-
Notifications
You must be signed in to change notification settings - Fork 122
/
BackToBackExponential.cpp
152 lines (132 loc) · 4.27 KB
/
BackToBackExponential.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
//----------------------------------------------------------------------
// Includes
//----------------------------------------------------------------------
#include "MantidCurveFitting/Functions/BackToBackExponential.h"
#include "MantidAPI/FunctionFactory.h"
#include <gsl/gsl_sf_erf.h>
#include <gsl/gsl_multifit_nlin.h>
#include <cmath>
#include <limits>
namespace Mantid {
namespace CurveFitting {
namespace Functions {
using namespace CurveFitting;
using namespace Kernel;
using namespace API;
DECLARE_FUNCTION(BackToBackExponential)
void BackToBackExponential::init() {
// Do not change the order of these parameters!
declareParameter("I", 0.0, "integrated intensity of the peak"); // 0
declareParameter("A", 1.0,
"exponential constant of rising part of neutron pulse"); // 1
declareParameter(
"B", 0.05, "exponential constant of decaying part of neutron pulse"); // 2
declareParameter("X0", 0.0, "peak position"); // 3
declareParameter(
"S", 1.0,
"standard deviation of gaussian part of peakshape function"); // 4
}
/**
* Get approximate height of the peak: function value at X0.
*/
double BackToBackExponential::height() const {
double x0 = getParameter(3);
std::vector<double> vec(1, x0);
FunctionDomain1DVector domain(vec);
FunctionValues values(domain);
function(domain, values);
return values[0];
}
/**
* Set new height of the peak. This method does this approximately.
* @param h :: New value for the height.
*/
void BackToBackExponential::setHeight(const double h) {
double h0 = height();
if (h0 == 0.0) {
setParameter(0, 1e-6);
h0 = height();
}
double area = getParameter(0); // == I
area *= h / h0;
if (area <= 0.0) {
area = 1e-6;
}
if (!std::isfinite(area)) {
area = std::numeric_limits<double>::max() / 2;
}
setParameter(0, area);
}
/**
* Get approximate peak width.
*/
double BackToBackExponential::fwhm() const { return 2 * getParameter("S"); }
/**
* Set new peak width approximately.
* @param w :: New value for the width.
*/
void BackToBackExponential::setFwhm(const double w) {
setParameter("S", w / 2.0);
}
void BackToBackExponential::function1D(double *out, const double *xValues,
const size_t nData) const {
/*
const double& I = getParameter("I");
const double& a = getParameter("A");
const double& b = getParameter("B");
const double& x0 = getParameter("X0");
const double& s = getParameter("S");
*/
const double I = getParameter(0);
const double a = getParameter(1);
const double b = getParameter(2);
const double x0 = getParameter(3);
const double s = getParameter(4);
// find the reasonable extent of the peak ~100 fwhm
double extent = expWidth();
if (s > extent)
extent = s;
extent *= 100;
double s2 = s * s;
double normFactor = a * b / (a + b) / 2;
// Needed for IntegratePeaksMD for cylinder profile fitted with b=0
if (normFactor == 0.0)
normFactor = 1.0;
for (size_t i = 0; i < nData; i++) {
double diff = xValues[i] - x0;
if (fabs(diff) < extent) {
double val = 0.0;
double arg1 = a / 2 * (a * s2 + 2 * diff);
val += exp(arg1 + gsl_sf_log_erfc((a * s2 + diff) /
sqrt(2 * s2))); // prevent overflow
double arg2 = b / 2 * (b * s2 - 2 * diff);
val += exp(arg2 + gsl_sf_log_erfc((b * s2 - diff) /
sqrt(2 * s2))); // prevent overflow
out[i] = I * val * normFactor;
} else
out[i] = 0.0;
}
}
/**
* Evaluate function derivatives numerically.
*/
void BackToBackExponential::functionDeriv1D(Jacobian *jacobian,
const double *xValues,
const size_t nData) {
FunctionDomain1DView domain(xValues, nData);
this->calNumericalDeriv(domain, *jacobian);
}
/**
* Calculate contribution to the width by the exponentials.
*/
double BackToBackExponential::expWidth() const {
const double a = getParameter(1);
const double b = getParameter(2);
// Needed for IntegratePeaksMD for cylinder profile fitted with b=0
if (a * b == 0.0)
return M_LN2;
return M_LN2 * (a + b) / (a * b);
}
} // namespace Functions
} // namespace CurveFitting
} // namespace Mantid