Mantid
Loading...
Searching...
No Matches
IkedaCarpenterPV.cpp
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2018 ISIS Rutherford Appleton Laboratory UKRI,
4// NScD Oak Ridge National Laboratory, European Spallation Source,
5// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
6// SPDX - License - Identifier: GPL - 3.0 +
7//----------------------------------------------------------------------
8// Includes
9//----------------------------------------------------------------------
11#include "MantidAPI/Axis.h"
24
25#include <cmath>
26#include <gsl/gsl_math.h>
27#include <gsl/gsl_multifit_nlin.h>
28#include <gsl/gsl_sf_erf.h>
29#include <limits>
30
32
33using namespace CurveFitting;
34
35namespace {
37Kernel::Logger g_log("IkedaCarpenterPV");
38} // namespace
39
40using namespace Kernel;
41
42using namespace CurveFitting::SpecialFunctionSupport;
43
44using namespace Geometry;
45using namespace Constraints;
46
47DECLARE_FUNCTION(IkedaCarpenterPV)
48
49double IkedaCarpenterPV::centre() const { return getParameter("X0"); }
50
51void IkedaCarpenterPV::setHeight(const double h) {
52 // calculate height of peakshape function corresponding to intensity = 1
53 setParameter("I", 1);
54 double h0 = height();
55
56 // to avoid devide by zero and to insane value for I to be set
57 double minCutOff = 100.0 * std::numeric_limits<double>::min();
58 if (h0 > 0 && h0 < minCutOff)
59 h0 = minCutOff;
60 if (h0 < 0 && h0 > -minCutOff)
61 h0 = -minCutOff;
62
63 // The intensity is then estimated to be h/h0
64 setParameter("I", h / h0);
65}
66
68 // return the function value at centre()
69 // using arrays - otherwise coverity warning
70 double h0[1];
71 double toCentre[1];
72 toCentre[0] = centre();
73 constFunction(h0, toCentre, 1);
74 return h0[0];
75}
76
77double IkedaCarpenterPV::fwhm() const {
78 double sigmaSquared = getParameter("SigmaSquared");
79 double gamma = getParameter("Gamma");
80
81 if (sigmaSquared < 0) {
82 g_log.debug() << "SigmaSquared NEGATIVE!.\n"
83 << "Likely due to a fit not converging properly\n"
84 << "If this is frequent problem please report to Mantid team.\n"
85 << "For now to calculate width force SigmaSquared positive.\n";
86 sigmaSquared = -sigmaSquared;
87 }
88 if (gamma < 0) {
89 g_log.debug() << "Gamma NEGATIVE!.\n"
90 << "Likely due to a fit not converging properly\n"
91 << "If this is frequent problem please report to Mantid team.\n"
92 << "For now to calculate width force Gamma positive.\n";
93 gamma = -gamma;
94 }
95
96 // Use the same Thompson-Cox-Hastings approximation that the peak shape
97 // itself uses (via convertVoigtToPseudo) so that the reported FWHM is
98 // consistent with the actual profile width.
99 double H = 1.0, eta = 0.5;
100 convertVoigtToPseudo(sigmaSquared, gamma, H, eta);
101 return H;
102}
103
104void IkedaCarpenterPV::setFwhm(const double w) {
105 setParameter("SigmaSquared", w * w / (32.0 * M_LN2)); // used 4.0 * 8.0 = 32.0
106 setParameter("Gamma", w / 2.0);
107}
108
109void IkedaCarpenterPV::setCentre(const double c) { setParameter("X0", c); }
110
112 declareParameter("I", 0.0,
113 "The integrated intensity of the peak. I.e. "
114 "approximately equal to HWHM times height of "
115 "peak");
116 this->lowerConstraint0("I");
117 declareParameter("Alpha0", 1.6, "Used to model fast decay constant");
118 this->lowerConstraint0("Alpha0");
119 declareParameter("Alpha1", 1.5, "Used to model fast decay constant");
120 this->lowerConstraint0("Alpha1");
121 declareParameter("Beta0", 31.9, "Inverse of slow decay constant");
122 this->lowerConstraint0("Beta0");
123 declareParameter("Kappa", 46.0, "Controls contribution of slow decay term");
124 this->lowerConstraint0("Kappa");
125 declareParameter("SigmaSquared", 1.0, "standard deviation squared (Voigt Guassian broadening)");
126 this->lowerConstraint0("SigmaSquared");
127 declareParameter("Gamma", 1.0, "Voigt Lorentzian broadening");
128 this->lowerConstraint0("Gamma");
129 declareParameter("X0", 0.0, "Peak position");
130 this->lowerConstraint0("X0");
131}
132
133void IkedaCarpenterPV::lowerConstraint0(const std::string &paramName) {
134 auto mixingConstraint = std::make_unique<BoundaryConstraint>(this, paramName, 0.0, true);
135 mixingConstraint->setPenaltyFactor(1e9);
136
137 addConstraint(std::move(mixingConstraint));
138}
139
148void IkedaCarpenterPV::calWavelengthAtEachDataPoint(const double *xValues, const size_t &nData) const {
149 // if wavelength vector already have the right size no need for resizing it
150 // further we make the assumption that no need to recalculate this vector if
151 // it already has the right size
152
153 if (m_waveLength.size() != nData) {
154 m_waveLength.resize(nData);
155
156 Mantid::Kernel::Unit_sptr wavelength = Mantid::Kernel::UnitFactory::Instance().create("Wavelength");
157 for (size_t i = 0; i < nData; i++) {
158 m_waveLength[i] = xValues[i];
159 }
160
161 // note if a version of convertValue was added which allows a double* as
162 // first argument then could avoid copying above plus only have to resize
163 // m_wavelength when its size smaller than nData
165 if (mws) {
166 Instrument_const_sptr instrument = mws->getInstrument();
167 Geometry::IComponent_const_sptr sample = instrument->getSample();
168 if (sample != nullptr) {
169 convertValue(m_waveLength, wavelength, mws, m_workspaceIndex);
170 } else {
171 g_log.warning() << "No sample set for instrument in workspace.\n"
172 << "Can't calculate wavelength in IkedaCarpenter.\n"
173 << "Default all wavelengths to one.\n"
174 << "Solution is to load appropriate instrument into workspace.\n";
175 for (size_t i = 0; i < nData; i++)
176 m_waveLength[i] = 1.0;
177 }
178 } else {
179 g_log.warning() << "Workspace not set.\n"
180 << "Can't calculate wavelength in IkedaCarpenter.\n"
181 << "Default all wavelengths to one.\n"
182 << "Solution call setMatrixWorkspace() for function.\n";
183 for (size_t i = 0; i < nData; i++)
184 m_waveLength[i] = 1.0;
185 }
186 }
187}
188
196void IkedaCarpenterPV::convertVoigtToPseudo(const double &voigtSigmaSq, const double &voigtGamma, double &H,
197 double &eta) const {
198 double fwhmGsq = 8.0 * M_LN2 * voigtSigmaSq;
199 double fwhmG = sqrt(fwhmGsq);
200 double fwhmG4 = fwhmGsq * fwhmGsq;
201 double fwhmL = voigtGamma;
202 double fwhmLsq = voigtGamma * voigtGamma;
203 double fwhmL4 = fwhmLsq * fwhmLsq;
204
205 H = pow(fwhmG4 * fwhmG + 2.69269 * fwhmG4 * fwhmL + 2.42843 * fwhmGsq * fwhmG * fwhmLsq +
206 4.47163 * fwhmGsq * fwhmLsq * fwhmL + 0.07842 * fwhmG * fwhmL4 + fwhmL4 * fwhmL,
207 0.2);
208
209 if (H == 0.0)
210 H = std::numeric_limits<double>::epsilon() * 1000.0;
211
212 double tmp = fwhmL / H;
213
214 eta = 1.36603 * tmp - 0.47719 * tmp * tmp + 0.11116 * tmp * tmp * tmp;
215}
216
217void IkedaCarpenterPV::constFunction(double *out, const double *xValues, const int &nData) const {
218 // Save the wavelength cache so that single-point evaluations (e.g. from
219 // height()) do not invalidate the full-spectrum cache used by functionLocal.
220 auto savedWaveLength = std::move(m_waveLength);
221 m_waveLength.clear();
222 functionLocal(out, xValues, static_cast<size_t>(nData));
223 m_waveLength = std::move(savedWaveLength);
224}
225
226void IkedaCarpenterPV::functionLocal(double *out, const double *xValues, const size_t nData) const {
227 const double I = getParameter("I");
228 const double alpha0 = getParameter("Alpha0");
229 const double alpha1 = getParameter("Alpha1");
230 const double beta0 = getParameter("Beta0");
231 const double kappa = getParameter("Kappa");
232 const double voigtsigmaSquared = getParameter("SigmaSquared");
233 const double voigtgamma = getParameter("Gamma");
234 const double X0 = getParameter("X0");
235
236 // cal pseudo voigt sigmaSq and gamma and eta
237 double gamma = 1.0; // dummy initialization
238 double eta = 0.5; // dummy initialization
239 convertVoigtToPseudo(voigtsigmaSquared, voigtgamma, gamma, eta);
240 const double eps = std::numeric_limits<double>::epsilon() * 1000.0;
241 double sigmaSquared = std::max(gamma * gamma / (8.0 * M_LN2), eps); // pseudo voigt sigma^2
242
243 const double beta = 1 / beta0;
244
245 // equations taken from Fullprof manual
246
247 const double k = 0.05;
248
249 const double someConst = 1.0 / sqrt(2.0 * sigmaSquared);
250
251 // update wavelength vector
252 calWavelengthAtEachDataPoint(xValues, nData);
253
254 for (size_t i = 0; i < nData; i++) {
255 double diff = xValues[i] - X0;
256
257 double R = exp(-81.799 / (m_waveLength[i] * m_waveLength[i] * kappa));
258 double alpha = 1.0 / (alpha0 + m_waveLength[i] * alpha1);
259
260 double a_minus = alpha * (1 - k);
261 double a_plus = alpha * (1 + k);
262 double x = a_minus - beta;
263 double y = alpha - beta;
264 double z = a_plus - beta;
265
266 double Nu = 1 - R * a_minus / x;
267 double Nv = 1 - R * a_plus / z;
268 double Ns = -2 * (1 - R * alpha / y);
269 double Nr = 2 * R * alpha * alpha * beta * k * k / (x * y * z);
270
271 double u = a_minus * (a_minus * sigmaSquared - 2 * diff) / 2.0;
272 double v = a_plus * (a_plus * sigmaSquared - 2 * diff) / 2.0;
273 double s = alpha * (alpha * sigmaSquared - 2 * diff) / 2.0;
274 double r = beta * (beta * sigmaSquared - 2 * diff) / 2.0;
275
276 double yu = (a_minus * sigmaSquared - diff) * someConst;
277 double yv = (a_plus * sigmaSquared - diff) * someConst;
278 double ys = (alpha * sigmaSquared - diff) * someConst;
279 double yr = (beta * sigmaSquared - diff) * someConst;
280
281 std::complex<double> zs = std::complex<double>(-alpha * diff, 0.5 * alpha * gamma);
282 std::complex<double> zu = (1 - k) * zs;
283 std::complex<double> zv = (1 + k) * zs;
284 std::complex<double> zr = std::complex<double>(-beta * diff, 0.5 * beta * gamma);
285
286 double N = 0.25 * alpha * (1 - k * k) / (k * k);
287
288 out[i] = I * N *
289 ((1 - eta) * (Nu * exp(u + gsl_sf_log_erfc(yu)) + Nv * exp(v + gsl_sf_log_erfc(yv)) +
290 Ns * exp(s + gsl_sf_log_erfc(ys)) + Nr * exp(r + gsl_sf_log_erfc(yr))) -
291 eta * 2.0 / M_PI *
292 (Nu * exponentialIntegral(zu).imag() + Nv * exponentialIntegral(zv).imag() +
293 Ns * exponentialIntegral(zs).imag() + Nr * exponentialIntegral(zr).imag()));
294 }
295}
296
297void IkedaCarpenterPV::functionDerivLocal(API::Jacobian *jacobian, const double *xValues, const size_t nData) {
298 // Parameter order: 0:I 1:Alpha0 2:Alpha1 3:Beta0 4:Kappa 5:SigmaSquared 6:Gamma 7:X0
299 // For more details on the origin of this derivative see
300 // docs/source/fitting/fitfunctions/IkedaCarpenterPV_analytical_derivative.rst
301
302 const double I = getParameter("I");
303 const double alpha0 = getParameter("Alpha0");
304 const double alpha1 = getParameter("Alpha1");
305 const double beta0 = getParameter("Beta0");
306 const double kappa = getParameter("Kappa");
307 const double voigtSigmaSq = getParameter("SigmaSquared");
308 const double voigtGamma = getParameter("Gamma");
309 const double X0 = getParameter("X0");
310
311 double H = 1.0, eta = 0.5;
312 convertVoigtToPseudo(voigtSigmaSq, voigtGamma, H, eta);
313
314 const double eps = std::numeric_limits<double>::epsilon() * 1000.0;
315 const double Hsafe = std::max(H, eps);
316 const double Hsafe2 = Hsafe * Hsafe;
317 const double S = std::max(Hsafe2 / (8.0 * M_LN2), eps);
318 const double beta = 1.0 / beta0;
319 const double beta2 = beta * beta;
320
321 constexpr double k = 0.05;
322 constexpr double kk = k * k;
323 constexpr double one_minus_k = 1.0 - k;
324 constexpr double one_plus_k = 1.0 + k;
325 const double one_minus_eta = 1.0 - eta;
326 constexpr double twoOverPi = 2.0 / M_PI;
327 const double etaTwoOverPi = eta * twoOverPi;
328 const double halfHsafe = 0.5 * Hsafe;
329 const double halfH2 = halfHsafe * halfHsafe;
330
331 const double invRoot2S = 1.0 / std::sqrt(2.0 * S);
332 const double rootSOver2 = std::sqrt(S / 2.0);
333 const double twoInvSqrtPi = 2.0 / std::sqrt(M_PI);
334 const double CexpSFactor = invRoot2S / (2.0 * S); // = 1/(2*S*sqrt(2S)), for dK/dS
335
336 // Voigt -> pseudo-Voigt chain derivatives
337 const double g2 = std::max(8.0 * M_LN2 * voigtSigmaSq, eps);
338 const double g = std::sqrt(g2);
339 const double g3 = g2 * g;
340 const double g4 = g2 * g2;
341
342 const double L = voigtGamma;
343 const double L2 = L * L;
344 const double L3 = L2 * L;
345 const double L4 = L2 * L2;
346
347 constexpr double c1 = 2.69269, c2 = 2.42843, c3 = 4.47163, c4 = 0.07842;
348
349 const double H4 = Hsafe2 * Hsafe2;
350 const double inv5H4 = 1.0 / (5.0 * std::max(H4, eps));
351
352 const double dPdg = 5.0 * g4 + 4.0 * c1 * g3 * L + 3.0 * c2 * g2 * L2 + 2.0 * c3 * g * L3 + c4 * L4;
353 const double dPdL = c1 * g4 + 2.0 * c2 * g3 * L + 3.0 * c3 * g2 * L2 + 4.0 * c4 * g * L3 + 5.0 * L4;
354
355 const double dgdV = (4.0 * M_LN2) / g;
356 const double dHdV = dPdg * dgdV * inv5H4;
357 const double dHdL = dPdL * inv5H4;
358
359 const double t = L / Hsafe;
360 const double deta_dt = 1.36603 - 0.95438 * t + 0.33348 * t * t;
361 const double invHsafe2 = 1.0 / Hsafe2;
362
363 const double deta_dV = deta_dt * (-L * dHdV * invHsafe2);
364 const double deta_dL = deta_dt * (1.0 / Hsafe - L * dHdL * invHsafe2);
365
366 const double inv4Ln2 = 1.0 / (4.0 * M_LN2);
367 const double dSdV = Hsafe * dHdV * inv4Ln2;
368 const double dSdL = Hsafe * dHdL * inv4Ln2;
369
370 const double kappa2_inv = 1.0 / (kappa * kappa);
371
372 calWavelengthAtEachDataPoint(xValues, nData);
373
374 for (size_t i = 0; i < nData; ++i) {
375 const double lambda = m_waveLength[i];
376 const double lambda2 = lambda * lambda;
377 const double diff = xValues[i] - X0;
378 const double diff2 = diff * diff;
379
380 const double R = std::exp(-81.799 / (lambda2 * kappa));
381 const double alpha = 1.0 / (alpha0 + lambda * alpha1);
382 const double alpha2 = alpha * alpha;
383
384 const double a_minus = alpha * one_minus_k;
385 const double a_plus = alpha * one_plus_k;
386
387 const double x = a_minus - beta;
388 const double y = alpha - beta;
389 const double z = a_plus - beta;
390
391 const double inv_x = 1.0 / x;
392 const double inv_y = 1.0 / y;
393 const double inv_z = 1.0 / z;
394 const double inv_x2 = inv_x * inv_x;
395 const double inv_y2 = inv_y * inv_y;
396 const double inv_z2 = inv_z * inv_z;
397
398 const double Nu = 1.0 - R * a_minus * inv_x;
399 const double Nv = 1.0 - R * a_plus * inv_z;
400 const double Ns = -2.0 * (1.0 - R * alpha * inv_y);
401 const double Nr = 2.0 * R * alpha2 * beta * kk * inv_x * inv_y * inv_z;
402
403 // Gaussian kernel values: K = exp(u + log_erfc(y))
404 const double amuS_d = a_minus * S - diff; // a_minus*S - diff (reused in dK)
405 const double apuS_d = a_plus * S - diff;
406 const double alS_d = alpha * S - diff;
407 const double beS_d = beta * S - diff;
408
409 const double u = 0.5 * a_minus * (amuS_d - diff);
410 const double v = 0.5 * a_plus * (apuS_d - diff);
411 const double s = 0.5 * alpha * (alS_d - diff);
412 const double r = 0.5 * beta * (beS_d - diff);
413
414 const double yu = amuS_d * invRoot2S;
415 const double yv = apuS_d * invRoot2S;
416 const double ys = alS_d * invRoot2S;
417 const double yr = beS_d * invRoot2S;
418
419 const double Ku = std::exp(u + gsl_sf_log_erfc(yu));
420 const double Kv = std::exp(v + gsl_sf_log_erfc(yv));
421 const double Ks = std::exp(s + gsl_sf_log_erfc(ys));
422 const double Kr = std::exp(r + gsl_sf_log_erfc(yr));
423
424 // Lorentzian kernel values
425 const std::complex<double> zu(-a_minus * diff, a_minus * halfHsafe);
426 const std::complex<double> zv(-a_plus * diff, a_plus * halfHsafe);
427 const std::complex<double> zs(-alpha * diff, alpha * halfHsafe);
428 const std::complex<double> zr(-beta * diff, beta * halfHsafe);
429
430 const std::complex<double> Fu = exponentialIntegral(zu);
431 const std::complex<double> Fv = exponentialIntegral(zv);
432 const std::complex<double> Fs = exponentialIntegral(zs);
433 const std::complex<double> Fr = exponentialIntegral(zr);
434
435 const double Ju = Fu.imag();
436 const double Jv = Fv.imag();
437 const double Js = Fs.imag();
438 const double Jr = Fr.imag();
439
440 // F'(z) = F(z) - 1/z, decomposed into real/imag to avoid complex division.
441 // For z = (-a*diff, a*halfH): 1/z = (-diff, -halfH) / (a*(diff^2+halfH^2))
442 // So F'_r = F_r + diff/(a*d2h2), F'_i = F_i + halfH/(a*d2h2)
443 const double inv_d2h2 = 1.0 / (diff2 + halfH2);
444 const double d_term = diff * inv_d2h2; // diff / (diff^2 + halfH^2)
445 const double h_term = halfHsafe * inv_d2h2; // halfH / (diff^2 + halfH^2)
446
447 const double Fpu_r = Fu.real() + d_term / a_minus;
448 const double Fpu_i = Fu.imag() + h_term / a_minus;
449 const double Fpv_r = Fv.real() + d_term / a_plus;
450 const double Fpv_i = Fv.imag() + h_term / a_plus;
451 const double Fps_r = Fs.real() + d_term / alpha;
452 const double Fps_i = Fs.imag() + h_term / alpha;
453 const double Fpr_r = Fr.real() + d_term / beta;
454 const double Fpr_i = Fr.imag() + h_term / beta;
455
456 // Function value components
457 const double Q = 0.25 * alpha * (1.0 - kk) / kk;
458 const double IQ = I * Q;
459
460 const double G = Nu * Ku + Nv * Kv + Ns * Ks + Nr * Kr;
461 const double Lmix = Nu * Ju + Nv * Jv + Ns * Js + Nr * Jr;
462 const double modelCore = one_minus_eta * G - etaTwoOverPi * Lmix;
463
464 // Common Gaussian derivative factor: (2/sqrt(pi)) * exp(-diff^2/(2S))
465 const double Cexp = twoInvSqrtPi * std::exp(-diff2 / (2.0 * S));
466 const double CexpRS2 = Cexp * rootSOver2; // for dK/da
467
468 // -- Inlined K kernel derivatives --
469 // dK/da = (a*S - diff)*K - Cexp*sqrt(S/2)
470 const double dKu_da = amuS_d * Ku - CexpRS2;
471 const double dKv_da = apuS_d * Kv - CexpRS2;
472 const double dKs_da = alS_d * Ks - CexpRS2;
473 const double dKr_da = beS_d * Kr - CexpRS2;
474
475 // dK/dS = 0.5*a^2*K - Cexp*(a*S+diff)/(2*S*sqrt(2S))
476 const double dKu_dS = 0.5 * a_minus * a_minus * Ku - Cexp * (a_minus * S + diff) * CexpSFactor;
477 const double dKv_dS = 0.5 * a_plus * a_plus * Kv - Cexp * (a_plus * S + diff) * CexpSFactor;
478 const double dKs_dS = 0.5 * alpha2 * Ks - Cexp * (alpha * S + diff) * CexpSFactor;
479 const double dKr_dS = 0.5 * beta2 * Kr - Cexp * (beta * S + diff) * CexpSFactor;
480
481 // dK/dd = -a*K + Cexp/sqrt(2S)
482 const double CexpIR2S = Cexp * invRoot2S;
483 const double dKu_dd = -a_minus * Ku + CexpIR2S;
484 const double dKv_dd = -a_plus * Kv + CexpIR2S;
485 const double dKs_dd = -alpha * Ks + CexpIR2S;
486 const double dKr_dd = -beta * Kr + CexpIR2S;
487
488 // -- J kernel derivatives --
489 // dJ/da = Im(Fp * (-diff, halfH)) = Fp_r*halfH - Fp_i*diff
490 const double dJu_da = Fpu_r * halfHsafe - Fpu_i * diff;
491 const double dJv_da = Fpv_r * halfHsafe - Fpv_i * diff;
492 const double dJs_da = Fps_r * halfHsafe - Fps_i * diff;
493 const double dJr_da = Fpr_r * halfHsafe - Fpr_i * diff;
494
495 // dJ/dH = 0.5*a*Re(Fp)
496 const double dJu_dH = 0.5 * a_minus * Fpu_r;
497 const double dJv_dH = 0.5 * a_plus * Fpv_r;
498 const double dJs_dH = 0.5 * alpha * Fps_r;
499 const double dJr_dH = 0.5 * beta * Fpr_r;
500
501 // dJ/dd = -a*Im(Fp)
502 const double dJu_dd = -a_minus * Fpu_i;
503 const double dJv_dd = -a_plus * Fpv_i;
504 const double dJs_dd = -alpha * Fps_i;
505 const double dJr_dd = -beta * Fpr_i;
506
507 // -- Coefficient derivatives wrt alpha, beta, R --
508 const double Rbeta = R * beta;
509 const double dNu_dalpha = Rbeta * one_minus_k * inv_x2;
510 const double dNv_dalpha = Rbeta * one_plus_k * inv_z2;
511 const double dNs_dalpha = -2.0 * Rbeta * inv_y2;
512 const double dNr_dalpha = Nr * (2.0 / alpha - one_minus_k * inv_x - inv_y - one_plus_k * inv_z);
513
514 const double dNu_dbeta = -R * a_minus * inv_x2;
515 const double dNv_dbeta = -R * a_plus * inv_z2;
516 const double dNs_dbeta = 2.0 * R * alpha * inv_y2;
517 const double dNr_dbeta = Nr * (1.0 / beta + inv_x + inv_y + inv_z);
518
519 const double dNu_dR = -a_minus * inv_x;
520 const double dNv_dR = -a_plus * inv_z;
521 const double dNs_dR = 2.0 * alpha * inv_y;
522 const double dNr_dR = 2.0 * alpha2 * beta * kk * inv_x * inv_y * inv_z;
523
524 // -- Mix derivatives --
525 // dG/dalpha: Kr doesn't depend on alpha, so dKr_dalpha term is zero
526 const double dG_dalpha = dNu_dalpha * Ku + Nu * one_minus_k * dKu_da + dNv_dalpha * Kv + Nv * one_plus_k * dKv_da +
527 dNs_dalpha * Ks + Ns * dKs_da + dNr_dalpha * Kr;
528
529 const double dL_dalpha = dNu_dalpha * Ju + Nu * one_minus_k * dJu_da + dNv_dalpha * Jv + Nv * one_plus_k * dJv_da +
530 dNs_dalpha * Js + Ns * dJs_da + dNr_dalpha * Jr;
531
532 // dG/dbeta: Ku,Kv,Ks don't depend on beta, so only Nr*dKr_da survives
533 const double dG_dbeta = dNu_dbeta * Ku + dNv_dbeta * Kv + dNs_dbeta * Ks + dNr_dbeta * Kr + Nr * dKr_da;
534
535 const double dL_dbeta = dNu_dbeta * Ju + dNv_dbeta * Jv + dNs_dbeta * Js + dNr_dbeta * Jr + Nr * dJr_da;
536
537 // dG/dR, dL/dR: R doesn't affect kernel arguments, only coefficients
538 const double dG_dR = dNu_dR * Ku + dNv_dR * Kv + dNs_dR * Ks + dNr_dR * Kr;
539 const double dL_dR = dNu_dR * Ju + dNv_dR * Jv + dNs_dR * Js + dNr_dR * Jr;
540
541 // dG/dS (K depends on S; J does not, so dL/dS = 0)
542 const double dG_dS = Nu * dKu_dS + Nv * dKv_dS + Ns * dKs_dS + Nr * dKr_dS;
543
544 // dL/dH (J depends on H; K does not, so dG/dH = 0)
545 const double dL_dH = Nu * dJu_dH + Nv * dJv_dH + Ns * dJs_dH + Nr * dJr_dH;
546
547 // dG/dd, dL/dd
548 const double dG_dd = Nu * dKu_dd + Nv * dKv_dd + Ns * dKs_dd + Nr * dKr_dd;
549 const double dL_dd = Nu * dJu_dd + Nv * dJv_dd + Ns * dJs_dd + Nr * dJr_dd;
550
551 // -- Derivatives wrt natural variables (noting dG_dH=0, dL_dS=0) --
552 const double df_dI = Q * modelCore;
553 const double df_dalpha = I * (Q / alpha * modelCore + Q * (one_minus_eta * dG_dalpha - etaTwoOverPi * dL_dalpha));
554 const double df_dbeta = IQ * (one_minus_eta * dG_dbeta - etaTwoOverPi * dL_dbeta);
555 const double df_dR = IQ * (one_minus_eta * dG_dR - etaTwoOverPi * dL_dR);
556 const double df_dS = IQ * one_minus_eta * dG_dS; // dL_dS = 0
557 const double df_dH = -IQ * etaTwoOverPi * dL_dH; // dG_dH = 0
558 const double df_deta = -IQ * (G + twoOverPi * Lmix);
559 const double df_ddiff = IQ * (one_minus_eta * dG_dd - etaTwoOverPi * dL_dd);
560
561 // -- Chain rule to Mantid fit parameters --
562 const double neg_alpha2 = -alpha2;
563 jacobian->set(i, 0, df_dI);
564 jacobian->set(i, 1, df_dalpha * neg_alpha2);
565 jacobian->set(i, 2, df_dalpha * neg_alpha2 * lambda);
566 jacobian->set(i, 3, df_dbeta * (-beta2));
567 jacobian->set(i, 4, df_dR * R * 81.799 * kappa2_inv / lambda2);
568 jacobian->set(i, 5, df_dS * dSdV + df_dH * dHdV + df_deta * deta_dV);
569 jacobian->set(i, 6, df_dS * dSdL + df_dH * dHdL + df_deta * deta_dL);
570 jacobian->set(i, 7, -df_ddiff);
571 }
572}
573
575 const auto *domain1D = dynamic_cast<const API::FunctionDomain1D *>(&domain);
576 if (!domain1D) {
577 calNumericalDeriv(domain, jacobian);
578 return;
579 }
580 functionDerivLocal(&jacobian, domain1D->getPointerAt(0), domain1D->size());
581}
582
585 auto interval = getDomainInterval(1e-2);
586
588 API::IntegrationResult result = integrator.integrate(*this, interval.first, interval.second);
589
590 return result.result;
591}
592
593void IkedaCarpenterPV::setMatrixWorkspace(std::shared_ptr<const API::MatrixWorkspace> workspace, size_t wi,
594 double startX, double endX) {
596
597 // Invalidate cached wavelengths so they are recalculated on the next
598 // function evaluation using the newly-set workspace and workspace index.
599 m_waveLength.clear();
600
601 if (workspace) {
602 // If the instrument definition provided IC parameters the base-class call above handles
603 // the formula evaluation and unit conversion.
604 if (isExplicitlySet(parameterIndex("Alpha0")))
605 return;
606
607 // This next section handles conversion when no instrument defaults are present
608 // in such a scenario the generic defaults will be used, which are in TOF units.
609
610 // This will scale these default time-dependent IC parameters so the initial plot guess
611 // is reasonable for non-TOF workspaces (e.g. d-spacing, wavelength).
612
613 // Only applies to units with a simple power-law relationship with TOF;
614 // units like DeltaE have a non-linear, energy-dependent conversion that
615 // may be kinematically forbidden at the peak position.
616 const auto wsUnit = workspace->getAxis(0)->unit();
617 auto tof = Mantid::Kernel::UnitFactory::Instance().create("TOF");
618 double factor, power;
619 if (!wsUnit->quickConversion(*tof, factor, power))
620 return;
621
622 const auto peakCentre = getParameter("X0");
623 if (peakCentre != 0.0) {
624 const auto tofCentre = convertValue(peakCentre, tof, workspace, wi);
625 if (tofCentre != 0.0) {
626 const double scaleFactor = peakCentre / tofCentre;
627 if (std::abs(scaleFactor - 1.0) > 1e-6) {
628 const auto scaleIfDefault = [&](const std::string &name, const double factor) {
629 const size_t idx = parameterIndex(name);
630 if (!isExplicitlySet(idx))
631 setParameter(idx, getParameter(idx) * factor, false);
632 };
633 scaleIfDefault("Alpha0", scaleFactor);
634 scaleIfDefault("Alpha1", scaleFactor);
635 scaleIfDefault("Beta0", scaleFactor);
636 scaleIfDefault("SigmaSquared", scaleFactor * scaleFactor);
637 scaleIfDefault("Gamma", scaleFactor);
638 }
639 }
640 }
641 }
642}
643
644} // namespace Mantid::CurveFitting::Functions
const std::vector< double > * lambda
gsl_vector * tmp
#define DECLARE_FUNCTION(classname)
Macro for declaring a new type of function to be used with the FunctionFactory.
IPeaksWorkspace_sptr workspace
Represent a domain for functions of one real argument.
Base class that represents the domain of a function.
size_t m_workspaceIndex
An index to a spectrum.
Definition IFunctionMW.h:44
void setMatrixWorkspace(std::shared_ptr< const API::MatrixWorkspace > workspace, size_t wi, double startX, double endX) override
Set MatrixWorkspace.
std::shared_ptr< const API::MatrixWorkspace > getMatrixWorkspace() const
Get shared pointer to the workspace.
virtual double getParameter(size_t i) const =0
Get i-th parameter.
double convertValue(double value, Kernel::Unit_sptr &outUnit, const std::shared_ptr< const MatrixWorkspace > &ws, size_t wsIndex) const
Convert a value from one unit (inUnit) to unit defined in workspace (ws)
void calNumericalDeriv(const FunctionDomain &domain, Jacobian &jacobian)
Calculate numerical derivatives.
virtual size_t parameterIndex(const std::string &name) const =0
Returns the index of parameter name.
virtual void addConstraint(std::unique_ptr< IConstraint > ic)
Add a constraint to function.
virtual bool isExplicitlySet(size_t i) const =0
Checks if a parameter has been set explicitly.
virtual void declareParameter(const std::string &name, double initValue=0, const std::string &description="")=0
Declare a new parameter.
void setParameter(const std::string &name, const double &value, bool explicitlySet=true) override
Override parent so that we may bust the cache when a parameter is set.
virtual std::pair< double, double > getDomainInterval(double level=DEFAULT_SEARCH_LEVEL) const
Get the interval on which the peak has all its values above a certain level.
Represents the Jacobian in IFitFunction::functionDeriv.
Definition Jacobian.h:22
virtual void set(size_t iY, size_t iP, double value)=0
Set a value to a Jacobian matrix element.
IntegrationResult integrate(const IPeakFunction &peakFunction, double lowerLimit, double upperLimit) const
Integration of peak function on the interval [a, b].
Provide Ikeda-Carpenter-pseudo-Voigt peak shape function interface to IPeakFunction.
void functionLocal(double *out, const double *xValues, const size_t nData) const override
Function evaluation method to be implemented in the inherited classes.
void constFunction(double *out, const double *xValues, const int &nData) const
calculate the const function
void setHeight(const double h) override
Sets the parameters such that height == h.
void setCentre(const double c) override
Sets the parameters such that centre == c.
double intensity() const override
Returns the integral intensity of the peak.
double centre() const override
overwrite IPeakFunction base class methods
void init() override
overwrite IFunction base class method, which declare function parameters
void setMatrixWorkspace(std::shared_ptr< const API::MatrixWorkspace > workspace, size_t wi, double startX, double endX) override
Set matrix workspace.
void setFwhm(const double w) override
Sets the parameters such that FWHM = w.
void lowerConstraint0(const std::string &paramName)
constrain all parameters to be non-negative
std::vector< double > m_waveLength
container for storing wavelength values for each data point
double height() const override
Returns the height of the function.
void convertVoigtToPseudo(const double &voigtSigmaSq, const double &voigtGamma, double &H, double &eta) const
convert voigt params to pseudo voigt params
void functionDerivLocal(API::Jacobian *out, const double *xValues, const size_t nData) override
Derivative evaluation method. Default is to calculate numerically.
std::string name() const override
overwrite IFunction base class methods
void functionDeriv(const API::FunctionDomain &domain, API::Jacobian &jacobian) override
Derivatives of function with respect to active parameters.
double fwhm() const override
Returns the peak FWHM.
void calWavelengthAtEachDataPoint(const double *xValues, const size_t &nData) const
method for updating m_waveLength
void debug(const std::string &msg)
Logs at debug level.
Definition Logger.cpp:145
void warning(const std::string &msg)
Logs at warning level.
Definition Logger.cpp:117
Kernel::Logger g_log("ExperimentInfo")
static logger object
std::shared_ptr< const MatrixWorkspace > MatrixWorkspace_const_sptr
shared pointer to the matrix workspace base class (const version)
std::complex< double > MANTID_CURVEFITTING_DLL exponentialIntegral(const std::complex< double > &z)
Compute exp(z)*E1(z) where z is complex and E1(z) is the Exponential Integral.
std::shared_ptr< const IComponent > IComponent_const_sptr
Typdef of a shared pointer to a const IComponent.
Definition IComponent.h:167
std::shared_ptr< const Instrument > Instrument_const_sptr
Shared pointer to an const instrument object.
std::shared_ptr< Unit > Unit_sptr
Shared pointer to the Unit base class.
Definition Unit.h:194