Mantid
Loading...
Searching...
No Matches
ThermalNeutronDtoTOFFunction.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 +
10#include <cmath>
11#include <gsl/gsl_sf_erf.h>
12
13using namespace Mantid::API;
14
15using namespace std;
16
18
19using namespace CurveFitting;
20
21//----------------------------------------------------------------------------------------------
22DECLARE_FUNCTION(ThermalNeutronDtoTOFFunction)
23
24
28
30 declareParameter("Dtt1", 1.0); // 0
31 declareParameter("Dtt1t", 1.0); // 1
32 declareParameter("Dtt2t", 1.0); // 2
33 declareParameter("Zero", 0.0); // 3
34 declareParameter("Zerot", 0.0); // 4
35
36 declareParameter("Width", 1.0); // 5
37 declareParameter("Tcross", 1.0); // 6
38}
39
43void ThermalNeutronDtoTOFFunction::function1D(double *out, const double *xValues, const size_t nData) const {
44 double dtt1 = getParameter("Dtt1");
45 double dtt1t = getParameter("Dtt1t");
46 double dtt2t = getParameter("Dtt2t");
47 double zero = getParameter("Zero");
48 double zerot = getParameter("Zerot");
49 double width = getParameter("Width");
50 double tcross = getParameter("Tcross");
51
52 for (size_t i = 0; i < nData; ++i) {
53 // out[i] = corefunction(xValues[i], dtt1, dtt1t, dtt2t, zero, zerot, width,
54 // tcross);
55 out[i] = calThermalNeutronTOF(xValues[i], dtt1, dtt1t, dtt2t, zero, zerot, width, tcross);
56 }
57}
58
59//------------------------------------------------------------------------------------------------
63void ThermalNeutronDtoTOFFunction::function1D(vector<double> &out, const vector<double> &xValues) const {
64 double dtt1 = getParameter(0);
65 double dtt1t = getParameter(1);
66 double dtt2t = getParameter(2);
67 double zero = getParameter(3);
68 double zerot = getParameter(4);
69 double width = getParameter(5);
70 double tcross = getParameter(6);
71
72 size_t nData = out.size();
73
74 for (size_t i = 0; i < nData; ++i) {
75 out[i] = calThermalNeutronTOF(xValues[i], dtt1, dtt1t, dtt2t, zero, zerot, width, tcross);
76 }
77}
78
79void ThermalNeutronDtoTOFFunction::functionDeriv1D(Jacobian *out, const double *xValues, const size_t nData) {
80 // 1. Get hold all parameters
81 const double dtt1 = getParameter("Dtt1");
82 const double dtt1t = getParameter("Dtt1t");
83 const double dtt2t = getParameter("Dtt2t");
84 const double zero = getParameter("Zero");
85 const double zerot = getParameter("Zerot");
86 const double width = getParameter("Width");
87 const double tcross = getParameter("Tcross");
88
89 // 2. Calcualtion
90 for (size_t i = 0; i < nData; ++i) {
91 // a) Some calcualtion
92 double x = xValues[i];
93 double n = 0.5 * gsl_sf_erfc(width * (tcross - 1 / x));
94 double u = width * (tcross - 1 / x);
95
96 double deriv_dtt1 = n * x;
97 double deriv_dtt1t = (1 - n) * x;
98 double deriv_dtt2t = (n - 1) / x;
99 double deriv_zero = n;
100 double deriv_zerot = (1 - n);
101 double deriv_width =
102 -(zero + dtt1 * x - zerot - dtt1t * x + dtt2t / x) * exp(-u * u) / sqrt(M_PI) * (tcross - 1 / x);
103 double deriv_tcross = -(zero + dtt1 * x - zerot - dtt1t * x + dtt2t / x) * exp(-u * u) / sqrt(M_PI) * width;
104
105 // b) Set
106 out->set(i, 0, deriv_dtt1);
107 out->set(i, 1, deriv_dtt1t);
108 out->set(i, 2, deriv_dtt2t);
109 out->set(i, 3, deriv_zero);
110 out->set(i, 4, deriv_zerot);
111 out->set(i, 5, deriv_width);
112 out->set(i, 6, deriv_tcross);
113 }
114}
115
118void ThermalNeutronDtoTOFFunction::functionDerivLocal(API::Jacobian * /*unused*/, const double * /*unused*/,
119 const size_t /*unused*/) {
120 throw Mantid::Kernel::Exception::NotImplementedError("functionDerivLocal is not implemented for "
121 "ThermalNeutronDtoTOFFunction.");
122}
123
124} // namespace Mantid::CurveFitting::Functions
#define DECLARE_FUNCTION(classname)
Macro for declaring a new type of function to be used with the FunctionFactory.
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.
double getParameter(size_t i) const override
Get i-th parameter.
void functionDeriv1D(API::Jacobian *out, const double *xValues, const size_t nData) override
Derivative.
void function1D(double *out, const double *xValues, const size_t nData) const override
Override.
void functionDerivLocal(API::Jacobian *, const double *, const size_t)
Derivative.
Marks code as not implemented yet.
Definition Exception.h:138
double calThermalNeutronTOF(double dh, double dtt1, double dtt1t, double dtt2t, double zero, double zerot, double width, double tcross)
Calcualte TOF from d-spacing value for thermal neutron.
STL namespace.