Mantid
Loading...
Searching...
No Matches
CubicSpline.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
14#include <boost/scoped_array.hpp>
15#include <gsl/gsl_errno.h>
16#include <gsl/gsl_spline.h>
17#include <valarray>
18
19namespace Mantid {
20namespace CurveFitting {
21namespace Functions {
31class MANTID_CURVEFITTING_DLL CubicSpline : public BackgroundFunction {
32
33public:
36
38 std::string name() const override { return "CubicSpline"; }
39 const std::string category() const override { return "Background"; }
40 void function1D(double *out, const double *xValues, const size_t nData) const override;
41 void derivative1D(double *out, const double *xValues, size_t nData, const size_t order) const override;
42 void setParameter(size_t i, const double &value, bool explicitlySet = true) override;
43 using ParamFunction::setParameter;
44
46 void setAttribute(const std::string &attName, const Attribute &) override;
47
49 void setXAttribute(const size_t index, double x);
50
51private:
53 const int m_min_points;
54
56 struct GSLFree {
57 void operator()(gsl_spline *spline) { gsl_spline_free(spline); }
58 void operator()(gsl_interp_accel *acc) { gsl_interp_accel_free(acc); }
59 } m_gslFree;
60
62 std::shared_ptr<gsl_interp_accel> m_acc;
63
65 std::shared_ptr<gsl_spline> m_spline;
66
68 mutable bool m_recalculateSpline;
69
71 void reallocGSLObjects(const int n);
72
74 void setupInput(boost::scoped_array<double> &x, boost::scoped_array<double> &y, int n) const;
75
77 void calculateSpline(double *out, const double *xValues, const size_t nData) const;
78
80 void calculateDerivative(double *out, const double *xValues, const size_t nData, const size_t order) const;
81
83 void initGSLObjects(boost::scoped_array<double> &x, boost::scoped_array<double> &y, int n) const;
84
86 void checkGSLError(const int status, const int errorType) const;
87
89 bool checkXInRange(double x) const;
90
92 double splineEval(const double x) const;
93};
94
95using CubicSpline_sptr = std::shared_ptr<CubicSpline>;
96using CubicSpline_const_sptr = const std::shared_ptr<CubicSpline>;
97
98} // namespace Functions
99} // namespace CurveFitting
100} // namespace Mantid
double value
The value of the point.
Definition: FitMW.cpp:51
std::map< DeltaEMode::Type, std::string > index
Definition: DeltaEMode.cpp:19
Attribute is a non-fitting parameter.
Definition: IFunction.h:282
A wrapper around GSL functions implementing cubic spline interpolation.
Definition: CubicSpline.h:31
bool m_recalculateSpline
Flag for checking if the spline needs recalculating.
Definition: CubicSpline.h:68
std::string name() const override
overwrite IFunction base class methods
Definition: CubicSpline.h:38
const std::string category() const override
overwrite IFunction base class methods
Definition: CubicSpline.h:39
std::shared_ptr< gsl_spline > m_spline
GSL data structure used to calculate spline.
Definition: CubicSpline.h:65
const int m_min_points
Minimum number of data points in spline.
Definition: CubicSpline.h:53
std::shared_ptr< gsl_interp_accel > m_acc
GSL interpolation accelerator object.
Definition: CubicSpline.h:62
const std::shared_ptr< CubicSpline > CubicSpline_const_sptr
Definition: CubicSpline.h:96
std::shared_ptr< CubicSpline > CubicSpline_sptr
Definition: CubicSpline.h:95
Helper class which provides the Collimation Length for SANS instruments.
Functor to free a GSL objects in a shared pointer.
Definition: CubicSpline.h:56