Mantid
Loading...
Searching...
No Matches
SplineBackground.h
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2008 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//----------------------------------------------------------------------
12#include "MantidAPI/Algorithm.h"
14#include "MantidCurveFitting/DllConfig.h"
15
16#include <gsl/gsl_bspline.h>
17#include <gsl/gsl_multifit.h>
18#include <gsl/gsl_statistics.h>
19
20namespace Mantid {
21namespace CurveFitting {
22namespace Algorithms {
23
29class MANTID_CURVEFITTING_DLL SplineBackground final : public API::Algorithm {
30public:
32 const std::string name() const override { return "SplineBackground"; }
34 int version() const override { return 1; }
35 const std::vector<std::string> seeAlso() const override { return {"Fit", "SplineInterpolation", "SplineSmoothing"}; }
37 const std::string category() const override { return "Optimization;CorrectionFunctions\\BackgroundCorrections"; }
39 const std::string summary() const override { return "Fit spectra background using b-splines."; }
40
41private:
42 // Overridden Algorithm methods
43 void init() override;
44 void exec() override;
45
47 void addWsDataToSpline(const API::MatrixWorkspace *ws, const size_t specNum, int expectedNumBins);
48
50 void allocateBSplinePointers(int numBins, int ncoeffs);
51
53 double calculateBinWeight(double errValue);
54
56 void freeBSplinePointers();
57
59 size_t calculateNumBinsToProcess(const API::MatrixWorkspace *ws);
60
63 API::MatrixWorkspace_sptr saveSplineOutput(const API::MatrixWorkspace_sptr &ws, const size_t spec);
64
66 void setupSpline(double xMin, double xMax, int numBins, int ncoeff);
67
70 gsl_bspline_workspace *splineToProcess{nullptr};
71 gsl_vector *inputSplineWs{nullptr};
72 gsl_vector *xData{nullptr}, *yData{nullptr};
73 gsl_vector *coefficients{nullptr}, *binWeights{nullptr};
74 gsl_matrix *fittedWs{nullptr}, *covariance{nullptr};
75 gsl_multifit_linear_workspace *weightedLinearFitWs{nullptr};
76 };
77
78 bSplinePointers m_splinePointers{};
79};
80
81} // namespace Algorithms
82} // namespace CurveFitting
83} // namespace Mantid
Base class from which all concrete algorithm classes should be derived.
Definition: Algorithm.h:85
Base MatrixWorkspace Abstract Class.
const std::string category() const override
Algorithm's category for identification overriding a virtual method.
const std::string name() const override
Algorithm's name for identification overriding a virtual method.
const std::string summary() const override
Summary of algorithms purpose.
const std::vector< std::string > seeAlso() const override
Function to return all of the seeAlso (these are not validated) algorithms related to this algorithm....
int version() const override
Algorithm's version for identification overriding a virtual method.
std::shared_ptr< MatrixWorkspace > MatrixWorkspace_sptr
shared pointer to the matrix workspace base class
Helper class which provides the Collimation Length for SANS instruments.
Struct holding various pointers required by GSL.