Mantid
Loading...
Searching...
No Matches
FitPeaks.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//----------------------------------------------------------------------
8// Includes
9//----------------------------------------------------------------------
11#include "MantidAPI/Axis.h"
19#include "MantidAPI/TableRow.h"
27#include "MantidHistogramData/EstimatePolynomial.h"
28#include "MantidHistogramData/Histogram.h"
29#include "MantidHistogramData/HistogramBuilder.h"
30#include "MantidHistogramData/HistogramIterator.h"
38
39#include "boost/algorithm/string.hpp"
40#include "boost/algorithm/string/trim.hpp"
41#include <limits>
42#include <utility>
43
44using namespace Mantid;
45using namespace Algorithms::PeakParameterHelper;
46using namespace Mantid::API;
47using namespace Mantid::DataObjects;
48using namespace Mantid::HistogramData;
49using namespace Mantid::Kernel;
50using namespace Mantid::Geometry;
51using Mantid::HistogramData::Histogram;
52using namespace std;
53
54namespace Mantid::Algorithms {
55
56namespace {
57namespace PropertyNames {
58const std::string INPUT_WKSP("InputWorkspace");
59const std::string OUTPUT_WKSP("OutputWorkspace");
60const std::string START_WKSP_INDEX("StartWorkspaceIndex");
61const std::string STOP_WKSP_INDEX("StopWorkspaceIndex");
62const std::string PEAK_CENTERS("PeakCenters");
63const std::string PEAK_CENTERS_WKSP("PeakCentersWorkspace");
64const std::string PEAK_FUNC("PeakFunction");
65const std::string BACK_FUNC("BackgroundType");
66const std::string FIT_WINDOW_LIST("FitWindowBoundaryList");
67const std::string FIT_WINDOW_WKSP("FitPeakWindowWorkspace");
68const std::string PEAK_WIDTH_PERCENT("PeakWidthPercent");
69const std::string PEAK_PARAM_NAMES("PeakParameterNames");
70const std::string PEAK_PARAM_VALUES("PeakParameterValues");
71const std::string PEAK_PARAM_TABLE("PeakParameterValueTable");
72const std::string FIT_FROM_RIGHT("FitFromRight");
73const std::string MINIMIZER("Minimizer");
74const std::string COST_FUNC("CostFunction");
75const std::string MAX_FIT_ITER("MaxFitIterations");
76const std::string BACKGROUND_Z_SCORE("FindBackgroundSigma");
77const std::string HIGH_BACKGROUND("HighBackground");
78const std::string POSITION_TOL("PositionTolerance");
79const std::string PEAK_MIN_HEIGHT("MinimumPeakHeight");
80const std::string CONSTRAIN_PEAK_POS("ConstrainPeakPositions");
81const std::string COPY_LAST_GOOD_PEAK_PARAMS("CopyLastGoodPeakParameters");
82const std::string OUTPUT_WKSP_MODEL("FittedPeaksWorkspace");
83const std::string OUTPUT_WKSP_PARAMS("OutputPeakParametersWorkspace");
84const std::string OUTPUT_WKSP_PARAM_ERRS("OutputParameterFitErrorsWorkspace");
85const std::string RAW_PARAMS("RawPeakParameters");
86const std::string PEAK_MIN_SIGNAL_TO_NOISE_RATIO("MinimumSignalToNoiseRatio");
87const std::string PEAK_MIN_TOTAL_COUNT("MinimumPeakTotalCount");
88const std::string PEAK_MIN_SIGNAL_TO_SIGMA_RATIO("MinimumSignalToSigmaRatio");
89} // namespace PropertyNames
90} // namespace
91
92namespace FitPeaksAlgorithm {
93
94//----------------------------------------------------------------------------------------------
96PeakFitResult::PeakFitResult(size_t num_peaks, size_t num_params) : m_function_parameters_number(num_params) {
97 // check input
98 if (num_peaks == 0 || num_params == 0)
99 throw std::runtime_error("No peak or no parameter error.");
100
101 //
102 m_fitted_peak_positions.resize(num_peaks, std::numeric_limits<double>::quiet_NaN());
103 m_costs.resize(num_peaks, DBL_MAX);
104 m_function_parameters_vector.resize(num_peaks);
105 m_function_errors_vector.resize(num_peaks);
106 for (size_t ipeak = 0; ipeak < num_peaks; ++ipeak) {
107 m_function_parameters_vector[ipeak].resize(num_params, std::numeric_limits<double>::quiet_NaN());
108 m_function_errors_vector[ipeak].resize(num_params, std::numeric_limits<double>::quiet_NaN());
109 }
110
111 return;
112}
113
114//----------------------------------------------------------------------------------------------
116
118
119//----------------------------------------------------------------------------------------------
126double PeakFitResult::getParameterError(size_t ipeak, size_t iparam) const {
127 return m_function_errors_vector[ipeak][iparam];
128}
129
130//----------------------------------------------------------------------------------------------
137double PeakFitResult::getParameterValue(size_t ipeak, size_t iparam) const {
138 return m_function_parameters_vector[ipeak][iparam];
139}
140
141//----------------------------------------------------------------------------------------------
142double PeakFitResult::getPeakPosition(size_t ipeak) const { return m_fitted_peak_positions[ipeak]; }
143
144//----------------------------------------------------------------------------------------------
145double PeakFitResult::getCost(size_t ipeak) const { return m_costs[ipeak]; }
146
147//----------------------------------------------------------------------------------------------
149void PeakFitResult::setRecord(size_t ipeak, const double cost, const double peak_position,
150 const FitFunction &fit_functions) {
151 // check input
152 if (ipeak >= m_costs.size())
153 throw std::runtime_error("Peak index is out of range.");
154
155 // set the values
156 m_costs[ipeak] = cost;
157
158 // set peak position
159 m_fitted_peak_positions[ipeak] = peak_position;
160
161 // transfer from peak function to vector
162 size_t peak_num_params = fit_functions.peakfunction->nParams();
163 for (size_t ipar = 0; ipar < peak_num_params; ++ipar) {
164 // peak function
165 m_function_parameters_vector[ipeak][ipar] = fit_functions.peakfunction->getParameter(ipar);
166 m_function_errors_vector[ipeak][ipar] = fit_functions.peakfunction->getError(ipar);
167 }
168 for (size_t ipar = 0; ipar < fit_functions.bkgdfunction->nParams(); ++ipar) {
169 // background function
170 m_function_parameters_vector[ipeak][ipar + peak_num_params] = fit_functions.bkgdfunction->getParameter(ipar);
171 m_function_errors_vector[ipeak][ipar + peak_num_params] = fit_functions.bkgdfunction->getError(ipar);
172 }
173}
174
175//----------------------------------------------------------------------------------------------
180void PeakFitResult::setBadRecord(size_t ipeak, const double peak_position) {
181 // check input
182 if (ipeak >= m_costs.size())
183 throw std::runtime_error("Peak index is out of range");
184 if (peak_position >= 0.)
185 throw std::runtime_error("Can only set negative postion for bad record");
186
187 // set the values
188 m_costs[ipeak] = DBL_MAX;
189
190 // set peak position
191 m_fitted_peak_positions[ipeak] = peak_position;
192
193 // transfer from peak function to vector
194 for (size_t ipar = 0; ipar < m_function_parameters_number; ++ipar) {
195 m_function_parameters_vector[ipeak][ipar] = 0.;
196 m_function_errors_vector[ipeak][ipar] = std::numeric_limits<double>::quiet_NaN();
197 }
198}
199
211
213
215
217
219
221
223
225
227 // the method should be used on an individual peak, not on a spectrum
228 assert(m_submitted_spectrum_peaks == 0);
229 assert(m_submitted_individual_peaks == 1);
230
231 // if a peak is rejected, it is rejected based on the very first check it fails
232 size_t individual_rejection_count = m_low_count_individual + m_not_enough_datapoints + m_low_snr;
233 assert(individual_rejection_count <= 1);
234
235 return individual_rejection_count == 1;
236}
237
240
241 // if no peaks were rejected by the pre-check, keep quiet
243 return "";
244
245 std::ostringstream os;
246 os << "Total number of peaks pre-checked before fitting: " << m_submitted_spectrum_peaks << "\n";
247 if (m_low_count_spectrum > 0)
248 os << m_low_count_spectrum << " peak(s) rejected: low signal count (whole spectrum).\n";
249 if (m_out_of_range > 0)
250 os << m_out_of_range << " peak(s) rejected: out of range.\n";
252 os << m_not_enough_datapoints << " peak(s) rejected: not enough X(Y) datapoints.\n";
254 os << m_low_count_individual << " peak(s) rejected: low signal count (individual peak).\n";
255 if (m_low_snr > 0)
256 os << m_low_snr << " peak(s) rejected: low signal-to-noise ratio.\n";
257
258 return os.str();
259}
260} // namespace FitPeaksAlgorithm
261
262//----------------------------------------------------------------------------------------------
264 : m_fitPeaksFromRight(true), m_fitIterations(50), m_numPeaksToFit(0), m_minPeakHeight(0.),
265 m_minSignalToNoiseRatio(0.), m_minPeakTotalCount(0.), m_peakPosTolCase234(false) {}
266
267//----------------------------------------------------------------------------------------------
272 "Name of the input workspace for peak fitting.");
275 "Name of the output workspace containing peak centers for "
276 "fitting offset."
277 "The output workspace is point data."
278 "Each workspace index corresponds to a spectrum. "
279 "Each X value ranges from 0 to N-1, where N is the number of "
280 "peaks to fit. "
281 "Each Y value is the peak position obtained by peak fitting. "
282 "Negative value is used for error signals. "
283 "-1 for data is zero; -2 for maximum value is smaller than "
284 "specified minimum value."
285 "and -3 for non-converged fitting.");
286
287 // properties about fitting range and criteria
288 auto mustBePositive = std::make_shared<BoundedValidator<int>>();
289 mustBePositive->setLower(0);
290 declareProperty(PropertyNames::START_WKSP_INDEX, 0, mustBePositive, "Starting workspace index for fit");
292 PropertyNames::STOP_WKSP_INDEX, EMPTY_INT(),
293 "Last workspace index for fit is the smaller of this value and the workspace index of last spectrum.");
294 // properties about peak positions to fit
295 declareProperty(std::make_unique<ArrayProperty<double>>(PropertyNames::PEAK_CENTERS),
296 "List of peak centers to use as initial guess for fit.");
298 std::make_unique<WorkspaceProperty<MatrixWorkspace>>(PropertyNames::PEAK_CENTERS_WKSP, "", Direction::Input,
300 "MatrixWorkspace containing referent peak centers for each spectrum, defined at the same workspace indices.");
301
302 const std::string peakcentergrp("Peak Positions");
303 setPropertyGroup(PropertyNames::PEAK_CENTERS, peakcentergrp);
304 setPropertyGroup(PropertyNames::PEAK_CENTERS_WKSP, peakcentergrp);
305
306 // properties about peak profile
307 const std::vector<std::string> peakNames = FunctionFactory::Instance().getFunctionNames<API::IPeakFunction>();
308 declareProperty(PropertyNames::PEAK_FUNC, "Gaussian", std::make_shared<StringListValidator>(peakNames),
309 "Use of a BackToBackExponential profile is only reccomended if the "
310 "coeficients to calculate A and B are defined in the instrument "
311 "Parameters.xml file.");
312 const vector<string> bkgdtypes{"Flat", "Linear", "Quadratic"};
313 declareProperty(PropertyNames::BACK_FUNC, "Linear", std::make_shared<StringListValidator>(bkgdtypes),
314 "Type of Background.");
315
316 const std::string funcgroup("Function Types");
317 setPropertyGroup(PropertyNames::PEAK_FUNC, funcgroup);
318 setPropertyGroup(PropertyNames::BACK_FUNC, funcgroup);
319
320 // properties about peak range including fitting window and peak width
321 // (percentage)
322 declareProperty(std::make_unique<ArrayProperty<double>>(PropertyNames::FIT_WINDOW_LIST),
323 "List of boundaries of the peak fitting window corresponding to "
324 "PeakCenters.");
325
326 declareProperty(std::make_unique<WorkspaceProperty<MatrixWorkspace>>(PropertyNames::FIT_WINDOW_WKSP, "",
328 "MatrixWorkspace containing peak windows for each peak center in each spectrum, defined at the same "
329 "workspace indices.");
330
331 auto min = std::make_shared<BoundedValidator<double>>();
332 min->setLower(1e-3);
333 // min->setUpper(1.); TODO make this a limit
334 declareProperty(PropertyNames::PEAK_WIDTH_PERCENT, EMPTY_DBL(), min,
335 "The estimated peak width as a "
336 "percentage of the d-spacing "
337 "of the center of the peak. Value must be less than 1.");
338
339 const std::string fitrangeegrp("Peak Range Setup");
340 setPropertyGroup(PropertyNames::PEAK_WIDTH_PERCENT, fitrangeegrp);
341 setPropertyGroup(PropertyNames::FIT_WINDOW_LIST, fitrangeegrp);
342 setPropertyGroup(PropertyNames::FIT_WINDOW_WKSP, fitrangeegrp);
343
344 // properties about peak parameters' names and value
345 declareProperty(std::make_unique<ArrayProperty<std::string>>(PropertyNames::PEAK_PARAM_NAMES),
346 "List of peak parameters' names");
347 declareProperty(std::make_unique<ArrayProperty<double>>(PropertyNames::PEAK_PARAM_VALUES),
348 "List of peak parameters' value");
349 declareProperty(std::make_unique<WorkspaceProperty<TableWorkspace>>(PropertyNames::PEAK_PARAM_TABLE, "",
351 "Name of the an optional workspace, whose each column "
352 "corresponds to given peak parameter names, "
353 "and each row corresponds to a subset of spectra.");
354
355 const std::string startvaluegrp("Starting Parameters Setup");
356 setPropertyGroup(PropertyNames::PEAK_PARAM_NAMES, startvaluegrp);
357 setPropertyGroup(PropertyNames::PEAK_PARAM_VALUES, startvaluegrp);
358 setPropertyGroup(PropertyNames::PEAK_PARAM_TABLE, startvaluegrp);
359
360 // optimization setup
361 declareProperty(PropertyNames::FIT_FROM_RIGHT, true,
362 "Flag for the order to fit peaks. If true, peaks are fitted "
363 "from rightmost;"
364 "Otherwise peaks are fitted from leftmost.");
365
366 const std::vector<std::string> minimizerOptions = API::FuncMinimizerFactory::Instance().getKeys();
367 declareProperty(PropertyNames::MINIMIZER, "Levenberg-Marquardt",
369 "Minimizer to use for fitting.");
370
371 const std::array<string, 2> costFuncOptions = {{"Least squares", "Rwp"}};
372 declareProperty(PropertyNames::COST_FUNC, "Least squares",
373 Kernel::IValidator_sptr(new Kernel::ListValidator<std::string>(costFuncOptions)), "Cost functions");
374
375 auto min_max_iter = std::make_shared<BoundedValidator<int>>();
376 min_max_iter->setLower(49);
377 declareProperty(PropertyNames::MAX_FIT_ITER, 50, min_max_iter, "Maximum number of function fitting iterations.");
378
379 const std::string optimizergrp("Optimization Setup");
380 setPropertyGroup(PropertyNames::MINIMIZER, optimizergrp);
381 setPropertyGroup(PropertyNames::COST_FUNC, optimizergrp);
382
383 // other helping information
384 std::ostringstream os;
385 os << "Deprecated property. Use " << PropertyNames::PEAK_MIN_SIGNAL_TO_NOISE_RATIO << " instead.";
386 declareProperty(PropertyNames::BACKGROUND_Z_SCORE, EMPTY_DBL(), os.str());
387
388 declareProperty(PropertyNames::HIGH_BACKGROUND, true,
389 "Flag whether the input data has high background compared to peak heights.");
390
391 declareProperty(std::make_unique<ArrayProperty<double>>(PropertyNames::POSITION_TOL),
392 "List of tolerance on fitted peak positions against given peak positions."
393 "If there is only one value given, then ");
394
395 declareProperty(PropertyNames::PEAK_MIN_HEIGHT, 0.,
396 "Used for validating peaks before and after fitting. If a peak's observed/estimated or "
397 "fitted height is under this value, the peak will be marked as error.");
398
399 declareProperty(PropertyNames::CONSTRAIN_PEAK_POS, true,
400 "If true peak position will be constrained by estimated positions "
401 "(highest Y value position) and "
402 "the peak width either estimted by observation or calculate.");
403
404 declareProperty(PropertyNames::COPY_LAST_GOOD_PEAK_PARAMS, true,
405 "If true, initial peak parameters (with the exception of peak centre) "
406 "may be copied from the last successfully fit peak in the spectra.");
407
408 // additional output for reviewing
409 declareProperty(std::make_unique<WorkspaceProperty<MatrixWorkspace>>(PropertyNames::OUTPUT_WKSP_MODEL, "",
411 "Name of the output matrix workspace with fitted peak. "
412 "This output workspace has the same dimension as the input workspace."
413 "The Y values belonged to peaks to fit are replaced by fitted value. "
414 "Values of estimated background are used if peak fails to be fit.");
415
416 declareProperty(std::make_unique<WorkspaceProperty<API::ITableWorkspace>>(PropertyNames::OUTPUT_WKSP_PARAMS, "",
418 "Name of table workspace containing all fitted peak parameters.");
419
420 // Optional output table workspace for each individual parameter's fitting
421 // error
422 declareProperty(std::make_unique<WorkspaceProperty<API::ITableWorkspace>>(PropertyNames::OUTPUT_WKSP_PARAM_ERRS, "",
424 "Name of workspace containing all fitted peak parameters' fitting error."
425 "It must be used along with FittedPeaksWorkspace and RawPeakParameters "
426 "(True)");
427
428 declareProperty(PropertyNames::RAW_PARAMS, true,
429 "false generates table with effective centre/width/height "
430 "parameters. true generates a table with peak function "
431 "parameters");
432
434 PropertyNames::PEAK_MIN_SIGNAL_TO_NOISE_RATIO, 0.,
435 "Used for validating peaks before fitting. If the signal-to-noise ratio is under this value, "
436 "the peak will be marked as error. This does not apply to peaks for which the noise cannot be estimated.");
437
438 declareProperty(PropertyNames::PEAK_MIN_TOTAL_COUNT, EMPTY_DBL(),
439 "Used for validating peaks before fitting. If the total peak window Y-value count "
440 "is under this value, the peak will be excluded from fitting and calibration.");
441
442 declareProperty(PropertyNames::PEAK_MIN_SIGNAL_TO_SIGMA_RATIO, 0.,
443 "Used for validating peaks after fitting. If the signal-to-sigma ratio is under this value, "
444 "the peak will be excluded from fitting and calibration.");
445
446 const std::string addoutgrp("Analysis");
447 setPropertyGroup(PropertyNames::OUTPUT_WKSP_PARAMS, addoutgrp);
448 setPropertyGroup(PropertyNames::OUTPUT_WKSP_MODEL, addoutgrp);
449 setPropertyGroup(PropertyNames::OUTPUT_WKSP_PARAM_ERRS, addoutgrp);
450 setPropertyGroup(PropertyNames::RAW_PARAMS, addoutgrp);
451}
452
453//----------------------------------------------------------------------------------------------
456std::map<std::string, std::string> FitPeaks::validateInputs() {
457 map<std::string, std::string> issues;
458
459 // check that min/max spectra indices make sense - only matters if both are specified
460 if (!(isDefault(PropertyNames::START_WKSP_INDEX) && isDefault(PropertyNames::STOP_WKSP_INDEX))) {
461 const int startIndex = getProperty(PropertyNames::START_WKSP_INDEX);
462 const int stopIndex = getProperty(PropertyNames::STOP_WKSP_INDEX);
463 if (startIndex > stopIndex) {
464 const std::string msg =
465 PropertyNames::START_WKSP_INDEX + " must be less than or equal to " + PropertyNames::STOP_WKSP_INDEX;
466 issues[PropertyNames::START_WKSP_INDEX] = msg;
467 issues[PropertyNames::STOP_WKSP_INDEX] = msg;
468 }
469 }
470
471 // check that the peak parameters are in parallel properties
472 bool haveCommonPeakParameters(false);
473 std::vector<string> suppliedParameterNames = getProperty(PropertyNames::PEAK_PARAM_NAMES);
474 std::vector<double> peakParamValues = getProperty(PropertyNames::PEAK_PARAM_VALUES);
475 if ((!suppliedParameterNames.empty()) || (!peakParamValues.empty())) {
476 haveCommonPeakParameters = true;
477 if (suppliedParameterNames.size() != peakParamValues.size()) {
478 issues[PropertyNames::PEAK_PARAM_NAMES] = "must have same number of values as PeakParameterValues";
479 issues[PropertyNames::PEAK_PARAM_VALUES] = "must have same number of values as PeakParameterNames";
480 }
481 }
482
483 // get the information out of the table
484 std::string partablename = getPropertyValue(PropertyNames::PEAK_PARAM_TABLE);
485 if (!partablename.empty()) {
486 if (haveCommonPeakParameters) {
487 const std::string msg = "Parameter value table and initial parameter "
488 "name/value vectors cannot be given "
489 "simultanenously.";
490 issues[PropertyNames::PEAK_PARAM_TABLE] = msg;
491 issues[PropertyNames::PEAK_PARAM_NAMES] = msg;
492 issues[PropertyNames::PEAK_PARAM_VALUES] = msg;
493 } else {
494 m_profileStartingValueTable = getProperty(PropertyNames::PEAK_PARAM_TABLE);
495 suppliedParameterNames = m_profileStartingValueTable->getColumnNames();
496 }
497 }
498
499 // check that the suggested peak parameter names exist in the peak function
500 if (!suppliedParameterNames.empty()) {
501 std::string peakfunctiontype = getPropertyValue(PropertyNames::PEAK_FUNC);
503 std::dynamic_pointer_cast<IPeakFunction>(API::FunctionFactory::Instance().createFunction(peakfunctiontype));
504
505 // put the names in a vector
506 std::vector<string> functionParameterNames;
507 for (size_t i = 0; i < m_peakFunction->nParams(); ++i)
508 functionParameterNames.emplace_back(m_peakFunction->parameterName(i));
509 // check that the supplied names are in the function
510 // it is acceptable to be missing parameters
511 const bool failed = std::any_of(suppliedParameterNames.cbegin(), suppliedParameterNames.cend(),
512 [&functionParameterNames](const auto &parName) {
513 return std::find(functionParameterNames.begin(), functionParameterNames.end(),
514 parName) == functionParameterNames.end();
515 });
516 if (failed) {
517 std::string msg = "Specified invalid parameter for peak function";
518 if (haveCommonPeakParameters)
519 issues[PropertyNames::PEAK_PARAM_NAMES] = msg;
520 else
521 issues[PropertyNames::PEAK_PARAM_TABLE] = msg;
522 }
523 }
524
525 // check inputs for uncertainty (fitting error)
526 const std::string error_table_name = getPropertyValue(PropertyNames::OUTPUT_WKSP_PARAM_ERRS);
527 if (!error_table_name.empty()) {
528 const bool use_raw_params = getProperty(PropertyNames::RAW_PARAMS);
529 if (!use_raw_params) {
530 issues[PropertyNames::OUTPUT_WKSP_PARAM_ERRS] = "Cannot be used with " + PropertyNames::RAW_PARAMS + "=False";
531 issues[PropertyNames::RAW_PARAMS] =
532 "Cannot be False with " + PropertyNames::OUTPUT_WKSP_PARAM_ERRS + " specified";
533 }
534 }
535
536 return issues;
537}
538
539//----------------------------------------------------------------------------------------------
541 // process inputs
543
544 // create output workspace: fitted peak positions
546
547 // create output workspace: fitted peaks' parameters values
549
550 // create output workspace: calculated from fitted peak and background
552
553 // fit peaks
554 auto fit_results = fitPeaks();
555
556 // set the output workspaces to properites
557 processOutputs(fit_results);
558}
559
560//----------------------------------------------------------------------------------------------
562 // input workspaces
564
565 if (m_inputMatrixWS->getAxis(0)->unit()->unitID() == "dSpacing")
566 m_inputIsDSpace = true;
567 else
568 m_inputIsDSpace = false;
569
570 // spectra to fit
571 int start_wi = getProperty(PropertyNames::START_WKSP_INDEX);
572 m_startWorkspaceIndex = static_cast<size_t>(start_wi);
573
574 // last spectrum's workspace index, which is included
575 int stop_wi = getProperty(PropertyNames::STOP_WKSP_INDEX);
576 if (isEmpty(stop_wi))
577 m_stopWorkspaceIndex = m_inputMatrixWS->getNumberHistograms() - 1;
578 else {
579 m_stopWorkspaceIndex = static_cast<size_t>(stop_wi);
580 if (m_stopWorkspaceIndex > m_inputMatrixWS->getNumberHistograms() - 1)
581 m_stopWorkspaceIndex = m_inputMatrixWS->getNumberHistograms() - 1;
582 }
583
584 // total number of spectra to be fit
586
587 // optimizer, cost function and fitting scheme
588 m_minimizer = getPropertyValue(PropertyNames::MINIMIZER);
589 m_costFunction = getPropertyValue(PropertyNames::COST_FUNC);
590 m_fitPeaksFromRight = getProperty(PropertyNames::FIT_FROM_RIGHT);
591 m_constrainPeaksPosition = getProperty(PropertyNames::CONSTRAIN_PEAK_POS);
592 m_fitIterations = getProperty(PropertyNames::MAX_FIT_ITER);
593 m_copyLastGoodPeakParameters = getProperty(PropertyNames::COPY_LAST_GOOD_PEAK_PARAMS);
594
595 // Peak centers, tolerance and fitting range
597 // check
598 if (m_numPeaksToFit == 0)
599 throw std::runtime_error("number of peaks to fit is zero.");
600 // about how to estimate the peak width
601 m_peakWidthPercentage = getProperty(PropertyNames::PEAK_WIDTH_PERCENT);
604 if (m_peakWidthPercentage >= 1.) // TODO
605 throw std::runtime_error("PeakWidthPercent must be less than 1");
606 g_log.debug() << "peak width/value = " << m_peakWidthPercentage << "\n";
607
608 // set up background
609 m_highBackground = getProperty(PropertyNames::HIGH_BACKGROUND);
610 double temp = getProperty(PropertyNames::BACKGROUND_Z_SCORE);
611 if (!isEmpty(temp)) {
612 std::ostringstream os;
613 os << "FitPeaks property \"" << PropertyNames::BACKGROUND_Z_SCORE << "\" is deprecated and will be ignored."
614 << "\n";
615 logNoOffset(4 /*warning*/, os.str());
616 }
617
618 // Set up peak and background functions
620
621 // about peak width and other peak parameter estimating method
622 if (m_peakWidthPercentage > 0.)
623 m_peakWidthEstimateApproach = EstimatePeakWidth::InstrumentResolution;
624 else if (isObservablePeakProfile((m_peakFunction->name())))
625 m_peakWidthEstimateApproach = EstimatePeakWidth::Observation;
626 else
627 m_peakWidthEstimateApproach = EstimatePeakWidth::NoEstimation;
628 // m_peakWidthEstimateApproach = EstimatePeakWidth::NoEstimation;
629 g_log.debug() << "Process inputs [3] peak type: " << m_peakFunction->name()
630 << ", background type: " << m_bkgdFunction->name() << "\n";
631
634
635 return;
636}
637
638//----------------------------------------------------------------------------------------------
642 // peak functions
643 std::string peakfunctiontype = getPropertyValue(PropertyNames::PEAK_FUNC);
645 std::dynamic_pointer_cast<IPeakFunction>(API::FunctionFactory::Instance().createFunction(peakfunctiontype));
646
647 // background functions
648 std::string bkgdfunctiontype = getPropertyValue(PropertyNames::BACK_FUNC);
649 std::string bkgdname;
650 if (bkgdfunctiontype == "Linear")
651 bkgdname = "LinearBackground";
652 else if (bkgdfunctiontype == "Flat") {
653 g_log.warning("There may be problems with Flat background");
654 bkgdname = "FlatBackground";
655 } else
656 bkgdname = bkgdfunctiontype;
658 std::dynamic_pointer_cast<IBackgroundFunction>(API::FunctionFactory::Instance().createFunction(bkgdname));
660 m_linearBackgroundFunction = std::dynamic_pointer_cast<IBackgroundFunction>(
661 API::FunctionFactory::Instance().createFunction("LinearBackground"));
662 else
664
665 // TODO check that both parameter names and values exist
666 // input peak parameters
667 std::string partablename = getPropertyValue(PropertyNames::PEAK_PARAM_TABLE);
668 m_peakParamNames = getProperty(PropertyNames::PEAK_PARAM_NAMES);
669
671 if (partablename.empty() && (!m_peakParamNames.empty())) {
672 // use uniform starting value of peak parameters
673 m_initParamValues = getProperty(PropertyNames::PEAK_PARAM_VALUES);
674 // convert the parameter name in string to parameter name in integer index
676 // m_uniformProfileStartingValue = true;
677 } else if ((!partablename.empty()) && m_peakParamNames.empty()) {
678 // use non-uniform starting value of peak parameters
680 } else if (peakfunctiontype != "Gaussian") {
681 // user specifies nothing
682 g_log.warning("Neither parameter value table nor initial "
683 "parameter name/value vectors is specified. Fitting might "
684 "not be reliable for peak profile other than Gaussian");
685 }
686
687 return;
688}
689
690//----------------------------------------------------------------------------------------------
695 // get peak fit window
696 std::vector<double> peakwindow = getProperty(PropertyNames::FIT_WINDOW_LIST);
697 std::string peakwindowname = getPropertyValue(PropertyNames::FIT_WINDOW_WKSP);
698 API::MatrixWorkspace_const_sptr peakwindowws = getProperty(PropertyNames::FIT_WINDOW_WKSP);
699
700 // in most case, calculate window by instrument resolution is False
701
702 if ((!peakwindow.empty()) && peakwindowname.empty()) {
703 // Peak windows are uniform among spectra: use vector for peak windows
704
705 // check peak positions
707 throw std::invalid_argument(
708 "Specifying peak windows with a list requires also specifying peak positions with a list.");
709 // check size
710 if (peakwindow.size() != m_numPeaksToFit * 2)
711 throw std::invalid_argument("Peak window vector must be twice as large as number of peaks.");
712
713 // set up window to m_peakWindowVector
715 for (size_t i = 0; i < m_numPeaksToFit; ++i) {
716 std::vector<double> peakranges(2);
717 peakranges[0] = peakwindow[i * 2];
718 peakranges[1] = peakwindow[i * 2 + 1];
719 // check peak window (range) against peak centers
720 if ((peakranges[0] < m_peakCenters[i]) && (m_peakCenters[i] < peakranges[1])) {
721 // pass check: set
722 m_peakWindowVector[i] = peakranges;
723 } else {
724 // failed
725 std::stringstream errss;
726 errss << "Peak " << i << ": user specifies an invalid range and peak center against " << peakranges[0] << " < "
727 << m_peakCenters[i] << " < " << peakranges[1];
728 throw std::invalid_argument(errss.str());
729 }
730 } // END-FOR
731 m_getPeakFitWindow = [this](std::size_t wi, std::size_t ipeak) -> std::pair<double, double> {
732 this->checkWorkspaceIndices(wi);
733 this->checkPeakIndices(wi, ipeak);
734 double left = this->m_peakWindowVector[ipeak][0];
735 double right = this->m_peakWindowVector[ipeak][1];
736 this->checkPeakWindowEdgeOrder(left, right);
737 return std::make_pair(left, right);
738 };
739 // END if list peak windows
740 } else if (peakwindow.empty() && peakwindowws != nullptr) {
741 // use matrix workspace for non-uniform peak windows
742 m_peakWindowWorkspace = getProperty(PropertyNames::FIT_WINDOW_WKSP);
743
744 // check each spectrum whether the window is defined with the correct size
745 for (std::size_t wi = m_startWorkspaceIndex; wi <= m_stopWorkspaceIndex; wi++) {
746 const auto &peakWindowX = m_peakWindowWorkspace->x(wi);
747 const auto &peakCenterX = m_peakCenterWorkspace->x(wi);
748 if (peakWindowX.empty()) {
749 std::stringstream errss;
750 errss << "Peak window required at workspace index " << wi << " "
751 << "which is undefined in the peak window workspace. "
752 << "Ensure workspace indices correspond in peak window workspace and input workspace "
753 << "when using start and stop indices.";
754 throw std::invalid_argument(errss.str());
755 }
756 // check size
757 if (peakWindowX.size() % 2 != 0) {
758 throw std::invalid_argument("The peak window vector must be even, with two edges for each peak center.");
759 }
760 if (peakWindowX.size() != peakCenterX.size() * 2) {
761 std::stringstream errss;
762 errss << "Peak window workspace index " << wi << " has incompatible number of fit windows "
763 << peakWindowX.size() / 2 << " with the number of peaks " << peakCenterX.size() << " to fit.";
764 throw std::invalid_argument(errss.str());
765 }
766
767 for (size_t ipeak = 0; ipeak < peakCenterX.size(); ++ipeak) {
768 double left_w_bound = peakWindowX[ipeak * 2];
769 double right_w_bound = peakWindowX[ipeak * 2 + 1];
770 double center = peakCenterX[ipeak];
771
772 if (!(left_w_bound < center && center < right_w_bound)) {
773 std::stringstream errss;
774 errss << "Workspace index " << wi << " has incompatible peak window "
775 << "(" << left_w_bound << ", " << right_w_bound << ") "
776 << "with " << ipeak << "-th expected peak's center " << center;
777 throw std::runtime_error(errss.str());
778 }
779 }
780 }
781 m_getPeakFitWindow = [this](std::size_t wi, std::size_t ipeak) -> std::pair<double, double> {
782 this->checkWorkspaceIndices(wi);
783 this->checkPeakIndices(wi, ipeak);
784 double left = m_peakWindowWorkspace->x(wi)[ipeak * 2];
785 double right = m_peakWindowWorkspace->x(wi)[ipeak * 2 + 1];
786 this->checkPeakWindowEdgeOrder(left, right);
787 return std::make_pair(left, right);
788 };
789 // END if workspace peak windows
790 } else if (peakwindow.empty()) {
791 // no peak window is defined, then the peak window will be estimated by
792 // delta(D)/D
794 // m_peakWindowMethod = PeakWindowMethod::TOLERANCE;
795 // m_calculateWindowInstrument = true;
796 m_getPeakFitWindow = [this](std::size_t wi, std::size_t ipeak) -> std::pair<double, double> {
797 this->checkWorkspaceIndices(wi);
798 this->checkPeakIndices(wi, ipeak);
799 // calcualte peak window by delta(d)/d
800 double peak_pos = this->m_getExpectedPeakPositions(wi)[ipeak];
801 // calcalate expected peak width
802 double estimate_peak_width = peak_pos * m_peakWidthPercentage;
803 // using the NUMBER THREE to estimate the peak window
804 double THREE = 3.0;
805 double left = peak_pos - estimate_peak_width * THREE;
806 double right = peak_pos + estimate_peak_width * THREE;
807 this->checkPeakWindowEdgeOrder(left, right);
808 return std::make_pair(left, right);
809 };
810 } else {
811 throw std::invalid_argument("Without definition of peak window, the "
812 "input workspace must be in unit of dSpacing "
813 "and Delta(D)/D must be given!");
814 }
815 } else {
816 // non-supported situation
817 throw std::invalid_argument("One and only one of peak window array and "
818 "peak window workspace can be specified.");
819 }
820
821 return;
822}
823
824//----------------------------------------------------------------------------------------------
834 // peak centers
835 m_peakCenters = getProperty(PropertyNames::PEAK_CENTERS);
836 API::MatrixWorkspace_const_sptr peakcenterws = getProperty(PropertyNames::PEAK_CENTERS_WKSP);
837 if (!peakcenterws)
838 g_log.notice("Peak centers are not specified by peak center workspace");
839
840 std::string peakpswsname = getPropertyValue(PropertyNames::PEAK_CENTERS_WKSP);
841 if ((!m_peakCenters.empty()) && peakcenterws == nullptr) {
842 // peak positions are uniform among all spectra
844 // number of peaks to fit!
846 m_getExpectedPeakPositions = [this](std::size_t wi) -> std::vector<double> {
847 this->checkWorkspaceIndices(wi);
848 return this->m_peakCenters;
849 };
850 } else if (m_peakCenters.empty() && peakcenterws != nullptr) {
851 // peak positions can be different among spectra
853 m_peakCenterWorkspace = getProperty(PropertyNames::PEAK_CENTERS_WKSP);
854 // number of peaks to fit must correspond to largest number of reference peaks
855 m_numPeaksToFit = 0;
856 g_log.debug() << "Input peak center workspace: " << m_peakCenterWorkspace->x(0).size() << ", "
857 << m_peakCenterWorkspace->y(0).size() << "\n";
858 for (std::size_t wi = m_startWorkspaceIndex; wi <= m_stopWorkspaceIndex; wi++) {
859 if (m_peakCenterWorkspace->x(wi).empty()) {
860 std::stringstream errss;
861 errss << "Fit peaks was asked to fit from workspace index " << m_startWorkspaceIndex << " "
862 << "until workspace index " << m_stopWorkspaceIndex << ". "
863 << "However, the peak center workspace does not have values defined "
864 << "at workspace index " << wi << ". "
865 << "Make sure the workspace indices between input and peak center workspaces correspond.";
866 g_log.error() << errss.str();
867 throw std::invalid_argument(errss.str());
868 }
869 // the number of peaks to try to fit should be the max number of peaks across spectra
870 m_numPeaksToFit = std::max(m_numPeaksToFit, m_peakCenterWorkspace->x(wi).size());
871 }
872 m_getExpectedPeakPositions = [this](std::size_t wi) -> std::vector<double> {
873 this->checkWorkspaceIndices(wi);
874 return this->m_peakCenterWorkspace->x(wi).rawData();
875 };
876 } else {
877 std::stringstream errss;
878 errss << "One and only one in 'PeakCenters' (vector) and "
879 "'PeakCentersWorkspace' shall be given. "
880 << "'PeakCenters' has size " << m_peakCenters.size() << ", and name of peak center workspace "
881 << "is " << peakpswsname;
882 throw std::invalid_argument(errss.str());
883 }
884
885 return;
886}
887
888//----------------------------------------------------------------------------------------------
895 // check code integrity
896 if (m_numPeaksToFit == 0)
897 throw std::runtime_error("ProcessInputPeakTolerance() must be called after "
898 "ProcessInputPeakCenters()");
899
900 // peak tolerance
901 m_peakPosTolerances = getProperty(PropertyNames::POSITION_TOL);
902
903 if (m_peakPosTolerances.empty()) {
904 // case 2, 3, 4
905 m_peakPosTolerances.clear();
906 m_peakPosTolCase234 = true;
907 } else if (m_peakPosTolerances.size() == 1) {
908 // only 1 uniform peak position tolerance is defined: expand to all peaks
909 double peak_tol = m_peakPosTolerances[0];
910 m_peakPosTolerances.resize(m_numPeaksToFit, peak_tol);
911 } else if (m_peakPosTolerances.size() != m_numPeaksToFit) {
912 // not uniform but number of peaks does not match
913 g_log.error() << "number of peak position tolerance " << m_peakPosTolerances.size()
914 << " is not same as number of peaks " << m_numPeaksToFit << "\n";
915 throw std::runtime_error("Number of peak position tolerances and number of "
916 "peaks to fit are inconsistent.");
917 }
918
919 // set the minimum peak height to 0 (default value) if not specified or invalid
920 m_minPeakHeight = getProperty(PropertyNames::PEAK_MIN_HEIGHT);
922 m_minPeakHeight = 0.;
923
924 // PEAK_MIN_HEIGHT used to function as both "peak height" and "total count" checker.
925 // Now the "total count" is checked by PEAK_MIN_TOTAL_COUNT, so set it accordingly.
926 m_minPeakTotalCount = getProperty(PropertyNames::PEAK_MIN_TOTAL_COUNT);
929 else {
930 // set the minimum peak total count to 0 if not specified or invalid
933 }
934
935 // set the signal-to-noise threshold to zero (default value) if not specified or invalid
936 m_minSignalToNoiseRatio = getProperty(PropertyNames::PEAK_MIN_SIGNAL_TO_NOISE_RATIO);
939
940 // set the signal-to-sigma threshold to zero (default value) if not specified or invalid
941 m_minSignalToSigmaRatio = getProperty(PropertyNames::PEAK_MIN_SIGNAL_TO_SIGMA_RATIO);
944}
945
946//----------------------------------------------------------------------------------------------
953 // get a map for peak profile parameter name and parameter index
954 std::map<std::string, size_t> parname_index_map;
955 for (size_t iparam = 0; iparam < m_peakFunction->nParams(); ++iparam)
956 parname_index_map.insert(std::make_pair(m_peakFunction->parameterName(iparam), iparam));
957
958 // define peak parameter names (class variable) if using table
961
962 // map the input parameter names to parameter indexes
963 for (const auto &paramName : m_peakParamNames) {
964 auto locator = parname_index_map.find(paramName);
965 if (locator != parname_index_map.end()) {
966 m_initParamIndexes.emplace_back(locator->second);
967 } else {
968 // a parameter name that is not defined in the peak profile function. An
969 // out-of-range index is thus set to this
970 g_log.warning() << "Given peak parameter " << paramName
971 << " is not an allowed parameter of peak "
972 "function "
973 << m_peakFunction->name() << "\n";
974 m_initParamIndexes.emplace_back(m_peakFunction->nParams() * 10);
975 }
976 }
977
978 return;
979}
980
981//----------------------------------------------------------------------------------------------
984std::vector<std::shared_ptr<FitPeaksAlgorithm::PeakFitResult>> FitPeaks::fitPeaks() {
985 API::Progress prog(this, 0., 1., m_numPeaksToFit - 1);
986
989 std::vector<std::shared_ptr<FitPeaksAlgorithm::PeakFitResult>> fit_result_vector(m_numSpectraToFit);
990
991 const int nThreads = FrameworkManager::Instance().getNumOMPThreads();
992 size_t chunkSize = m_numSpectraToFit / nThreads;
993
994 std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> pre_check_result =
995 std::make_shared<FitPeaksAlgorithm::PeakFitPreCheckResult>();
996
997 PRAGMA_OMP(parallel for schedule(dynamic, 1) )
998 for (int ithread = 0; ithread < nThreads; ithread++) {
1000 auto iws_begin = m_startWorkspaceIndex + chunkSize * static_cast<size_t>(ithread);
1001 auto iws_end = (ithread == nThreads - 1) ? m_stopWorkspaceIndex + 1 : iws_begin + chunkSize;
1002
1003 // vector to store fit params for last good fit to each peak
1004 std::vector<std::vector<double>> lastGoodPeakParameters(m_numPeaksToFit,
1005 std::vector<double>(m_peakFunction->nParams(), 0.0));
1006 // track which spectrum index last successfully fitted each peak
1007 std::vector<size_t> lastGoodPeakSpectra(m_numPeaksToFit, 0);
1008
1009 for (auto wi = iws_begin; wi < iws_end; ++wi) {
1010 // peaks to fit
1011 std::vector<double> expected_peak_centers = m_getExpectedPeakPositions(static_cast<size_t>(wi));
1012
1013 // initialize output for this
1014 size_t numfuncparams = m_peakFunction->nParams() + m_bkgdFunction->nParams();
1015 std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> fit_result =
1016 std::make_shared<FitPeaksAlgorithm::PeakFitResult>(m_numPeaksToFit, numfuncparams);
1017
1018 std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> spectrum_pre_check_result =
1019 std::make_shared<FitPeaksAlgorithm::PeakFitPreCheckResult>();
1020
1021 fitSpectrumPeaks(static_cast<size_t>(wi), expected_peak_centers, fit_result, lastGoodPeakParameters,
1022 lastGoodPeakSpectra, spectrum_pre_check_result);
1023
1024 PARALLEL_CRITICAL(FindPeaks_WriteOutput) {
1025 writeFitResult(static_cast<size_t>(wi), expected_peak_centers, fit_result);
1026 fit_result_vector[wi - m_startWorkspaceIndex] = fit_result;
1027 *pre_check_result += *spectrum_pre_check_result;
1028 }
1029 prog.report();
1030 }
1032 }
1034 logNoOffset(5 /*notice*/, pre_check_result->getReport());
1035 return fit_result_vector;
1036}
1037
1038namespace {
1039// Forward declarations
1040bool estimateBackgroundParameters(const Histogram &histogram, const std::pair<size_t, size_t> &peak_window,
1041 const API::IBackgroundFunction_sptr &bkgd_function);
1042void reduceByBackground(const API::IBackgroundFunction_sptr &bkgd_func, const std::vector<double> &vec_x,
1043 std::vector<double> &vec_y);
1044template <typename vector_like>
1045void rangeToIndexBounds(const vector_like &vecx, const double range_left, const double range_right, size_t &left_index,
1046 size_t &right_index);
1047
1049std::vector<std::string> supported_peak_profiles{"Gaussian", "Lorentzian", "PseudoVoigt", "Voigt",
1050 "BackToBackExponential"};
1051
1052//----------------------------------------------------------------------------------------------
1057double estimateBackgroundNoise(const std::vector<double> &vec_y) {
1058 // peak window must have a certain minimum number of data points necessary to do the statistics
1059 size_t half_number_of_bkg_datapoints{5};
1060 if (vec_y.size() < 2 * half_number_of_bkg_datapoints + 3 /*a magic number*/)
1061 return DBL_MIN; // can't estimate the noise
1062
1063 // the specified number of left-most and right-most data points in the peak window are assumed to represent
1064 // background. Combine these data points into a single vector
1065 std::vector<double> vec_bkg;
1066 vec_bkg.resize(2 * half_number_of_bkg_datapoints);
1067 std::copy(vec_y.begin(), vec_y.begin() + half_number_of_bkg_datapoints, vec_bkg.begin());
1068 std::copy(vec_y.end() - half_number_of_bkg_datapoints, vec_y.end(), vec_bkg.begin() + half_number_of_bkg_datapoints);
1069
1070 // estimate the noise as the standard deviation of the combined background vector, but without outliers
1071 std::vector<double> zscore_vec = Kernel::getZscore(vec_bkg);
1072 std::vector<double> vec_bkg_no_outliers;
1073 vec_bkg_no_outliers.resize(vec_bkg.size());
1074 double zscore_crit = 3.; // using three-sigma rule
1075 for (size_t ii = 0; ii < vec_bkg.size(); ii++) {
1076 if (zscore_vec[ii] <= zscore_crit)
1077 vec_bkg_no_outliers.push_back(vec_bkg[ii]);
1078 }
1079
1080 if (vec_bkg_no_outliers.size() < half_number_of_bkg_datapoints)
1081 return DBL_MIN; // can't estimate the noise
1082
1083 auto intensityStatistics = Kernel::getStatistics(vec_bkg_no_outliers, StatOptions::CorrectedStdDev);
1084 return intensityStatistics.standard_deviation;
1085}
1086
1087//----------------------------------------------------------------------------------------------
1095template <typename vector_like>
1096void rangeToIndexBounds(const vector_like &elems, const double range_left, const double range_right, size_t &left_index,
1097 size_t &right_index) {
1098 const auto left_iter = std::lower_bound(elems.cbegin(), elems.cend(), range_left);
1099 const auto right_iter = std::upper_bound(elems.cbegin(), elems.cend(), range_right);
1100
1101 left_index = std::distance(elems.cbegin(), left_iter);
1102 right_index = std::distance(elems.cbegin(), right_iter);
1103 right_index = std::min(right_index, elems.size() - 1);
1104}
1105
1106//----------------------------------------------------------------------------------------------
1112void reduceByBackground(const API::IBackgroundFunction_sptr &bkgd_func, const std::vector<double> &vec_x,
1113 std::vector<double> &vec_y) {
1114 // calculate the background
1115 FunctionDomain1DVector vectorx(vec_x.begin(), vec_x.end());
1116 FunctionValues vector_bkgd(vectorx);
1117 bkgd_func->function(vectorx, vector_bkgd);
1118
1119 // subtract the background from the supplied data
1120 for (size_t i = 0; i < vec_y.size(); ++i) {
1121 (vec_y)[i] -= vector_bkgd[i];
1122 // Note, E is not changed here
1123 }
1124}
1125
1126//----------------------------------------------------------------------------------------------
1130class LoggingOffsetSentry {
1131public:
1132 LoggingOffsetSentry(Algorithm *const alg) : m_alg(alg) {
1135 }
1136 ~LoggingOffsetSentry() { m_alg->setLoggingOffset(m_loggingOffset); }
1137
1138private:
1141};
1142} // namespace
1143
1144//----------------------------------------------------------------------------------------------
1147void FitPeaks::fitSpectrumPeaks(size_t wi, const std::vector<double> &expected_peak_centers,
1148 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result,
1149 std::vector<std::vector<double>> &lastGoodPeakParameters,
1150 std::vector<size_t> &lastGoodPeakSpectra,
1151 const std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> &pre_check_result) {
1152 assert(fit_result->getNumberPeaks() == m_numPeaksToFit);
1153 pre_check_result->setNumberOfSubmittedSpectrumPeaks(m_numPeaksToFit);
1154 // if the whole spectrum has low count, do not fit any peaks for that spectrum
1156 for (size_t i = 0; i < m_numPeaksToFit; ++i)
1157 fit_result->setBadRecord(i, -1.);
1158 pre_check_result->setNumberOfSpectrumPeaksWithLowCount(m_numPeaksToFit);
1159 return;
1160 }
1161
1162 // Set up sub algorithm Fit for peak and background
1163 IAlgorithm_sptr peak_fitter; // both peak and background (combo)
1164 try {
1165 peak_fitter = createChildAlgorithm("Fit", -1, -1, false);
1166 } catch (Exception::NotFoundError &) {
1167 std::stringstream errss;
1168 errss << "The FitPeak algorithm requires the CurveFitting library";
1169 g_log.error(errss.str());
1170 throw std::runtime_error(errss.str());
1171 }
1172
1173 // Clone background function
1174 IBackgroundFunction_sptr bkgdfunction = std::dynamic_pointer_cast<API::IBackgroundFunction>(m_bkgdFunction->clone());
1175
1176 // set up properties of algorithm (reference) 'Fit'
1177 peak_fitter->setProperty("Minimizer", m_minimizer);
1178 peak_fitter->setProperty("CostFunction", m_costFunction);
1179 peak_fitter->setProperty("CalcErrors", true);
1180
1181 const double x0 = m_inputMatrixWS->histogram(wi).x().front();
1182 const double xf = m_inputMatrixWS->histogram(wi).x().back();
1183
1184 // index of previous peak in same spectrum (initially invalid)
1185 size_t prev_peak_index = m_numPeaksToFit;
1186 bool neighborPeakSameSpectrum = false;
1187 size_t number_of_out_of_range_peaks{0};
1188 for (size_t fit_index = 0; fit_index < m_numPeaksToFit; ++fit_index) {
1189 // convert fit index to peak index (in ascending order)
1190 size_t peak_index(fit_index);
1192 peak_index = m_numPeaksToFit - fit_index - 1;
1193
1194 // reset the background function
1195 for (size_t i = 0; i < bkgdfunction->nParams(); ++i)
1196 bkgdfunction->setParameter(i, 0.);
1197
1198 double expected_peak_pos = expected_peak_centers[peak_index];
1199 std::pair<double, double> peak_window_i = m_getPeakFitWindow(wi, peak_index);
1200
1201 // clone peak function for each peak (need to do this so can
1202 // set center and calc any parameters from xml)
1203 auto peakfunction = std::dynamic_pointer_cast<API::IPeakFunction>(m_peakFunction->clone());
1204 peakfunction->setCentre(expected_peak_pos);
1205 peakfunction->setMatrixWorkspace(m_inputMatrixWS, wi, peak_window_i.first, peak_window_i.second);
1206
1207 std::map<size_t, double> keep_values;
1208 for (size_t ipar = 0; ipar < peakfunction->nParams(); ++ipar) {
1209 if (peakfunction->isFixed(ipar)) {
1210 // save value of these parameters which have just been calculated
1211 // if they were set to be fixed (e.g. for the B2Bexp this would
1212 // typically be A and B but not Sigma)
1213 keep_values[ipar] = peakfunction->getParameter(ipar);
1214 // let them be free to fit as these are typically refined from a
1215 // focussed bank
1216 peakfunction->unfix(ipar);
1217 }
1218 }
1219
1220 // Determine whether to set starting parameter from fitted value
1221 // of same peak but different spectrum
1222 bool samePeakCrossSpectrum = (lastGoodPeakParameters[peak_index].size() >
1223 static_cast<size_t>(std::count_if(lastGoodPeakParameters[peak_index].begin(),
1224 lastGoodPeakParameters[peak_index].end(),
1225 [&](auto const &val) { return val <= 1e-10; })));
1226
1227 // Check whether current spectrum's pixel (detector ID) is close to the
1228 // spectrum that last successfully fitted this peak.
1229 try {
1230 if (wi > 0 && samePeakCrossSpectrum) {
1231 size_t lastGoodWi = lastGoodPeakSpectra[peak_index];
1232 std::shared_ptr<const Geometry::Detector> pdetector =
1233 std::dynamic_pointer_cast<const Geometry::Detector>(m_inputMatrixWS->getDetector(lastGoodWi));
1234 std::shared_ptr<const Geometry::Detector> cdetector =
1235 std::dynamic_pointer_cast<const Geometry::Detector>(m_inputMatrixWS->getDetector(wi));
1236
1237 // If they do have detector ID
1238 if (pdetector && cdetector) {
1239 auto prev_id = pdetector->getID();
1240 auto curr_id = cdetector->getID();
1241 if (prev_id + 1 != curr_id)
1242 samePeakCrossSpectrum = false;
1243 } else {
1244 samePeakCrossSpectrum = false;
1245 }
1246
1247 } else {
1248 // first spectrum in the workspace: no peak's fitting result to copy
1249 // from
1250 samePeakCrossSpectrum = false;
1251 }
1252 } catch (const std::runtime_error &) {
1253 // workspace does not have detector ID set: there is no guarantee that the
1254 // adjacent spectra can have similar peak profiles
1255 samePeakCrossSpectrum = false;
1256 }
1257
1258 // Set starting values of the peak function
1259 if (samePeakCrossSpectrum) { // somePeakFit
1260 // Get from local best result
1261 for (size_t i = 0; i < peakfunction->nParams(); ++i) {
1262 peakfunction->setParameter(i, lastGoodPeakParameters[peak_index][i]);
1263 }
1264 } else if (neighborPeakSameSpectrum && m_copyLastGoodPeakParameters) {
1265 // set the peak parameters from last good fit from ANY peak in the spectrum
1266 for (size_t i = 0; i < peakfunction->nParams(); ++i) {
1267 peakfunction->setParameter(i, lastGoodPeakParameters[prev_peak_index][i]);
1268 }
1269 }
1270
1271 // reset center though - don't know before hand which element this is
1272 peakfunction->setCentre(expected_peak_pos);
1273
1274 // reset value of parameters that were fixed (but are now free to vary)
1275 for (const auto &[ipar, value] : keep_values) {
1276 peakfunction->setParameter(ipar, value);
1277 }
1278
1279 double cost(DBL_MAX);
1280 if (expected_peak_pos <= x0 || expected_peak_pos >= xf) {
1281 // out of range and there won't be any fit
1282 peakfunction->setIntensity(0);
1283 number_of_out_of_range_peaks++;
1284 } else {
1285 // Decide whether to estimate peak width by observation
1286 // If no peaks fitted in the same or cross spectrum then the user supplied
1287 // parameters will be used if present and the width will not be estimated
1288 // (note this will overwrite parameter values caluclated from
1289 // Parameters.xml)
1290 // When CopyLastGoodPeakParameters is disabled, treat each peak as if it
1291 // has no fitted neighbour so that user-specified initial values are
1292 // re-applied as the baseline rather than falling back to function defaults.
1293 auto useUserSpecifedIfGiven =
1294 !(samePeakCrossSpectrum || (neighborPeakSameSpectrum && m_copyLastGoodPeakParameters));
1295 bool observe_peak_width = decideToEstimatePeakParams(useUserSpecifedIfGiven, peakfunction);
1296
1297 if (observe_peak_width && m_peakWidthEstimateApproach == EstimatePeakWidth::NoEstimation) {
1298 g_log.warning("Peak width can be estimated as ZERO. The result can be wrong");
1299 }
1300
1301 // do fitting with peak and background function (no analysis at this point)
1302 std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> peak_pre_check_result =
1303 std::make_shared<FitPeaksAlgorithm::PeakFitPreCheckResult>();
1304 cost = fitIndividualPeak(wi, peak_fitter, expected_peak_pos, peak_window_i, observe_peak_width, peakfunction,
1305 bkgdfunction, peak_pre_check_result);
1306 if (peak_pre_check_result->isIndividualPeakRejected())
1307 fit_result->setBadRecord(peak_index, -1.);
1308
1309 if (m_minSignalToSigmaRatio > 0) {
1310 if (calculateSignalToSigmaRatio(wi, peak_window_i, peakfunction) < m_minSignalToSigmaRatio) {
1311 fit_result->setBadRecord(peak_index, -1.);
1312 cost = DBL_MAX;
1313 }
1314 }
1315
1316 *pre_check_result += *peak_pre_check_result; // keep track of the rejection count within the spectrum
1317 }
1318 pre_check_result->setNumberOfOutOfRangePeaks(number_of_out_of_range_peaks);
1319
1320 // process fitting result
1321 FitPeaksAlgorithm::FitFunction fit_function;
1322 fit_function.peakfunction = peakfunction;
1323 fit_function.bkgdfunction = bkgdfunction;
1324
1325 auto good_fit = processSinglePeakFitResult(wi, peak_index, cost, expected_peak_centers, fit_function,
1326 fit_result); // sets the record
1327
1328 if (good_fit) {
1329 // reset the flag such that there is at a peak fit in this spectrum
1330 neighborPeakSameSpectrum = true;
1331 prev_peak_index = peak_index;
1332 // copy values and record which spectrum they came from
1333 for (size_t i = 0; i < lastGoodPeakParameters[peak_index].size(); ++i) {
1334 lastGoodPeakParameters[peak_index][i] = peakfunction->getParameter(i);
1335 }
1336 lastGoodPeakSpectra[peak_index] = wi;
1337 }
1338 }
1339
1340 return;
1341}
1342
1343//----------------------------------------------------------------------------------------------
1352bool FitPeaks::decideToEstimatePeakParams(const bool firstPeakInSpectrum,
1353 const API::IPeakFunction_sptr &peak_function) {
1354 // should observe the peak width if the user didn't supply all of the peak
1355 // function parameters
1356 bool observe_peak_shape(m_initParamIndexes.size() != peak_function->nParams());
1357
1358 if (!m_initParamIndexes.empty()) {
1359 // user specifies starting value of peak parameters
1360 if (firstPeakInSpectrum) {
1361 // set the parameter values in a vector and loop over it
1362 // first peak. using the user-specified value
1363 for (size_t i = 0; i < m_initParamIndexes.size(); ++i) {
1364 const size_t param_index = m_initParamIndexes[i];
1365 const double param_value = m_initParamValues[i];
1366 peak_function->setParameter(param_index, param_value);
1367 }
1368 } else {
1369 // using the fitted paramters from the previous fitting result
1370 // do noting
1371 }
1372 } else {
1373 // no previously defined peak parameters: observation is thus required
1374 observe_peak_shape = true;
1375 }
1376
1377 return observe_peak_shape;
1378}
1379
1380//----------------------------------------------------------------------------------------------
1392bool FitPeaks::processSinglePeakFitResult(size_t wsindex, size_t peakindex, const double cost,
1393 const std::vector<double> &expected_peak_positions,
1394 const FitPeaksAlgorithm::FitFunction &fitfunction,
1395 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result) {
1396 // determine peak position tolerance
1397 double postol(DBL_MAX);
1398 bool case23(false);
1399 if (m_peakPosTolCase234) {
1400 // peak tolerance is not defined
1401 if (m_numPeaksToFit == 1) {
1402 // case (d) one peak only
1403 postol = m_inputMatrixWS->histogram(wsindex).x().back() - m_inputMatrixWS->histogram(wsindex).x().front();
1404 } else {
1405 // case b and c: more than 1 peaks without defined peak tolerance
1406 case23 = true;
1407 }
1408 } else {
1409 // user explicitly specified
1410 if (peakindex >= m_peakPosTolerances.size())
1411 throw std::runtime_error("Peak tolerance out of index");
1412 postol = m_peakPosTolerances[peakindex];
1413 }
1414
1415 // get peak position and analyze the fitting is good or not by various
1416 // criteria
1417 auto peak_pos = fitfunction.peakfunction->centre();
1418 auto peak_fwhm = fitfunction.peakfunction->fwhm();
1419 bool good_fit(false);
1420 if ((cost < 0) || (cost >= DBL_MAX - 1.) || std::isnan(cost)) {
1421 // unphysical cost function value
1422 peak_pos = -4;
1423 } else if (fitfunction.peakfunction->height() < m_minPeakHeight) {
1424 // peak height is under minimum request
1425 peak_pos = -3;
1426 } else if (case23) {
1427 // case b and c to check peak position without defined peak tolerance
1428 std::pair<double, double> fitwindow = m_getPeakFitWindow(wsindex, peakindex);
1429 if (fitwindow.first < fitwindow.second) {
1430 // peak fit window is specified or calculated: use peak window as position
1431 // tolerance
1432 if (peak_pos < fitwindow.first || peak_pos > fitwindow.second) {
1433 // peak is out of fit window
1434 peak_pos = -2;
1435 g_log.debug() << "Peak position " << peak_pos << " is out of fit "
1436 << "window boundary " << fitwindow.first << ", " << fitwindow.second << "\n";
1437 } else if (peak_fwhm > (fitwindow.second - fitwindow.first)) {
1438 // peak is too wide or window is too small
1439 peak_pos = -2.25;
1440 g_log.debug() << "Peak position " << peak_pos << " has fwhm "
1441 << "wider than the fit window " << fitwindow.second - fitwindow.first << "\n";
1442 } else {
1443 good_fit = true;
1444 }
1445 } else {
1446 // use the 1/2 distance to neiboring peak without defined peak window
1447 double left_bound(-1);
1448 if (peakindex > 0)
1449 left_bound = 0.5 * (expected_peak_positions[peakindex] - expected_peak_positions[peakindex - 1]);
1450 double right_bound(-1);
1451 if (peakindex < m_numPeaksToFit - 1)
1452 right_bound = 0.5 * (expected_peak_positions[peakindex + 1] - expected_peak_positions[peakindex]);
1453 if (left_bound < 0)
1454 left_bound = right_bound;
1455 if (right_bound < left_bound)
1456 right_bound = left_bound;
1457 if (left_bound < 0 || right_bound < 0)
1458 throw std::runtime_error("Code logic error such that left or right "
1459 "boundary of peak position is negative.");
1460 if (peak_pos < left_bound || peak_pos > right_bound) {
1461 peak_pos = -2.5;
1462 } else if (peak_fwhm > (right_bound - left_bound)) {
1463 // peak is too wide or window is too small
1464 peak_pos = -2.75;
1465 g_log.debug() << "Peak position " << peak_pos << " has fwhm "
1466 << "wider than the fit window " << right_bound - left_bound << "\n";
1467 } else {
1468 good_fit = true;
1469 }
1470 }
1471 } else if (fabs(fitfunction.peakfunction->centre() - expected_peak_positions[peakindex]) > postol) {
1472 // peak center is not within tolerance
1473 peak_pos = -5;
1474 g_log.debug() << "Peak position difference "
1475 << fabs(fitfunction.peakfunction->centre() - expected_peak_positions[peakindex])
1476 << " is out of range of tolerance: " << postol << "\n";
1477 } else {
1478 // all criteria are passed
1479 good_fit = true;
1480 }
1481
1482 // set cost function to DBL_MAX if fitting is bad
1483 double adjust_cost(cost);
1484 if (!good_fit) {
1485 // set the cost function value to DBL_MAX
1486 adjust_cost = DBL_MAX;
1487 }
1488
1489 // reset cost
1490 if (adjust_cost > DBL_MAX - 1) {
1491 fitfunction.peakfunction->setIntensity(0);
1492 }
1493
1494 // chi2
1495 fit_result->setRecord(peakindex, adjust_cost, peak_pos, fitfunction);
1496
1497 return good_fit;
1498}
1499
1500//----------------------------------------------------------------------------------------------
1506void FitPeaks::calculateFittedPeaks(const std::vector<std::shared_ptr<FitPeaksAlgorithm::PeakFitResult>> &fit_results) {
1507 // check
1508 if (!m_fittedParamTable)
1509 throw std::runtime_error("No parameters");
1510
1511 const size_t num_peakfunc_params = m_peakFunction->nParams();
1512 const size_t num_bkgdfunc_params = m_bkgdFunction->nParams();
1513
1515 for (int64_t iiws = m_startWorkspaceIndex; iiws <= static_cast<int64_t>(m_stopWorkspaceIndex); ++iiws) {
1517 std::size_t iws = static_cast<std::size_t>(iiws);
1518 // get a copy of peak function and background function
1519 IPeakFunction_sptr peak_function = std::dynamic_pointer_cast<IPeakFunction>(m_peakFunction->clone());
1520 IBackgroundFunction_sptr bkgd_function = std::dynamic_pointer_cast<IBackgroundFunction>(m_bkgdFunction->clone());
1521 std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> fit_result_i = fit_results[iws - m_startWorkspaceIndex];
1522 // FIXME - This is a just a pure check
1523 if (!fit_result_i)
1524 throw std::runtime_error("There is something wroing with PeakFitResult vector!");
1525
1526 for (size_t ipeak = 0; ipeak < m_numPeaksToFit; ++ipeak) {
1527 // get and set the peak function parameters
1528 const double chi2 = fit_result_i->getCost(ipeak);
1529 if (chi2 > 10.e10)
1530 continue;
1531
1532 for (size_t iparam = 0; iparam < num_peakfunc_params; ++iparam)
1533 peak_function->setParameter(iparam, fit_result_i->getParameterValue(ipeak, iparam));
1534 for (size_t iparam = 0; iparam < num_bkgdfunc_params; ++iparam)
1535 bkgd_function->setParameter(iparam, fit_result_i->getParameterValue(ipeak, num_peakfunc_params + iparam));
1536 // use domain and function to calcualte
1537 // get the range of start and stop to construct a function domain
1538 const auto &vec_x = m_fittedPeakWS->points(iws);
1539 std::pair<double, double> peakwindow = m_getPeakFitWindow(iws, ipeak);
1540 auto start_x_iter = std::lower_bound(vec_x.begin(), vec_x.end(), peakwindow.first);
1541 auto stop_x_iter = std::lower_bound(vec_x.begin(), vec_x.end(), peakwindow.second);
1542
1543 if (start_x_iter == stop_x_iter)
1544 throw std::runtime_error("Range size is zero in calculateFittedPeaks");
1545
1546 FunctionDomain1DVector domain(start_x_iter, stop_x_iter);
1547 FunctionValues values(domain);
1548 CompositeFunction_sptr comp_func = std::make_shared<API::CompositeFunction>();
1549 comp_func->addFunction(peak_function);
1550 comp_func->addFunction(bkgd_function);
1551 comp_func->function(domain, values);
1552
1553 // copy over the values
1554 std::size_t istart = static_cast<size_t>(start_x_iter - vec_x.begin());
1555 std::size_t istop = static_cast<size_t>(stop_x_iter - vec_x.begin());
1556 for (std::size_t yindex = istart; yindex < istop; ++yindex) {
1557 m_fittedPeakWS->dataY(iws)[yindex] = values.getCalculated(yindex - istart);
1558 }
1559 } // END-FOR (ipeak)
1561 } // END-FOR (iws)
1563
1564 return;
1565}
1566
1567double FitPeaks::calculateSignalToSigmaRatio(const size_t &iws, const std::pair<double, double> &peakWindow,
1568 const API::IPeakFunction_sptr &peakFunction) {
1569 const auto &vecX = m_inputMatrixWS->points(iws);
1570 auto startX = std::lower_bound(vecX.begin(), vecX.end(), peakWindow.first);
1571 auto stopX = std::lower_bound(vecX.begin(), vecX.end(), peakWindow.second);
1572
1573 FunctionDomain1DVector domain(startX, stopX);
1574 FunctionValues values(domain);
1575
1576 peakFunction->function(domain, values);
1577 auto peakValues = values.toVector();
1578
1579 const auto &errors = m_inputMatrixWS->readE(iws);
1580 auto startE = errors.begin() + (startX - vecX.begin());
1581 auto stopE = errors.begin() + (stopX - vecX.begin());
1582 std::vector<double> peakErrors(startE, stopE);
1583
1584 double peakSum = std::accumulate(peakValues.cbegin(), peakValues.cend(), 0.0);
1585 double sigma = sqrt(std::accumulate(peakErrors.cbegin(), peakErrors.cend(), 0.0, VectorHelper::SumSquares<double>()));
1586
1587 return peakSum / ((sigma == 0) ? 1 : sigma);
1588}
1589
1590namespace {
1591bool estimateBackgroundParameters(const Histogram &histogram, const std::pair<size_t, size_t> &peak_window,
1592 const API::IBackgroundFunction_sptr &bkgd_function) {
1593 // for estimating background parameters
1594 // 0 = constant, 1 = linear
1595 const auto POLYNOMIAL_ORDER = std::min<size_t>(1, bkgd_function->nParams());
1596
1597 if (peak_window.first >= peak_window.second)
1598 throw std::runtime_error("Invalid peak window");
1599
1600 // reset the background function
1601 const auto nParams = bkgd_function->nParams();
1602 for (size_t i = 0; i < nParams; ++i)
1603 bkgd_function->setParameter(i, 0.);
1604
1605 // 10 is a magic number that worked in a variety of situations
1606 const size_t iback_start = peak_window.first + 10;
1607 const size_t iback_stop = peak_window.second - 10;
1608
1609 // use the simple way to find linear background
1610 // there aren't enough bins in the window to try to estimate so just leave the
1611 // estimate at zero
1612 if (iback_start < iback_stop) {
1613 double bkgd_a0{0.}; // will be fit
1614 double bkgd_a1{0.}; // may be fit
1615 double bkgd_a2{0.}; // will be ignored
1616 double chisq{DBL_MAX}; // how well the fit worked
1617 HistogramData::estimateBackground(POLYNOMIAL_ORDER, histogram, peak_window.first, peak_window.second, iback_start,
1618 iback_stop, bkgd_a0, bkgd_a1, bkgd_a2, chisq);
1619 // update the background function with the result
1620 bkgd_function->setParameter(0, bkgd_a0);
1621 if (nParams > 1)
1622 bkgd_function->setParameter(1, bkgd_a1);
1623 // quadratic term is always estimated to be zero
1624
1625 // TODO: return false if chisq is too large
1626 return true;
1627 }
1628
1629 return false; // too few data points for the fit
1630}
1631} // anonymous namespace
1632
1633//----------------------------------------------------------------------------------------------
1639bool FitPeaks::isObservablePeakProfile(const std::string &peakprofile) {
1640 return (std::find(supported_peak_profiles.begin(), supported_peak_profiles.end(), peakprofile) !=
1641 supported_peak_profiles.end());
1642}
1643
1644//----------------------------------------------------------------------------------------------
1647bool FitPeaks::fitBackground(const size_t &ws_index, const std::pair<double, double> &fit_window,
1648 const double &expected_peak_pos, const API::IBackgroundFunction_sptr &bkgd_func) {
1649 constexpr size_t MIN_POINTS{10}; // TODO explain why 10
1650
1651 // find out how to fit background
1652 const auto &points = m_inputMatrixWS->histogram(ws_index).points();
1653 size_t start_index = findXIndex(points.rawData(), fit_window.first);
1654 size_t expected_peak_index = findXIndex(points.rawData(), expected_peak_pos, start_index);
1655 size_t stop_index = findXIndex(points.rawData(), fit_window.second, expected_peak_index);
1656
1657 // treat 5 as a magic number - TODO explain why
1658 bool good_fit(false);
1659 if (expected_peak_index - start_index > MIN_POINTS && stop_index - expected_peak_index > MIN_POINTS) {
1660 // enough data points left for multi-domain fitting
1661 // set a smaller fit window
1662 const std::pair<double, double> vec_min{fit_window.first, points[expected_peak_index + 5]};
1663 const std::pair<double, double> vec_max{points[expected_peak_index - 5], fit_window.second};
1664
1665 // reset background function value
1666 for (size_t n = 0; n < bkgd_func->nParams(); ++n)
1667 bkgd_func->setParameter(n, 0);
1668
1669 double chi2 = fitFunctionMD(bkgd_func, m_inputMatrixWS, ws_index, vec_min, vec_max);
1670
1671 // process
1672 if (chi2 < DBL_MAX - 1) {
1673 good_fit = true;
1674 }
1675
1676 } else {
1677 // fit as a single domain function. check whether the result is good or bad
1678
1679 // TODO FROM HERE!
1680 g_log.debug() << "Don't know what to do with background fitting with single "
1681 << "domain function! " << (expected_peak_index - start_index) << " points to the left "
1682 << (stop_index - expected_peak_index) << " points to the right\n";
1683 }
1684
1685 return good_fit;
1686}
1687
1688//----------------------------------------------------------------------------------------------
1691double FitPeaks::fitIndividualPeak(size_t wi, const API::IAlgorithm_sptr &fitter, const double expected_peak_center,
1692 const std::pair<double, double> &fitwindow, const bool estimate_peak_width,
1693 const API::IPeakFunction_sptr &peakfunction,
1694 const API::IBackgroundFunction_sptr &bkgdfunc,
1695 const std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> &pre_check_result) {
1696 pre_check_result->setNumberOfSubmittedIndividualPeaks(1);
1697 double cost(DBL_MAX);
1698
1699 // make sure the number of data points satisfies the number of fitting parameters plus a magic cushion of 2.
1700 size_t min_required_datapoints{peakfunction->nParams() + bkgdfunc->nParams() + 2};
1701 size_t number_of_datapoints = histRangeToDataPointCount(wi, fitwindow);
1702 if (number_of_datapoints < min_required_datapoints) {
1703 pre_check_result->setNumberOfPeaksWithNotEnoughDataPoints(1);
1704 return cost;
1705 }
1706
1707 // check the number of counts in the peak window
1708 if (m_minPeakTotalCount >= 0.0 && numberCounts(wi, fitwindow) <= m_minPeakTotalCount) {
1709 pre_check_result->setNumberOfIndividualPeaksWithLowCount(1);
1710 return cost;
1711 }
1712
1713 // exclude a peak with a low signal-to-noise ratio
1714 if (m_minSignalToNoiseRatio > 0.0 && calculateSignalToNoiseRatio(wi, fitwindow, bkgdfunc) < m_minSignalToNoiseRatio) {
1715 pre_check_result->setNumberOfPeaksWithLowSignalToNoise(1);
1716 return cost;
1717 }
1718
1719 if (m_highBackground) {
1720 // fit peak with high background!
1721 cost = fitFunctionHighBackground(fitter, fitwindow, wi, expected_peak_center, estimate_peak_width, peakfunction,
1722 bkgdfunc);
1723 } else {
1724 // fit peak and background
1725 cost = fitFunctionSD(fitter, peakfunction, bkgdfunc, m_inputMatrixWS, wi, fitwindow, expected_peak_center,
1726 estimate_peak_width, true);
1727 }
1728
1729 return cost;
1730}
1731
1732//----------------------------------------------------------------------------------------------
1739 const API::IBackgroundFunction_sptr &bkgd_function,
1740 const API::MatrixWorkspace_sptr &dataws, size_t wsindex,
1741 const std::pair<double, double> &peak_range, const double &expected_peak_center,
1742 bool estimate_peak_width, bool estimate_background) {
1743 std::stringstream errorid;
1744 errorid << "(WorkspaceIndex=" << wsindex << " PeakCentre=" << expected_peak_center << ")";
1745
1746 // validate peak window
1747 if (peak_range.first >= peak_range.second) {
1748 std::stringstream msg;
1749 msg << "Invalid peak window: xmin>xmax (" << peak_range.first << ", " << peak_range.second << ")" << errorid.str();
1750 throw std::runtime_error(msg.str());
1751 }
1752
1753 // determine the peak window in terms of vector indexes
1754 const auto &histogram = dataws->histogram(wsindex);
1755 const auto &vector_x = histogram.points();
1756 const auto start_index = findXIndex(vector_x, peak_range.first);
1757 const auto stop_index = findXIndex(vector_x, peak_range.second, start_index);
1758 if (start_index == stop_index)
1759 throw std::runtime_error("Range size is zero in fitFunctionSD");
1760 std::pair<size_t, size_t> peak_index_window = std::make_pair(start_index, stop_index);
1761
1762 // Estimate background
1763 if (estimate_background) {
1764 if (!estimateBackgroundParameters(histogram, peak_index_window, bkgd_function)) {
1765 return DBL_MAX;
1766 }
1767 }
1768
1769 // Estimate peak profile parameter
1770 peak_function->setCentre(expected_peak_center); // set expected position first
1771 int result = estimatePeakParameters(histogram, peak_index_window, peak_function, bkgd_function, estimate_peak_width,
1773
1774 if (result != GOOD) {
1775 peak_function->setCentre(expected_peak_center);
1776 if (result == NOSIGNAL || result == LOWPEAK) {
1777 return DBL_MAX; // exit early - don't fit
1778 }
1779 }
1780
1781 // Create the composition function
1782 CompositeFunction_sptr comp_func = std::make_shared<API::CompositeFunction>();
1783 comp_func->addFunction(peak_function);
1784 comp_func->addFunction(bkgd_function);
1785 IFunction_sptr fitfunc = std::dynamic_pointer_cast<IFunction>(comp_func);
1786
1787 // Set the properties
1788 fit->setProperty("Function", fitfunc);
1789 fit->setProperty("InputWorkspace", dataws);
1790 fit->setProperty("WorkspaceIndex", static_cast<int>(wsindex));
1791 fit->setProperty("MaxIterations", m_fitIterations); // magic number
1792 fit->setProperty("StartX", peak_range.first);
1793 fit->setProperty("EndX", peak_range.second);
1794 fit->setProperty("IgnoreInvalidData", true);
1795
1797 // set up a constraint on peak position
1798 double peak_center = peak_function->centre();
1799 double peak_width = peak_function->fwhm();
1800 std::stringstream peak_center_constraint;
1801 peak_center_constraint << (peak_center - 0.5 * peak_width) << " < f0." << peak_function->getCentreParameterName()
1802 << " < " << (peak_center + 0.5 * peak_width);
1803 fit->setProperty("Constraints", peak_center_constraint.str());
1804 }
1805
1806 // Execute fit and get result of fitting background
1807 g_log.debug() << "[E1201] FitSingleDomain Before fitting, Fit function: " << fit->asString() << "\n";
1808 errorid << " starting function [" << comp_func->asString() << "]";
1809 try {
1810 fit->execute();
1811 g_log.debug() << "[E1202] FitSingleDomain After fitting, Fit function: " << fit->asString() << "\n";
1812
1813 if (!fit->isExecuted()) {
1814 g_log.warning() << "Fitting peak SD (single domain) failed to execute. " + errorid.str();
1815 return DBL_MAX;
1816 }
1817 } catch (std::invalid_argument &e) {
1818 errorid << ": " << e.what();
1819 g_log.warning() << "\nWhile fitting " + errorid.str();
1820 return DBL_MAX; // probably the wrong thing to do
1821 }
1822
1823 // Retrieve result
1824 std::string fitStatus = fit->getProperty("OutputStatus");
1825 double chi2{std::numeric_limits<double>::max()};
1826 if (fitStatus == "success") {
1827 chi2 = fit->getProperty("OutputChi2overDoF");
1828 }
1829
1830 return chi2;
1831}
1832
1833//----------------------------------------------------------------------------------------------
1835 const size_t wsindex, const std::pair<double, double> &vec_xmin,
1836 const std::pair<double, double> &vec_xmax) {
1837 // Note: after testing it is found that multi-domain Fit cannot be reused
1839 try {
1840 fit = createChildAlgorithm("Fit", -1, -1, false);
1841 } catch (Exception::NotFoundError &) {
1842 std::stringstream errss;
1843 errss << "The FitPeak algorithm requires the CurveFitting library";
1844 throw std::runtime_error(errss.str());
1845 }
1846 // set up background fit instance
1847 fit->setProperty("Minimizer", m_minimizer);
1848 fit->setProperty("CostFunction", m_costFunction);
1849 fit->setProperty("CalcErrors", true);
1850
1851 // This use multi-domain; but does not know how to set up IFunction_sptr
1852 // fitfunc,
1853 std::shared_ptr<MultiDomainFunction> md_function = std::make_shared<MultiDomainFunction>();
1854
1855 // Set function first
1856 md_function->addFunction(std::move(fit_function));
1857
1858 // set domain for function with index 0 covering both sides
1859 md_function->clearDomainIndices();
1860 md_function->setDomainIndices(0, {0, 1});
1861
1862 // Set the properties
1863 fit->setProperty("Function", std::dynamic_pointer_cast<IFunction>(md_function));
1864 fit->setProperty("InputWorkspace", dataws);
1865 fit->setProperty("WorkspaceIndex", static_cast<int>(wsindex));
1866 fit->setProperty("StartX", vec_xmin.first);
1867 fit->setProperty("EndX", vec_xmax.first);
1868 fit->setProperty("InputWorkspace_1", dataws);
1869 fit->setProperty("WorkspaceIndex_1", static_cast<int>(wsindex));
1870 fit->setProperty("StartX_1", vec_xmin.second);
1871 fit->setProperty("EndX_1", vec_xmax.second);
1872 fit->setProperty("MaxIterations", m_fitIterations);
1873 fit->setProperty("IgnoreInvalidData", true);
1874
1875 // Execute
1876 fit->execute();
1877 if (!fit->isExecuted()) {
1878 throw runtime_error("Fit is not executed on multi-domain function/data. ");
1879 }
1880
1881 // Retrieve result
1882 std::string fitStatus = fit->getProperty("OutputStatus");
1883
1884 double chi2 = DBL_MAX;
1885 if (fitStatus == "success") {
1886 chi2 = fit->getProperty("OutputChi2overDoF");
1887 }
1888
1889 return chi2;
1890}
1891
1892//----------------------------------------------------------------------------------------------
1894double FitPeaks::fitFunctionHighBackground(const IAlgorithm_sptr &fit, const std::pair<double, double> &fit_window,
1895 const size_t &ws_index, const double &expected_peak_center,
1896 bool observe_peak_shape, const API::IPeakFunction_sptr &peakfunction,
1897 const API::IBackgroundFunction_sptr &bkgdfunc) {
1899
1900 // high background to reduce
1901 API::IBackgroundFunction_sptr high_bkgd_function =
1902 std::dynamic_pointer_cast<API::IBackgroundFunction>(m_linearBackgroundFunction->clone());
1903
1904 // Fit the background first if there is enough data points
1905 fitBackground(ws_index, fit_window, expected_peak_center, high_bkgd_function);
1906
1907 // Get partial of the data
1908 std::vector<double> vec_x, vec_y, vec_e;
1909 getRangeData(ws_index, fit_window, vec_x, vec_y, vec_e);
1910
1911 // Reduce the background
1912 reduceByBackground(high_bkgd_function, vec_x, vec_y);
1913 for (std::size_t n = 0; n < bkgdfunc->nParams(); ++n)
1914 bkgdfunc->setParameter(n, 0);
1915
1916 // Create a new workspace
1917 API::MatrixWorkspace_sptr reduced_bkgd_ws = createMatrixWorkspace(vec_x, vec_y, vec_e);
1918
1919 // Fit peak with background
1920 fitFunctionSD(fit, peakfunction, bkgdfunc, reduced_bkgd_ws, 0, {vec_x.front(), vec_x.back()}, expected_peak_center,
1921 observe_peak_shape, false);
1922
1923 // add the reduced background back
1924 bkgdfunc->setParameter(0, bkgdfunc->getParameter(0) + high_bkgd_function->getParameter(0));
1925 bkgdfunc->setParameter(1, bkgdfunc->getParameter(1) + // TODO doesn't work for flat background
1926 high_bkgd_function->getParameter(1));
1927
1928 double cost = fitFunctionSD(fit, peakfunction, bkgdfunc, m_inputMatrixWS, ws_index, {vec_x.front(), vec_x.back()},
1929 expected_peak_center, false, false);
1930
1931 return cost;
1932}
1933
1934//----------------------------------------------------------------------------------------------
1937 const std::vector<double> &vec_y,
1938 const std::vector<double> &vec_e) {
1939 std::size_t size = vec_x.size();
1940 std::size_t ysize = vec_y.size();
1941
1942 HistogramBuilder builder;
1943 builder.setX(size);
1944 builder.setY(ysize);
1945 MatrixWorkspace_sptr matrix_ws = create<Workspace2D>(1, builder.build());
1946
1947 auto &dataX = matrix_ws->mutableX(0);
1948 auto &dataY = matrix_ws->mutableY(0);
1949 auto &dataE = matrix_ws->mutableE(0);
1950
1951 dataX.assign(vec_x.cbegin(), vec_x.cend());
1952 dataY.assign(vec_y.cbegin(), vec_y.cend());
1953 dataE.assign(vec_e.cbegin(), vec_e.cend());
1954 return matrix_ws;
1955}
1956
1957//----------------------------------------------------------------------------------------------
1961 // create output workspace for peak positions: can be partial spectra to input
1962 // workspace
1963 m_outputPeakPositionWorkspace = create<Workspace2D>(m_numSpectraToFit, Points(m_numPeaksToFit));
1964 // set default
1965 for (std::size_t wi = 0; wi < m_numSpectraToFit; ++wi) {
1966 // convert to workspace index of input data workspace
1967 std::size_t inp_wi = wi + m_startWorkspaceIndex;
1968 std::vector<double> expected_position = m_getExpectedPeakPositions(inp_wi);
1969 for (std::size_t ipeak = 0; ipeak < expected_position.size(); ++ipeak) {
1970 m_outputPeakPositionWorkspace->dataX(wi)[ipeak] = expected_position[ipeak];
1971 }
1972 }
1973
1974 return;
1975}
1976
1977//----------------------------------------------------------------------------------------------
1985 const std::vector<std::string> &param_names, bool with_chi2) {
1986 // add columns
1987 table_ws->addColumn("int", "wsindex");
1988 table_ws->addColumn("int", "peakindex");
1989 for (const auto &param_name : param_names)
1990 table_ws->addColumn("double", param_name);
1991 if (with_chi2)
1992 table_ws->addColumn("double", "chi2");
1993
1994 // add rows
1995 const size_t numParam = m_fittedParamTable->columnCount() - 3;
1996 for (size_t iws = m_startWorkspaceIndex; iws <= m_stopWorkspaceIndex; ++iws) {
1997 for (size_t ipeak = 0; ipeak < m_numPeaksToFit; ++ipeak) {
1998 API::TableRow newRow = table_ws->appendRow();
1999 newRow << static_cast<int>(iws); // workspace index
2000 newRow << static_cast<int>(ipeak); // peak number
2001 for (size_t iparam = 0; iparam < numParam; ++iparam)
2002 newRow << 0.; // parameters for each peak
2003 if (with_chi2)
2004 newRow << DBL_MAX; // chisq
2005 }
2006 }
2007
2008 return;
2009}
2010
2011//----------------------------------------------------------------------------------------------
2017 // peak parameter workspace
2018 m_rawPeaksTable = getProperty(PropertyNames::RAW_PARAMS);
2019
2020 // create parameters
2021 // peak
2022 std::vector<std::string> param_vec;
2023 if (m_rawPeaksTable) {
2024 param_vec = m_peakFunction->getParameterNames();
2025 } else {
2026 param_vec.emplace_back("centre");
2027 param_vec.emplace_back("width");
2028 param_vec.emplace_back("height");
2029 param_vec.emplace_back("intensity");
2030 }
2031 // background
2032 for (size_t iparam = 0; iparam < m_bkgdFunction->nParams(); ++iparam)
2033 param_vec.emplace_back(m_bkgdFunction->parameterName(iparam));
2034
2035 // parameter value table
2036 m_fittedParamTable = std::make_shared<TableWorkspace>();
2038
2039 // for error workspace
2040 std::string fiterror_table_name = getPropertyValue(PropertyNames::OUTPUT_WKSP_PARAM_ERRS);
2041 // do nothing if user does not specifiy
2042 if (fiterror_table_name.empty()) {
2043 // not specified
2044 m_fitErrorTable = nullptr;
2045 } else {
2046 // create table and set up parameter table
2047 m_fitErrorTable = std::make_shared<TableWorkspace>();
2049 }
2050
2051 return;
2052}
2053
2054//----------------------------------------------------------------------------------------------
2059 // matrix workspace contained calculated peaks from fitting
2060 std::string fit_ws_name = getPropertyValue(PropertyNames::OUTPUT_WKSP_MODEL);
2061 if (fit_ws_name.size() == 0) {
2062 // skip if user does not specify
2063 m_fittedPeakWS = nullptr;
2064 return;
2065 }
2066
2067 // create a wokspace with same size as in the input matrix workspace
2068 m_fittedPeakWS = create<Workspace2D>(*m_inputMatrixWS);
2069}
2070
2071//----------------------------------------------------------------------------------------------
2073void FitPeaks::processOutputs(std::vector<std::shared_ptr<FitPeaksAlgorithm::PeakFitResult>> fit_result_vec) {
2075 setProperty(PropertyNames::OUTPUT_WKSP_PARAMS, m_fittedParamTable);
2076
2077 if (m_fitErrorTable) {
2078 g_log.warning("Output error table workspace");
2079 setProperty(PropertyNames::OUTPUT_WKSP_PARAM_ERRS, m_fitErrorTable);
2080 } else {
2081 g_log.warning("No error table output");
2082 }
2083
2084 // optional
2086 g_log.debug("about to calcualte fitted peaks");
2087 calculateFittedPeaks(std::move(fit_result_vec));
2088 setProperty(PropertyNames::OUTPUT_WKSP_MODEL, m_fittedPeakWS);
2089 }
2090}
2091
2092//----------------------------------------------------------------------------------------------
2097double FitPeaks::numberCounts(size_t iws) {
2098 const std::vector<double> &vec_y = m_inputMatrixWS->histogram(iws).y().rawData();
2099 double total = std::accumulate(vec_y.begin(), vec_y.end(), 0.);
2100 return total;
2101}
2102
2103//----------------------------------------------------------------------------------------------
2109double FitPeaks::numberCounts(size_t iws, const std::pair<double, double> &range) {
2110 // get data range
2111 std::vector<double> vec_x, vec_y, vec_e;
2112 getRangeData(iws, range, vec_x, vec_y, vec_e);
2113 // sum up all counts
2114 double total = std::accumulate(vec_y.begin(), vec_y.end(), 0.);
2115 return total;
2116}
2117
2118//----------------------------------------------------------------------------------------------
2124size_t FitPeaks::histRangeToDataPointCount(size_t iws, const std::pair<double, double> &range) {
2125 size_t left_index, right_index;
2126 histRangeToIndexBounds(iws, range, left_index, right_index);
2127 size_t number_dp = right_index - left_index + 1;
2128 if (m_inputMatrixWS->isHistogramData())
2129 number_dp -= 1;
2130 assert(number_dp > 0);
2131 return number_dp;
2132}
2133
2134GNU_DIAG_OFF("dangling-reference")
2135
2136//----------------------------------------------------------------------------------------------
2143void FitPeaks::histRangeToIndexBounds(size_t iws, const std::pair<double, double> &range, size_t &left_index,
2144 size_t &right_index) {
2145 const auto &orig_x = m_inputMatrixWS->histogram(iws).x();
2146 rangeToIndexBounds(orig_x, range.first, range.second, left_index, right_index);
2147
2148 // handle an invalid range case. For the histogram point data, make sure the number of data points is non-zero as
2149 // well.
2150 if (left_index >= right_index || (m_inputMatrixWS->isHistogramData() && left_index == right_index - 1)) {
2151 std::stringstream err_ss;
2152 err_ss << "Unable to get a valid subset of histogram from given fit window. "
2153 << "Histogram X: " << orig_x.front() << "," << orig_x.back() << "; Range: " << range.first << ","
2154 << range.second;
2155 throw std::runtime_error(err_ss.str());
2156 }
2157}
2158
2159//----------------------------------------------------------------------------------------------
2167void FitPeaks::getRangeData(size_t iws, const std::pair<double, double> &range, std::vector<double> &vec_x,
2168 std::vector<double> &vec_y, std::vector<double> &vec_e) {
2169 // convert range to index boundaries
2170 size_t left_index, right_index;
2171 histRangeToIndexBounds(iws, range, left_index, right_index);
2172
2173 // copy X, Y and E
2174 size_t num_elements_x = right_index - left_index;
2175
2176 vec_x.resize(num_elements_x);
2177 const auto &orig_x = m_inputMatrixWS->histogram(iws).x();
2178 std::copy(orig_x.begin() + left_index, orig_x.begin() + right_index, vec_x.begin());
2179
2180 size_t num_datapoints = m_inputMatrixWS->isHistogramData() ? num_elements_x - 1 : num_elements_x;
2181
2182 const std::vector<double> orig_y = m_inputMatrixWS->histogram(iws).y().rawData();
2183 const std::vector<double> orig_e = m_inputMatrixWS->histogram(iws).e().rawData();
2184 vec_y.resize(num_datapoints);
2185 vec_e.resize(num_datapoints);
2186 std::copy(orig_y.begin() + left_index, orig_y.begin() + left_index + num_datapoints, vec_y.begin());
2187 std::copy(orig_e.begin() + left_index, orig_e.begin() + left_index + num_datapoints, vec_e.begin());
2188}
2189
2190GNU_DIAG_ON("dangling-reference")
2191
2192//----------------------------------------------------------------------------------------------
2199double FitPeaks::calculateSignalToNoiseRatio(size_t iws, const std::pair<double, double> &range,
2200 const API::IBackgroundFunction_sptr &bkgd_function) {
2201 // convert range to index boundaries
2202 size_t left_index, right_index;
2203 histRangeToIndexBounds(iws, range, left_index, right_index);
2204
2205 // estimate background level by Y(X) fitting
2206 if (!estimateBackgroundParameters(m_inputMatrixWS->histogram(iws), std::pair<size_t, size_t>(left_index, right_index),
2207 bkgd_function))
2208 return 0.0; // failed to estimate background parameters
2209
2210 // get X,Y,and E for the data range
2211 std::vector<double> vec_x, vec_y, vec_e;
2212 getRangeData(iws, range, vec_x, vec_y, vec_e);
2213 if (vec_x.empty())
2214 return 0.0;
2215
2216 // subtract background from Y-values
2217 reduceByBackground(bkgd_function, vec_x, vec_y);
2218
2219 // estimate the signal as the highest Y-value in the data range
2220 auto it_max = std::max_element(vec_y.begin(), vec_y.end());
2221 double signal = vec_y[it_max - vec_y.begin()];
2222 if (signal <= DBL_MIN)
2223 return 0.0;
2224
2225 // estimate noise from background. If noise is zero, or impossible to estimate, return DBL_MAX so that the peak
2226 // won't be rejected.
2227 double noise = estimateBackgroundNoise(vec_y);
2228 if (noise <= DBL_MIN)
2229 return DBL_MAX;
2230
2231 // finally, calculate the signal-to-noise ratio
2232 return signal / noise;
2233}
2234
2235//----------------------------------------------------------------------------------------------
2237
2238void FitPeaks::checkWorkspaceIndices(std::size_t const &wi) {
2239 if (wi < m_startWorkspaceIndex || wi > m_stopWorkspaceIndex) {
2240 std::stringstream errss;
2241 errss << "Workspace index " << wi << " is out of range "
2242 << "[" << m_startWorkspaceIndex << ", " << m_stopWorkspaceIndex << "]";
2243 throw std::runtime_error(errss.str());
2244 }
2245}
2246
2247void FitPeaks::checkPeakIndices(std::size_t const &wi, std::size_t const &ipeak) {
2248 // check peak index
2249 if (ipeak >= m_getExpectedPeakPositions(wi).size()) {
2250 std::stringstream errss;
2251 errss << "Peak index " << ipeak << " is out of range (" << m_numPeaksToFit << ")";
2252 throw std::runtime_error(errss.str());
2253 }
2254}
2255
2256void FitPeaks::checkPeakWindowEdgeOrder(double const &left, double const &right) {
2257 if (left >= right) {
2258 std::stringstream errss;
2259 errss << "Peak window is inappropriate for workspace index: " << left << " >= " << right;
2260 throw std::runtime_error(errss.str());
2261 }
2262}
2263
2264//----------------------------------------------------------------------------------------------
2273void FitPeaks::writeFitResult(size_t wi, const std::vector<double> &expected_positions,
2274 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result) {
2275 // convert to
2276 size_t out_wi = wi - m_startWorkspaceIndex;
2277 if (out_wi >= m_outputPeakPositionWorkspace->getNumberHistograms()) {
2278 g_log.error() << "workspace index " << wi << " is out of output peak position workspace "
2279 << "range of spectra, which contains " << m_outputPeakPositionWorkspace->getNumberHistograms()
2280 << " spectra"
2281 << "\n";
2282 throw std::runtime_error("Out of boundary to set output peak position workspace");
2283 }
2284
2285 // Fill the output peak position workspace
2286 for (size_t ipeak = 0; ipeak < m_numPeaksToFit; ++ipeak) {
2287 double exp_peak_pos(expected_positions[ipeak]);
2288 double fitted_peak_pos = fit_result->getPeakPosition(ipeak);
2289 double peak_chi2 = fit_result->getCost(ipeak);
2290
2291 m_outputPeakPositionWorkspace->mutableX(out_wi)[ipeak] = exp_peak_pos;
2292 m_outputPeakPositionWorkspace->mutableY(out_wi)[ipeak] = fitted_peak_pos;
2293 m_outputPeakPositionWorkspace->mutableE(out_wi)[ipeak] = peak_chi2;
2294 }
2295
2296 // Output the peak parameters to the table workspace
2297 // check vector size
2298
2299 // last column of the table is for chi2
2300 size_t chi2_index = m_fittedParamTable->columnCount() - 1;
2301
2302 // check TableWorkspace and given FitResult
2303 if (m_rawPeaksTable) {
2304 // duplicate from FitPeakResult to table workspace
2305 // check again with the column size versus peak parameter values
2306 if (fit_result->getNumberParameters() != m_fittedParamTable->columnCount() - 3) {
2307 g_log.error() << "Peak of type (" << m_peakFunction->name() << ") has " << fit_result->getNumberParameters()
2308 << " parameters. Parameter table shall have 3 more "
2309 "columns. But not it has "
2310 << m_fittedParamTable->columnCount() << " columns\n";
2311 throw std::runtime_error("Peak parameter vector for one peak has different sizes to output "
2312 "table workspace");
2313 }
2314 } else {
2315 // effective peak profile parameters: need to re-construct the peak function
2316 if (4 + m_bkgdFunction->nParams() != m_fittedParamTable->columnCount() - 3) {
2317
2318 std::stringstream err_ss;
2319 err_ss << "Peak has 4 effective peak parameters and " << m_bkgdFunction->nParams() << " background parameters "
2320 << ". Parameter table shall have 3 more columns. But not it has " << m_fittedParamTable->columnCount()
2321 << " columns";
2322 throw std::runtime_error(err_ss.str());
2323 }
2324 }
2325
2326 // go through each peak
2327 // get a copy of peak function and background function
2328 IPeakFunction_sptr peak_function = std::dynamic_pointer_cast<IPeakFunction>(m_peakFunction->clone());
2329 size_t num_peakfunc_params = peak_function->nParams();
2330 size_t num_bkgd_params = m_bkgdFunction->nParams();
2331
2332 for (size_t ipeak = 0; ipeak < m_numPeaksToFit; ++ipeak) {
2333 // get row number
2334 size_t row_index = out_wi * m_numPeaksToFit + ipeak;
2335
2336 // treat as different cases for writing out raw or effective parametr
2337 if (m_rawPeaksTable) {
2338 // duplicate from FitPeakResult to table workspace
2339 for (size_t iparam = 0; iparam < num_peakfunc_params + num_bkgd_params; ++iparam) {
2340 size_t col_index = iparam + 2;
2341 // fitted parameter's value
2342 m_fittedParamTable->cell<double>(row_index, col_index) = fit_result->getParameterValue(ipeak, iparam);
2343 // fitted parameter's fitting error
2344 if (m_fitErrorTable) {
2345 m_fitErrorTable->cell<double>(row_index, col_index) = fit_result->getParameterError(ipeak, iparam);
2346 }
2347
2348 } // end for (iparam)
2349 } else {
2350 // effective peak profile parameter
2351 // construct the peak function
2352 for (size_t iparam = 0; iparam < num_peakfunc_params; ++iparam)
2353 peak_function->setParameter(iparam, fit_result->getParameterValue(ipeak, iparam));
2354
2355 // set the effective peak parameters
2356 m_fittedParamTable->cell<double>(row_index, 2) = peak_function->centre();
2357 m_fittedParamTable->cell<double>(row_index, 3) = peak_function->fwhm();
2358 m_fittedParamTable->cell<double>(row_index, 4) = peak_function->height();
2359 m_fittedParamTable->cell<double>(row_index, 5) = peak_function->intensity();
2360
2361 // background
2362 for (size_t iparam = 0; iparam < num_bkgd_params; ++iparam)
2363 m_fittedParamTable->cell<double>(row_index, 6 + iparam) =
2364 fit_result->getParameterValue(ipeak, num_peakfunc_params + iparam);
2365 }
2366
2367 // set chi2
2368 m_fittedParamTable->cell<double>(row_index, chi2_index) = fit_result->getCost(ipeak);
2369 }
2370
2371 return;
2372}
2373
2374//----------------------------------------------------------------------------------------------
2376 std::string height_name("");
2377
2378 std::vector<std::string> peak_parameters = peak_function->getParameterNames();
2379 for (const auto &parName : peak_parameters) {
2380 if (parName == "Height") {
2381 height_name = "Height";
2382 break;
2383 } else if (parName == "I") {
2384 height_name = "I";
2385 break;
2386 } else if (parName == "Intensity") {
2387 height_name = "Intensity";
2388 break;
2389 }
2390 }
2391
2392 if (height_name.empty())
2393 throw std::runtime_error("Peak height parameter name cannot be found.");
2394
2395 return height_name;
2396}
2397
2398// A client, like PDCalibration, may set a logging offset to make FitPeaks less "chatty".
2399// This method temporarily removes the logging offset and logs the message at its priority level.
2400void FitPeaks::logNoOffset(const size_t &priority, const std::string &msg) {
2401 LoggingOffsetSentry sentry(this);
2402
2403 switch (priority) {
2404 case 4: // warning
2405 g_log.warning() << msg;
2406 break;
2407 case 5: // notice
2408 g_log.notice() << msg;
2409 break;
2410 default:
2411 assert(false); // not implemented yet
2412 }
2413}
2414
2416
2417} // namespace Mantid::Algorithms
#define DECLARE_ALGORITHM(classname)
Definition Algorithm.h:538
double value
The value of the point.
Definition FitMW.cpp:51
Algorithm *const m_alg
int m_loggingOffset
double sigma
Definition GetAllEi.cpp:156
double left
double right
#define fabs(x)
Definition Matrix.cpp:22
#define PARALLEL_START_INTERRUPT_REGION
Begins a block to skip processing is the algorithm has been interupted Note the end of the block if n...
#define PARALLEL_CRITICAL(name)
#define PARALLEL_END_INTERRUPT_REGION
Ends a block to skip processing is the algorithm has been interupted Note the start of the block if n...
#define PARALLEL_FOR_IF(condition)
Empty definitions - to enable set your complier to enable openMP.
#define PRAGMA_OMP(expression)
#define PARALLEL_CHECK_INTERRUPT_REGION
Adds a check after a Parallel region to see if it was interupted.
#define GNU_DIAG_ON(x)
#define GNU_DIAG_OFF(x)
This is a collection of macros for turning compiler warnings off in a controlled manner.
Base class from which all concrete algorithm classes should be derived.
Definition Algorithm.h:76
void declareProperty(std::unique_ptr< Kernel::Property > p, const std::string &doc="") override
Add a property to the list of managed properties.
std::string getPropertyValue(const std::string &name) const override
Get the value of a property as a string.
TypedValue getProperty(const std::string &name) const override
Get the value of a property.
virtual std::shared_ptr< Algorithm > createChildAlgorithm(const std::string &name, const double startProgress=-1., const double endProgress=-1., const bool enableLogging=true, const int &version=-1)
Create a Child Algorithm.
void setLoggingOffset(const int value) override
gets the logging priority offset
Kernel::Logger & g_log
Definition Algorithm.h:422
int getLoggingOffset() const override
returns the logging priority offset
bool isDefault(const std::string &name) const
static bool isEmpty(const NumT toCheck)
checks that the value was not set by users, uses the value in empty double/int.
Implements FunctionDomain1D with its own storage in form of a std::vector.
A class to store values calculated by a function.
const std::vector< double > & toVector() const
Return the calculated values as a vector.
double getCalculated(size_t i) const
Get i-th calculated value.
An interface to a peak function, which extend the interface of IFunctionWithLocation by adding method...
Helper class for reporting progress from algorithms.
Definition Progress.h:25
TableRow represents a row in a TableWorkspace.
Definition TableRow.h:39
A property class for workspaces.
PeakFitPreCheckResult & operator+=(const PeakFitPreCheckResult &another)
Definition FitPeaks.cpp:200
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
PeakFitResult(size_t num_peaks, size_t num_params)
Holds all of the fitting information for a single spectrum.
Definition FitPeaks.cpp:96
double getParameterValue(size_t ipeak, size_t iparam) const
get the fitted value of a particular parameter
Definition FitPeaks.cpp:137
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:149
double getParameterError(size_t ipeak, size_t iparam) const
get the fitting error of a particular parameter
Definition FitPeaks.cpp:126
void setBadRecord(size_t ipeak, const double peak_position)
The peak postition should be negative and indicates what went wrong.
Definition FitPeaks.cpp:180
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:312
void calculateFittedPeaks(const std::vector< std::shared_ptr< FitPeaksAlgorithm::PeakFitResult > > &fit_results)
calculate peak+background for fitted
double calculateSignalToNoiseRatio(size_t iws, const std::pair< double, double > &range, const API::IBackgroundFunction_sptr &bkgd_function)
calculate signal-to-noise ratio in histogram range
API::MatrixWorkspace_const_sptr m_peakCenterWorkspace
Definition FitPeaks.h:286
void generateFittedParametersValueWorkspaces()
Generate output workspaces.
void setupParameterTableWorkspace(const API::ITableWorkspace_sptr &table_ws, const std::vector< std::string > &param_names, bool with_chi2)
Set up parameter table (parameter value or error)
std::vector< double > m_peakPosTolerances
tolerances for fitting peak positions
Definition FitPeaks.h:308
API::IPeakFunction_sptr m_peakFunction
Peak profile name.
Definition FitPeaks.h:263
bool fitBackground(const size_t &ws_index, const std::pair< double, double > &fit_window, const double &expected_peak_pos, const API::IBackgroundFunction_sptr &bkgd_func)
fit background
double fitFunctionSD(const API::IAlgorithm_sptr &fit, const API::IPeakFunction_sptr &peak_function, const API::IBackgroundFunction_sptr &bkgd_function, const API::MatrixWorkspace_sptr &dataws, size_t wsindex, const std::pair< double, double > &peak_range, const double &expected_peak_center, bool estimate_peak_width, bool estimate_background)
Methods to fit functions (general)
API::MatrixWorkspace_sptr m_outputPeakPositionWorkspace
output workspace for peak positions
Definition FitPeaks.h:247
API::MatrixWorkspace_sptr m_fittedPeakWS
matrix workspace contained calcalated peaks+background from fitted result it has same number of spect...
Definition FitPeaks.h:259
API::ITableWorkspace_const_sptr m_profileStartingValueTable
table workspace for profile parameters' starting value
Definition FitPeaks.h:325
std::string m_minimizer
Minimzer.
Definition FitPeaks.h:270
API::IBackgroundFunction_sptr m_linearBackgroundFunction
Linear background function for high background fitting.
Definition FitPeaks.h:267
bool m_peakPosTolCase234
peak positon tolerance case b, c and d
Definition FitPeaks.h:345
void fitSpectrumPeaks(size_t wi, const std::vector< double > &expected_peak_centers, const std::shared_ptr< FitPeaksAlgorithm::PeakFitResult > &fit_result, std::vector< std::vector< double > > &lastGoodPeakParameters, std::vector< size_t > &lastGoodPeakSpectra, const std::shared_ptr< FitPeaksAlgorithm::PeakFitPreCheckResult > &pre_check_result)
fit peaks in a same spectrum
std::map< std::string, std::string > validateInputs() override
Validate inputs.
Definition FitPeaks.cpp:456
double m_peakWidthPercentage
flag to estimate peak width from
Definition FitPeaks.h:298
API::IBackgroundFunction_sptr m_bkgdFunction
Background function.
Definition FitPeaks.h:265
size_t histRangeToDataPointCount(size_t iws, const std::pair< double, double > &range)
convert a histogram range to index boundaries
void checkPeakIndices(std::size_t const &, std::size_t const &)
double fitFunctionHighBackground(const API::IAlgorithm_sptr &fit, const std::pair< double, double > &fit_window, const size_t &ws_index, const double &expected_peak_center, bool observe_peak_shape, const API::IPeakFunction_sptr &peakfunction, const API::IBackgroundFunction_sptr &bkgdfunc)
fit a single peak with high background
void processInputPeakCenters()
peak centers
Definition FitPeaks.cpp:833
std::vector< std::shared_ptr< FitPeaksAlgorithm::PeakFitResult > > fitPeaks()
suites of method to fit peaks
Definition FitPeaks.cpp:984
double m_minPeakHeight
minimum peak height without background and it also serves as the criteria for observed peak parameter
Definition FitPeaks.h:332
std::string m_costFunction
Cost function.
Definition FitPeaks.h:272
void processInputFunctions()
process inputs for peak and background functions
Definition FitPeaks.cpp:641
void processInputPeakTolerance()
process inputs about fitted peak positions' tolerance
Definition FitPeaks.cpp:894
void writeFitResult(size_t wi, const std::vector< double > &expected_positions, const std::shared_ptr< FitPeaksAlgorithm::PeakFitResult > &fit_result)
Write result of peak fit per spectrum to output analysis workspaces.
void histRangeToIndexBounds(size_t iws, const std::pair< double, double > &range, size_t &left_index, size_t &right_index)
convert a histogram range to index boundaries
void exec() override
Main exec method.
Definition FitPeaks.cpp:540
void init() override
Init.
Definition FitPeaks.cpp:270
void logNoOffset(const size_t &priority, const std::string &msg)
std::size_t m_numSpectraToFit
total number of spectra to be fit
Definition FitPeaks.h:306
std::vector< std::string > m_peakParamNames
input peak parameters' names
Definition FitPeaks.h:320
double calculateSignalToSigmaRatio(const size_t &iws, const std::pair< double, double > &peakWindow, const API::IPeakFunction_sptr &peakFunction)
bool isObservablePeakProfile(const std::string &peakprofile)
check whether FitPeaks supports observation on a certain peak profile's parameters (width!...
void processInputFitRanges()
process inputs for peak fitting range
Definition FitPeaks.cpp:694
std::size_t m_numPeaksToFit
the number of peaks to fit in all spectra
Definition FitPeaks.h:288
std::size_t m_startWorkspaceIndex
start index
Definition FitPeaks.h:302
API::MatrixWorkspace_sptr createMatrixWorkspace(const std::vector< double > &vec_x, const std::vector< double > &vec_y, const std::vector< double > &vec_e)
Create a single spectrum workspace for fitting.
void generateOutputPeakPositionWS()
main method to create output workspaces
void generateCalculatedPeaksWS()
Generate workspace for calculated values.
double fitFunctionMD(API::IFunction_sptr fit_function, const API::MatrixWorkspace_sptr &dataws, const size_t wsindex, const std::pair< double, double > &vec_xmin, const std::pair< double, double > &vec_xmax)
std::string getPeakHeightParameterName(const API::IPeakFunction_const_sptr &peak_function)
Get the parameter name for peak height (I or height or etc)
API::ITableWorkspace_sptr m_fittedParamTable
output analysis workspaces table workspace for fitted parameters
Definition FitPeaks.h:250
bool decideToEstimatePeakParams(const bool firstPeakInSpectrum, const API::IPeakFunction_sptr &peak_function)
Decide whether to estimate peak parameters.
void convertParametersNameToIndex()
Convert peak function's parameter names to parameter index for fast access.
Definition FitPeaks.cpp:952
void processOutputs(std::vector< std::shared_ptr< FitPeaksAlgorithm::PeakFitResult > > fit_result_vec)
Set the workspaces and etc to output properties.
bool m_highBackground
flag for high background
Definition FitPeaks.h:341
std::vector< double > m_initParamValues
input peak parameters' starting values corresponding to above peak parameter names
Definition FitPeaks.h:323
std::vector< double > m_peakCenters
Designed peak positions and tolerance.
Definition FitPeaks.h:285
double fitIndividualPeak(size_t wi, const API::IAlgorithm_sptr &fitter, const double expected_peak_center, const std::pair< double, double > &fitwindow, const bool estimate_peak_width, const API::IPeakFunction_sptr &peakfunction, const API::IBackgroundFunction_sptr &bkgdfunc, const std::shared_ptr< FitPeaksAlgorithm::PeakFitPreCheckResult > &pre_check_result)
Fit an individual peak.
std::vector< std::vector< double > > m_peakWindowVector
peak windows
Definition FitPeaks.h:316
std::function< std::pair< double, double >(std::size_t const &, std::size_t const &)> m_getPeakFitWindow
Definition FitPeaks.h:292
API::MatrixWorkspace_sptr m_inputMatrixWS
mandatory input and output workspaces
Definition FitPeaks.h:242
std::vector< size_t > m_initParamIndexes
input starting parameters' indexes in peak function
Definition FitPeaks.h:282
bool m_fitPeaksFromRight
Fit from right or left.
Definition FitPeaks.h:274
std::function< std::vector< double >(std::size_t const &)> m_getExpectedPeakPositions
Definition FitPeaks.h:291
bool m_rawPeaksTable
flag to show that the pamarameters in table are raw parameters or effective parameters
Definition FitPeaks.h:255
void processInputs()
process inputs (main and child algorithms)
Definition FitPeaks.cpp:561
void checkWorkspaceIndices(std::size_t const &)
Get the expected peak's position.
void checkPeakWindowEdgeOrder(double const &, double const &)
int m_fitIterations
Fit iterations.
Definition FitPeaks.h:276
double numberCounts(size_t iws)
sum up all counts in histogram
API::MatrixWorkspace_const_sptr m_peakWindowWorkspace
Definition FitPeaks.h:317
bool m_uniformProfileStartingValue
flag for profile startng value being uniform or not
Definition FitPeaks.h:327
API::ITableWorkspace_sptr m_fitErrorTable
table workspace for fitted parameters' fitting error. This is optional
Definition FitPeaks.h:252
void getRangeData(size_t iws, const std::pair< double, double > &range, std::vector< double > &vec_x, std::vector< double > &vec_y, std::vector< double > &vec_e)
get vector X, Y and E in a given range
std::size_t m_stopWorkspaceIndex
stop index (workspace index of the last spectrum included)
Definition FitPeaks.h:304
bool processSinglePeakFitResult(size_t wsindex, size_t peakindex, const double cost, const std::vector< double > &expected_peak_positions, const FitPeaksAlgorithm::FitFunction &fitfunction, const std::shared_ptr< FitPeaksAlgorithm::PeakFitResult > &fit_result)
Process the result from fitting a single peak.
Support for a property that holds an array of values.
Exception for when an item is not found in a collection.
Definition Exception.h:145
IPropertyManager * setProperty(const std::string &name, const T &value)
Templated method to set the value of a PropertyWithValue.
void setPropertyGroup(const std::string &name, const std::string &group)
Set the group for a given property.
ListValidator is a validator that requires the value of a property to be one of a defined list of pos...
void debug(const std::string &msg)
Logs at debug level.
Definition Logger.cpp:145
void notice(const std::string &msg)
Logs at notice level.
Definition Logger.cpp:126
void error(const std::string &msg)
Logs at error level.
Definition Logger.cpp:108
void warning(const std::string &msg)
Logs at warning level.
Definition Logger.cpp:117
void report()
Increments the loop counter by 1, then sends the progress notification on behalf of its algorithm.
StartsWithValidator is a validator that requires the value of a property to start with one of the str...
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 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< CompositeFunction > CompositeFunction_sptr
shared pointer to the composite function base class
std::string const OUTPUT_WKSP("OutputWorkspace")
std::string const INPUT_WKSP("InputWorkspace")
MANTID_ALGORITHMS_DLL size_t findXIndex(const vector_like &vecx, const double x, const size_t startindex=0)
Get an index of a value in a sorted vector.
MANTID_ALGORITHMS_DLL int estimatePeakParameters(const HistogramData::Histogram &histogram, const std::pair< size_t, size_t > &peak_window, const API::IPeakFunction_sptr &peakfunction, const API::IBackgroundFunction_sptr &bkgdfunction, bool observe_peak_width, const EstimatePeakWidth peakWidthEstimateApproach, const double peakWidthPercentage, const double minPeakHeight)
Estimate peak parameters by 'observation'.
Statistics getStatistics(const std::vector< TYPE > &data, const unsigned int flags=StatOptions::AllStats)
Return a statistics object for the given data set.
std::vector< double > getZscore(const std::vector< TYPE > &data)
Return the Z score values for a dataset.
std::shared_ptr< IValidator > IValidator_sptr
A shared_ptr to an IValidator.
Definition IValidator.h:26
std::enable_if< std::is_pointer< Arg >::value, bool >::type threadSafe(Arg workspace)
Thread-safety check Checks the workspace to ensure it is suitable for multithreaded access.
const std::string OUTPUT_WKSP("OutputWorkspace")
const std::string INPUT_WKSP("InputWorkspace")
Helper class which provides the Collimation Length for SANS instruments.
constexpr int EMPTY_INT() noexcept
Returns what we consider an "empty" integer within a property.
Definition EmptyValues.h:24
constexpr double EMPTY_DBL() noexcept
Returns what we consider an "empty" double within a property.
Definition EmptyValues.h:42
STL namespace.
API::IBackgroundFunction_sptr bkgdfunction
Definition FitPeaks.h:35
@ Input
An input workspace.
Definition Property.h:53
@ Output
An output workspace.
Definition Property.h:54
Functor to accumulate a sum of squares.