-
Notifications
You must be signed in to change notification settings - Fork 122
/
ThermalNeutronDtoTOFFunction.cpp
125 lines (106 loc) · 4.32 KB
/
ThermalNeutronDtoTOFFunction.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
// Mantid Repository : https://github.com/mantidproject/mantid
//
// Copyright © 2018 ISIS Rutherford Appleton Laboratory UKRI,
// NScD Oak Ridge National Laboratory, European Spallation Source,
// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
// SPDX - License - Identifier: GPL - 3.0 +
#include "MantidCurveFitting/Functions/ThermalNeutronDtoTOFFunction.h"
#include "MantidAPI/FunctionDomain1D.h"
#include "MantidAPI/FunctionFactory.h"
#include "MantidKernel/System.h"
#include <cmath>
#include <gsl/gsl_sf_erf.h>
using namespace Mantid::API;
using namespace std;
namespace Mantid::CurveFitting::Functions {
using namespace CurveFitting;
//----------------------------------------------------------------------------------------------
DECLARE_FUNCTION(ThermalNeutronDtoTOFFunction)
/**
* Define the fittable parameters
*/
void ThermalNeutronDtoTOFFunction::init() {
/// Instrument geometry related
declareParameter("Dtt1", 1.0); // 0
declareParameter("Dtt1t", 1.0); // 1
declareParameter("Dtt2t", 1.0); // 2
declareParameter("Zero", 0.0); // 3
declareParameter("Zerot", 0.0); // 4
declareParameter("Width", 1.0); // 5
declareParameter("Tcross", 1.0); // 6
}
/** Main function
* xValues containing the d-space value of peaks centres
*/
void ThermalNeutronDtoTOFFunction::function1D(double *out, const double *xValues, const size_t nData) const {
double dtt1 = getParameter("Dtt1");
double dtt1t = getParameter("Dtt1t");
double dtt2t = getParameter("Dtt2t");
double zero = getParameter("Zero");
double zerot = getParameter("Zerot");
double width = getParameter("Width");
double tcross = getParameter("Tcross");
for (size_t i = 0; i < nData; ++i) {
// out[i] = corefunction(xValues[i], dtt1, dtt1t, dtt2t, zero, zerot, width,
// tcross);
out[i] = calThermalNeutronTOF(xValues[i], dtt1, dtt1t, dtt2t, zero, zerot, width, tcross);
}
}
//------------------------------------------------------------------------------------------------
/** Main function
* xValues containing the d-space value of peaks centres
*/
void ThermalNeutronDtoTOFFunction::function1D(vector<double> &out, const vector<double> &xValues) const {
double dtt1 = getParameter(0);
double dtt1t = getParameter(1);
double dtt2t = getParameter(2);
double zero = getParameter(3);
double zerot = getParameter(4);
double width = getParameter(5);
double tcross = getParameter(6);
size_t nData = out.size();
for (size_t i = 0; i < nData; ++i) {
out[i] = calThermalNeutronTOF(xValues[i], dtt1, dtt1t, dtt2t, zero, zerot, width, tcross);
}
}
void ThermalNeutronDtoTOFFunction::functionDeriv1D(Jacobian *out, const double *xValues, const size_t nData) {
// 1. Get hold all parameters
const double dtt1 = getParameter("Dtt1");
const double dtt1t = getParameter("Dtt1t");
const double dtt2t = getParameter("Dtt2t");
const double zero = getParameter("Zero");
const double zerot = getParameter("Zerot");
const double width = getParameter("Width");
const double tcross = getParameter("Tcross");
// 2. Calcualtion
for (size_t i = 0; i < nData; ++i) {
// a) Some calcualtion
double x = xValues[i];
double n = 0.5 * gsl_sf_erfc(width * (tcross - 1 / x));
double u = width * (tcross - 1 / x);
double deriv_dtt1 = n * x;
double deriv_dtt1t = (1 - n) * x;
double deriv_dtt2t = (n - 1) / x;
double deriv_zero = n;
double deriv_zerot = (1 - n);
double deriv_width =
-(zero + dtt1 * x - zerot - dtt1t * x + dtt2t / x) * exp(-u * u) / sqrt(M_PI) * (tcross - 1 / x);
double deriv_tcross = -(zero + dtt1 * x - zerot - dtt1t * x + dtt2t / x) * exp(-u * u) / sqrt(M_PI) * width;
// b) Set
out->set(i, 0, deriv_dtt1);
out->set(i, 1, deriv_dtt1t);
out->set(i, 2, deriv_dtt2t);
out->set(i, 3, deriv_zero);
out->set(i, 4, deriv_zerot);
out->set(i, 5, deriv_width);
out->set(i, 6, deriv_tcross);
}
}
/** Some forbidden function
*/
void ThermalNeutronDtoTOFFunction::functionDerivLocal(API::Jacobian * /*unused*/, const double * /*unused*/,
const size_t /*unused*/) {
throw Mantid::Kernel::Exception::NotImplementedError("functionDerivLocal is not implemented for "
"ThermalNeutronDtoTOFFunction.");
}
} // namespace Mantid::CurveFitting::Functions