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