Mantid
Loading...
Searching...
No Matches
BackToBackExponential.h
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2007 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#pragma once
8
9//----------------------------------------------------------------------
10// Includes
11//----------------------------------------------------------------------
13#include "MantidCurveFitting/DllConfig.h"
14
15namespace Mantid {
16namespace CurveFitting {
17namespace Functions {
38class MANTID_CURVEFITTING_DLL BackToBackExponential : public API::IPeakFunction {
39public:
42
44 double centre() const override { return getParameter("X0"); }
45 void setCentre(const double c) override { setParameter("X0", c); }
46 double height() const override;
47 void setHeight(const double h) override;
48 double fwhm() const override;
49 void setFwhm(const double w) override;
50 double intensity() const override { return getParameter("I"); }
51 double intensityError() const override { return getError("I"); }
52 void setIntensity(const double newIntensity) override { setParameter("I", newIntensity); }
53 std::string getWidthParameterName() const override { return "S"; }
54
56 std::string name() const override { return "BackToBackExponential"; }
57 const std::string category() const override { return "Peak"; }
58 void function1D(double *out, const double *xValues, const size_t nData) const override;
59 void functionDeriv1D(API::Jacobian *jacobian, const double *xValues, const size_t nData) override;
60
61protected:
63 void init() override;
65 void functionLocal(double *, const double *, const size_t) const override {}
67 void functionDerivLocal(API::Jacobian *, const double *, const size_t) override {}
68 double expWidth() const;
69};
70
71using BackToBackExponential_sptr = std::shared_ptr<BackToBackExponential>;
72
73} // namespace Functions
74} // namespace CurveFitting
75} // namespace Mantid
double height
Definition: GetAllEi.cpp:155
An interface to a peak function, which extend the interface of IFunctionWithLocation by adding method...
Definition: IPeakFunction.h:51
Represents the Jacobian in IFitFunction::functionDeriv.
Definition: Jacobian.h:22
Provide BackToBackExponential peak shape function interface to IPeakFunction.
void functionDerivLocal(API::Jacobian *, const double *, const size_t) override
Derivative evaluation method to be implemented in the inherited classes.
void functionLocal(double *, const double *, const size_t) const override
Function evaluation method to be implemented in the inherited classes.
double intensity() const override
Returns the integral intensity of the peak.
std::string name() const override
overwrite IFunction base class methods
double centre() const override
overwrite IPeakFunction base class methods
double intensityError() const override
Error in the integrated intensity of the peak due to uncertainties in the values of the fit parameter...
std::string getWidthParameterName() const override
Get the name of the parameter that is changed when the fwhm is changed.
const std::string category() const override
The categories the Fit function belong to.
void setCentre(const double c) override
Sets the parameters such that centre == c.
void setIntensity(const double newIntensity) override
Sets the integral intensity of the peak.
std::shared_ptr< BackToBackExponential > BackToBackExponential_sptr
Helper class which provides the Collimation Length for SANS instruments.