Mantid
Loading...
Searching...
No Matches
FitPeaks.h
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#pragma once
8
16#include "MantidAlgorithms/DllConfig.h"
21
22#include <utility>
23
24namespace Mantid {
25namespace HistogramData {
26class HistogramX;
27class HistogramY;
28} // namespace HistogramData
29
30namespace Algorithms {
31
32namespace FitPeaksAlgorithm {
37
39public:
40 PeakFitResult(size_t num_peaks, size_t num_params);
41 double getPeakPosition(size_t ipeak) const;
42 double getCost(size_t ipeak) const;
43 size_t getNumberParameters() const;
44 size_t getNumberPeaks() const;
45 double getParameterValue(size_t ipeak, size_t iparam) const;
46 double getParameterError(size_t ipeak, size_t iparam) const;
47 void setRecord(size_t ipeak, const double cost, const double peak_position, const FitFunction &fit_functions);
48 void setBadRecord(size_t ipeak, const double peak_position);
49 void setFunctionParameters(size_t ipeak, std::vector<double> &param_values);
50
51private:
54 // goodness of fitting
55 std::vector<double> m_costs;
56 // fitted peak positions
57 std::vector<double> m_fitted_peak_positions;
58 // fitted peak and background parameters
59 std::vector<std::vector<double>> m_function_parameters_vector;
61 std::vector<std::vector<double>> m_function_errors_vector;
62};
63
65public:
70
71public:
72 void setNumberOfSubmittedSpectrumPeaks(const size_t n);
73 void setNumberOfSubmittedIndividualPeaks(const size_t n);
74 void setNumberOfSpectrumPeaksWithLowCount(const size_t n);
75 void setNumberOfOutOfRangePeaks(const size_t n);
78 void setNumberOfPeaksWithLowSignalToNoise(const size_t n);
79 bool isIndividualPeakRejected() const;
80 std::string getReport() const;
81
82private:
83 // number of peaks submitted for spectrum fitting
85 // number of peaks submitted for individual fitting. Since some spectra might fail a pre-check, not all peaks might
86 // make it to the individual fitting
88 // number of peaks rejected as a whole spectrum due to its low signal count
90 // number of peaks rejected individually because their predicted position is out of range
92 // number of peaks rejected individually due to low signal count
94 // number of peask rejected due to not enough data points
96 // number of peaks rejected due to low signal-to-noise ratio
97 size_t m_low_snr;
98};
99} // namespace FitPeaksAlgorithm
100
101class MANTID_ALGORITHMS_DLL FitPeaks final : public API::Algorithm {
102public:
103 FitPeaks();
104
106 const std::string name() const override { return "FitPeaks"; }
107
109 const std::string summary() const override { return "Fit one or multiple peaks in all spectra of a given workspace"; }
110
112 int version() const override { return (1); }
113
115 const std::string category() const override { return "Optimization"; }
116
117 std::map<std::string, std::string> validateInputs() override;
118
119private:
121 void init() override;
123 void exec() override;
124
126 void processInputs();
128 void processInputPeakCenters();
130 void processInputPeakTolerance();
132 void processInputFunctions();
134 void processInputFitRanges();
135
137 void generateFittedParametersValueWorkspaces();
139 void generateOutputPeakPositionWS();
141 void generateCalculatedPeaksWS();
142
145 void convertParametersNameToIndex();
146
148 std::vector<std::shared_ptr<FitPeaksAlgorithm::PeakFitResult>> fitPeaks();
149
151 void fitSpectrumPeaks(size_t wi, const std::vector<double> &expected_peak_centers,
152 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result,
153 std::vector<std::vector<double>> &lastGoodPeakParameters,
154 const std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> &pre_check_result);
155
157 bool fitBackground(const size_t &ws_index, const std::pair<double, double> &fit_window,
158 const double &expected_peak_pos, const API::IBackgroundFunction_sptr &bkgd_func);
159
160 // Peak fitting suite
161 double fitIndividualPeak(size_t wi, const API::IAlgorithm_sptr &fitter, const double expected_peak_center,
162 const std::pair<double, double> &fitwindow, const bool estimate_peak_width,
163 const API::IPeakFunction_sptr &peakfunction, const API::IBackgroundFunction_sptr &bkgdfunc,
164 const std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> &pre_check_result);
165
167 double fitFunctionSD(const API::IAlgorithm_sptr &fit, const API::IPeakFunction_sptr &peak_function,
168 const API::IBackgroundFunction_sptr &bkgd_function, const API::MatrixWorkspace_sptr &dataws,
169 size_t wsindex, const std::pair<double, double> &peak_range, const double &expected_peak_center,
170 bool estimate_peak_width, bool estimate_background);
171
172 double fitFunctionMD(API::IFunction_sptr fit_function, const API::MatrixWorkspace_sptr &dataws, const size_t wsindex,
173 const std::pair<double, double> &vec_xmin, const std::pair<double, double> &vec_xmax);
174
176 double fitFunctionHighBackground(const API::IAlgorithm_sptr &fit, const std::pair<double, double> &fit_window,
177 const size_t &ws_index, const double &expected_peak_center, bool observe_peak_shape,
178 const API::IPeakFunction_sptr &peakfunction,
179 const API::IBackgroundFunction_sptr &bkgdfunc);
180
181 void setupParameterTableWorkspace(const API::ITableWorkspace_sptr &table_ws,
182 const std::vector<std::string> &param_names, bool with_chi2);
183
185 void histRangeToIndexBounds(size_t iws, const std::pair<double, double> &range, size_t &left_index,
186 size_t &right_index);
187
189 size_t histRangeToDataPointCount(size_t iws, const std::pair<double, double> &range);
190
192 void getRangeData(size_t iws, const std::pair<double, double> &range, std::vector<double> &vec_x,
193 std::vector<double> &vec_y, std::vector<double> &vec_e);
194
196 double numberCounts(size_t iws);
197
199 double numberCounts(size_t iws, const std::pair<double, double> &range);
200
202 double calculateSignalToNoiseRatio(size_t iws, const std::pair<double, double> &range,
203 const API::IBackgroundFunction_sptr &bkgd_function);
204
205 API::MatrixWorkspace_sptr createMatrixWorkspace(const std::vector<double> &vec_x, const std::vector<double> &vec_y,
206 const std::vector<double> &vec_e);
207
208 bool decideToEstimatePeakParams(const bool firstPeakInSpectrum, const API::IPeakFunction_sptr &peak_function);
209
211 bool processSinglePeakFitResult(size_t wsindex, size_t peakindex, const double cost,
212 const std::vector<double> &expected_peak_positions,
213 const FitPeaksAlgorithm::FitFunction &fitfunction,
214 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result);
215
217 void calculateFittedPeaks(const std::vector<std::shared_ptr<FitPeaksAlgorithm::PeakFitResult>> &fit_results);
218
219 double calculateSignalToSigmaRatio(const size_t &iws, const std::pair<double, double> &peakWindow,
220 const API::IPeakFunction_sptr &peakFunction);
221
223 std::string getPeakHeightParameterName(const API::IPeakFunction_const_sptr &peak_function);
224
226 void processOutputs(std::vector<std::shared_ptr<FitPeaksAlgorithm::PeakFitResult>> fit_result_vec);
227
229 void writeFitResult(size_t wi, const std::vector<double> &expected_positions,
230 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result);
231
234 bool isObservablePeakProfile(const std::string &peakprofile);
235
236 // log a message disregarding the current logging offset
237 void logNoOffset(const size_t &priority, const std::string &msg);
238
239 //------- Workspaces-------------------------------------
246 API::MatrixWorkspace_sptr m_outputPeakPositionWorkspace; // output workspace for peak positions
259
260 //-------- Functions ------------------------------------------------------
267
269 std::string m_minimizer;
271 std::string m_costFunction;
276
277 //-------- Input param init values --------------------------------
279 std::vector<size_t> m_initParamIndexes;
280
282 std::vector<double> m_peakCenters;
285 std::size_t m_numPeaksToFit;
287
288 std::function<std::vector<double>(std::size_t const &)> m_getExpectedPeakPositions;
289 std::function<std::pair<double, double>(std::size_t const &, std::size_t const &)> m_getPeakFitWindow;
290 void checkWorkspaceIndices(std::size_t const &);
291 void checkPeakIndices(std::size_t const &, std::size_t const &);
292 void checkPeakWindowEdgeOrder(double const &, double const &);
293
296
297 //--------- Fitting range -----------------------------------------
303 std::size_t m_numSpectraToFit;
305 std::vector<double> m_peakPosTolerances;
306
311
313 std::vector<std::vector<double>> m_peakWindowVector;
315
317 std::vector<std::string> m_peakParamNames;
320 std::vector<double> m_initParamValues;
325
326 // Criteria for fitting peaks
330
331 // Criteria for rejecting non-peaks or weak peaks from fitting
334
336
339
340 //----- Result criterias ---------------
343};
344
345} // namespace Algorithms
346} // namespace Mantid
Base class from which all concrete algorithm classes should be derived.
Definition Algorithm.h:76
PeakFitPreCheckResult & operator+=(const PeakFitPreCheckResult &another)
Definition FitPeaks.cpp:199
size_t m_function_parameters_number
number of function parameters
Definition FitPeaks.h:53
std::vector< std::vector< double > > m_function_parameters_vector
Definition FitPeaks.h:59
double getParameterValue(size_t ipeak, size_t iparam) const
get the fitted value of a particular parameter
Definition FitPeaks.cpp:136
void setFunctionParameters(size_t ipeak, std::vector< double > &param_values)
std::vector< std::vector< double > > m_function_errors_vector
fitted peak and background parameters' fitting error
Definition FitPeaks.h:61
void setRecord(size_t ipeak, const double cost, const double peak_position, const FitFunction &fit_functions)
set the peak fitting record/parameter for one peak
Definition FitPeaks.cpp:148
double getParameterError(size_t ipeak, size_t iparam) const
get the fitting error of a particular parameter
Definition FitPeaks.cpp:125
void setBadRecord(size_t ipeak, const double peak_position)
The peak postition should be negative and indicates what went wrong.
Definition FitPeaks.cpp:179
Algorithms::PeakParameterHelper::EstimatePeakWidth m_peakWidthEstimateApproach
Flag for observing peak width: there are 3 states (1) no estimation (2) from 'observation' (3) calcul...
Definition FitPeaks.h:309
API::MatrixWorkspace_const_sptr m_peakCenterWorkspace
Definition FitPeaks.h:283
std::vector< double > m_peakPosTolerances
tolerances for fitting peak positions
Definition FitPeaks.h:305
API::IPeakFunction_sptr m_peakFunction
Peak profile name.
Definition FitPeaks.h:262
API::MatrixWorkspace_sptr m_outputPeakPositionWorkspace
output workspace for peak positions
Definition FitPeaks.h:246
API::MatrixWorkspace_sptr m_fittedPeakWS
matrix workspace contained calcalated peaks+background from fitted result it has same number of spect...
Definition FitPeaks.h:258
API::ITableWorkspace_const_sptr m_profileStartingValueTable
table workspace for profile parameters' starting value
Definition FitPeaks.h:322
std::string m_minimizer
Minimzer.
Definition FitPeaks.h:269
API::IBackgroundFunction_sptr m_linearBackgroundFunction
Linear background function for high background fitting.
Definition FitPeaks.h:266
bool m_peakPosTolCase234
peak positon tolerance case b, c and d
Definition FitPeaks.h:342
double m_peakWidthPercentage
flag to estimate peak width from
Definition FitPeaks.h:295
API::IBackgroundFunction_sptr m_bkgdFunction
Background function.
Definition FitPeaks.h:264
int version() const override
Algorithm's version.
Definition FitPeaks.h:112
const std::string summary() const override
Summary of algorithms purpose.
Definition FitPeaks.h:109
const std::string category() const override
Algorithm's category for identification.
Definition FitPeaks.h:115
double m_minPeakHeight
minimum peak height without background and it also serves as the criteria for observed peak parameter
Definition FitPeaks.h:329
std::string m_costFunction
Cost function.
Definition FitPeaks.h:271
std::size_t m_numSpectraToFit
total number of spectra to be fit
Definition FitPeaks.h:303
std::vector< std::string > m_peakParamNames
input peak parameters' names
Definition FitPeaks.h:317
std::size_t m_numPeaksToFit
the number of peaks to fit in all spectra
Definition FitPeaks.h:285
std::size_t m_startWorkspaceIndex
start index
Definition FitPeaks.h:299
const std::string name() const override
Algorithm's name.
Definition FitPeaks.h:106
API::ITableWorkspace_sptr m_fittedParamTable
output analysis workspaces table workspace for fitted parameters
Definition FitPeaks.h:249
DataObjects::EventWorkspace_const_sptr m_inputEventWS
event workspace for input
Definition FitPeaks.h:244
bool m_highBackground
flag for high background
Definition FitPeaks.h:338
std::vector< double > m_initParamValues
input peak parameters' starting values corresponding to above peak parameter names
Definition FitPeaks.h:320
std::vector< double > m_peakCenters
Designed peak positions and tolerance.
Definition FitPeaks.h:282
std::vector< std::vector< double > > m_peakWindowVector
peak windows
Definition FitPeaks.h:313
std::function< std::pair< double, double >(std::size_t const &, std::size_t const &)> m_getPeakFitWindow
Definition FitPeaks.h:289
API::MatrixWorkspace_sptr m_inputMatrixWS
mandatory input and output workspaces
Definition FitPeaks.h:241
std::vector< size_t > m_initParamIndexes
input starting parameters' indexes in peak function
Definition FitPeaks.h:279
bool m_fitPeaksFromRight
Fit from right or left.
Definition FitPeaks.h:273
std::function< std::vector< double >(std::size_t const &)> m_getExpectedPeakPositions
Definition FitPeaks.h:288
bool m_rawPeaksTable
flag to show that the pamarameters in table are raw parameters or effective parameters
Definition FitPeaks.h:254
int m_fitIterations
Fit iterations.
Definition FitPeaks.h:275
API::MatrixWorkspace_const_sptr m_peakWindowWorkspace
Definition FitPeaks.h:314
bool m_uniformProfileStartingValue
flag for profile startng value being uniform or not
Definition FitPeaks.h:324
API::ITableWorkspace_sptr m_fitErrorTable
table workspace for fitted parameters' fitting error. This is optional
Definition FitPeaks.h:251
std::size_t m_stopWorkspaceIndex
stop index (workspace index of the last spectrum included)
Definition FitPeaks.h:301
std::shared_ptr< IAlgorithm > IAlgorithm_sptr
shared pointer to Mantid::API::IAlgorithm
std::shared_ptr< IBackgroundFunction > IBackgroundFunction_sptr
std::shared_ptr< IPeakFunction > IPeakFunction_sptr
std::shared_ptr< ITableWorkspace > ITableWorkspace_sptr
shared pointer to Mantid::API::ITableWorkspace
std::shared_ptr< const IPeakFunction > IPeakFunction_const_sptr
std::shared_ptr< const ITableWorkspace > ITableWorkspace_const_sptr
shared pointer to Mantid::API::ITableWorkspace (const version)
std::shared_ptr< const MatrixWorkspace > MatrixWorkspace_const_sptr
shared pointer to the matrix workspace base class (const version)
std::shared_ptr< IFunction > IFunction_sptr
shared pointer to the function base class
Definition IFunction.h:743
std::shared_ptr< MatrixWorkspace > MatrixWorkspace_sptr
shared pointer to the matrix workspace base class
std::shared_ptr< const EventWorkspace > EventWorkspace_const_sptr
shared pointer to a const Workspace2D
Helper class which provides the Collimation Length for SANS instruments.
API::IBackgroundFunction_sptr bkgdfunction
Definition FitPeaks.h:35