Mantid
Loading...
Searching...
No Matches
Framework
CurveFitting
src
Functions
ThermalNeutronBk2BkExpBeta.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
#include "
MantidCurveFitting/Functions/ThermalNeutronBk2BkExpBeta.h
"
8
#include "
MantidAPI/FunctionFactory.h
"
9
10
#include <cmath>
11
#include <gsl/gsl_sf_erf.h>
12
13
using namespace
std
;
14
15
using namespace
Mantid
;
16
17
using namespace
Mantid::CurveFitting
;
18
19
using namespace
Mantid::API
;
20
21
namespace
Mantid::CurveFitting::Functions
{
22
23
using namespace
CurveFitting;
24
25
DECLARE_FUNCTION
(
ThermalNeutronBk2BkExpBeta
)
26
27
//----------------------------------------------------------------------------------------------
30
void
ThermalNeutronBk2BkExpBeta
::init() {
31
declareParameter(
"Width"
, 1.0);
32
declareParameter(
"Tcross"
, 1.0);
33
34
declareParameter(
"Beta0"
, 0.0);
35
declareParameter(
"Beta1"
, 0.0);
36
declareParameter(
"Beta0t"
, 0.0);
37
declareParameter(
"Beta1t"
, 0.0);
38
}
39
40
//----------------------------------------------------------------------------------------------
43
void
ThermalNeutronBk2BkExpBeta::function1D
(
double
*out,
const
double
*xValues,
const
size_t
nData)
const
{
44
double
width =
getParameter
(
"Width"
);
45
double
tcross =
getParameter
(
"Tcross"
);
46
double
beta0 =
getParameter
(
"Beta0"
);
47
double
beta1 =
getParameter
(
"Beta1"
);
48
double
beta0t =
getParameter
(
"Beta0t"
);
49
double
beta1t =
getParameter
(
"Beta1t"
);
50
51
for
(
size_t
i = 0; i < nData; ++i) {
52
out[i] =
corefunction
(xValues[i], width, tcross, beta0, beta1, beta0t, beta1t);
53
}
54
}
55
58
void
ThermalNeutronBk2BkExpBeta::functionDeriv
(
const
FunctionDomain
&domain,
Jacobian
&jacobian) {
59
calNumericalDeriv
(domain, jacobian);
60
}
61
62
//----------------------------------------------------------------------------------------------
65
double
ThermalNeutronBk2BkExpBeta::corefunction
(
double
dh,
double
width,
double
tcross,
double
beta0,
double
beta1,
66
double
beta0t,
double
beta1t)
const
{
67
double
n
= 0.5 * gsl_sf_erfc(width * (tcross - 1.0 / dh));
68
double
beta = 1.0 / (
n
* (beta0 + beta1 * dh) + (1.0 -
n
) * (beta0t - beta1t / dh));
69
70
return
beta;
71
}
72
73
}
// namespace Mantid::CurveFitting::Functions
n
size_t n
Definition
AugmentedLagrangianOptimizer.cpp:42
FunctionFactory.h
DECLARE_FUNCTION
#define DECLARE_FUNCTION(classname)
Macro for declaring a new type of function to be used with the FunctionFactory.
Definition
FunctionFactory.h:151
ThermalNeutronBk2BkExpBeta.h
Mantid::API::FunctionDomain
Base class that represents the domain of a function.
Definition
FunctionDomain.h:31
Mantid::API::IFunction::calNumericalDeriv
void calNumericalDeriv(const FunctionDomain &domain, Jacobian &jacobian)
Calculate numerical derivatives.
Definition
IFunction.cpp:1106
Mantid::API::Jacobian
Represents the Jacobian in IFitFunction::functionDeriv.
Definition
Jacobian.h:22
Mantid::API::ParamFunction::getParameter
double getParameter(size_t i) const override
Get i-th parameter.
Definition
ParamFunction.cpp:62
Mantid::CurveFitting::Functions::ThermalNeutronBk2BkExpBeta
ThermalNeutronBk2BkExpBETA : Function to calculate Beta of Bk2Bk Exponential function from Thermal Ne...
Definition
ThermalNeutronBk2BkExpBeta.h:21
Mantid::CurveFitting::Functions::ThermalNeutronBk2BkExpBeta::function1D
void function1D(double *out, const double *xValues, const size_t nData) const override
Override.
Definition
ThermalNeutronBk2BkExpBeta.cpp:43
Mantid::CurveFitting::Functions::ThermalNeutronBk2BkExpBeta::corefunction
double corefunction(double dh, double width, double tcross, double beta0, double beta1, double beta0t, double beta1t) const
Core function (inline) to calcualte TOF_h from d-spacing.
Definition
ThermalNeutronBk2BkExpBeta.cpp:65
Mantid::CurveFitting::Functions::ThermalNeutronBk2BkExpBeta::functionDeriv
void functionDeriv(const API::FunctionDomain &domain, API::Jacobian &jacobian) override
Derivative to overwrite.
Definition
ThermalNeutronBk2BkExpBeta.cpp:58
Mantid::API
Definition
AbsorptionCorrection.h:20
Mantid::CurveFitting::Functions
Definition
VesuvioCalculateGammaBackground.h:23
Mantid::CurveFitting
Definition
IFunction1D.h:19
Mantid
Helper class which provides the Collimation Length for SANS instruments.
Definition
AbsorptionCorrection.h:18
std
STL namespace.
Generated by
1.9.8