27#include "MantidHistogramData/EstimatePolynomial.h"
28#include "MantidHistogramData/Histogram.h"
29#include "MantidHistogramData/HistogramBuilder.h"
30#include "MantidHistogramData/HistogramIterator.h"
39#include "boost/algorithm/string.hpp"
40#include "boost/algorithm/string/trim.hpp"
45using namespace Algorithms::PeakParameterHelper;
51using Mantid::HistogramData::Histogram;
60const std::string START_WKSP_INDEX(
"StartWorkspaceIndex");
61const std::string STOP_WKSP_INDEX(
"StopWorkspaceIndex");
62const std::string PEAK_CENTERS(
"PeakCenters");
63const std::string PEAK_CENTERS_WKSP(
"PeakCentersWorkspace");
64const std::string PEAK_FUNC(
"PeakFunction");
65const std::string BACK_FUNC(
"BackgroundType");
66const std::string FIT_WINDOW_LIST(
"FitWindowBoundaryList");
67const std::string FIT_WINDOW_WKSP(
"FitPeakWindowWorkspace");
68const std::string PEAK_WIDTH_PERCENT(
"PeakWidthPercent");
69const std::string PEAK_PARAM_NAMES(
"PeakParameterNames");
70const std::string PEAK_PARAM_VALUES(
"PeakParameterValues");
71const std::string PEAK_PARAM_TABLE(
"PeakParameterValueTable");
72const std::string FIT_FROM_RIGHT(
"FitFromRight");
73const std::string MINIMIZER(
"Minimizer");
74const std::string COST_FUNC(
"CostFunction");
75const std::string MAX_FIT_ITER(
"MaxFitIterations");
76const std::string BACKGROUND_Z_SCORE(
"FindBackgroundSigma");
77const std::string HIGH_BACKGROUND(
"HighBackground");
78const std::string POSITION_TOL(
"PositionTolerance");
79const std::string PEAK_MIN_HEIGHT(
"MinimumPeakHeight");
80const std::string CONSTRAIN_PEAK_POS(
"ConstrainPeakPositions");
81const std::string COPY_LAST_GOOD_PEAK_PARAMS(
"CopyLastGoodPeakParameters");
82const std::string OUTPUT_WKSP_MODEL(
"FittedPeaksWorkspace");
83const std::string OUTPUT_WKSP_PARAMS(
"OutputPeakParametersWorkspace");
84const std::string OUTPUT_WKSP_PARAM_ERRS(
"OutputParameterFitErrorsWorkspace");
85const std::string RAW_PARAMS(
"RawPeakParameters");
86const std::string PEAK_MIN_SIGNAL_TO_NOISE_RATIO(
"MinimumSignalToNoiseRatio");
87const std::string PEAK_MIN_TOTAL_COUNT(
"MinimumPeakTotalCount");
88const std::string PEAK_MIN_SIGNAL_TO_SIGMA_RATIO(
"MinimumSignalToSigmaRatio");
92namespace FitPeaksAlgorithm {
98 if (num_peaks == 0 || num_params == 0)
99 throw std::runtime_error(
"No peak or no parameter error.");
103 m_costs.resize(num_peaks, DBL_MAX);
106 for (
size_t ipeak = 0; ipeak < num_peaks; ++ipeak) {
153 throw std::runtime_error(
"Peak index is out of range.");
162 size_t peak_num_params = fit_functions.
peakfunction->nParams();
163 for (
size_t ipar = 0; ipar < peak_num_params; ++ipar) {
168 for (
size_t ipar = 0; ipar < fit_functions.
bkgdfunction->nParams(); ++ipar) {
183 throw std::runtime_error(
"Peak index is out of range");
184 if (peak_position >= 0.)
185 throw std::runtime_error(
"Can only set negative postion for bad record");
233 assert(individual_rejection_count <= 1);
235 return individual_rejection_count == 1;
245 std::ostringstream os;
256 os <<
m_low_snr <<
" peak(s) rejected: low signal-to-noise ratio.\n";
264 : m_fitPeaksFromRight(true), m_fitIterations(50), m_numPeaksToFit(0), m_minPeakHeight(0.),
265 m_minSignalToNoiseRatio(0.), m_minPeakTotalCount(0.), m_peakPosTolCase234(false) {}
272 "Name of the input workspace for peak fitting.");
275 "Name of the output workspace containing peak centers for "
277 "The output workspace is point data."
278 "Each workspace index corresponds to a spectrum. "
279 "Each X value ranges from 0 to N-1, where N is the number of "
281 "Each Y value is the peak position obtained by peak fitting. "
282 "Negative value is used for error signals. "
283 "-1 for data is zero; -2 for maximum value is smaller than "
284 "specified minimum value."
285 "and -3 for non-converged fitting.");
288 auto mustBePositive = std::make_shared<BoundedValidator<int>>();
289 mustBePositive->setLower(0);
290 declareProperty(PropertyNames::START_WKSP_INDEX, 0, mustBePositive,
"Starting workspace index for fit");
292 PropertyNames::STOP_WKSP_INDEX,
EMPTY_INT(),
293 "Last workspace index for fit is the smaller of this value and the workspace index of last spectrum.");
296 "List of peak centers to use as initial guess for fit.");
300 "MatrixWorkspace containing referent peak centers for each spectrum, defined at the same workspace indices.");
302 const std::string peakcentergrp(
"Peak Positions");
307 const std::vector<std::string> peakNames = FunctionFactory::Instance().getFunctionNames<
API::IPeakFunction>();
308 declareProperty(PropertyNames::PEAK_FUNC,
"Gaussian", std::make_shared<StringListValidator>(peakNames),
309 "Use of a BackToBackExponential profile is only reccomended if the "
310 "coeficients to calculate A and B are defined in the instrument "
311 "Parameters.xml file.");
312 const vector<string> bkgdtypes{
"Flat",
"Linear",
"Quadratic"};
313 declareProperty(PropertyNames::BACK_FUNC,
"Linear", std::make_shared<StringListValidator>(bkgdtypes),
314 "Type of Background.");
316 const std::string funcgroup(
"Function Types");
323 "List of boundaries of the peak fitting window corresponding to "
328 "MatrixWorkspace containing peak windows for each peak center in each spectrum, defined at the same "
329 "workspace indices.");
331 auto min = std::make_shared<BoundedValidator<double>>();
335 "The estimated peak width as a "
336 "percentage of the d-spacing "
337 "of the center of the peak. Value must be less than 1.");
339 const std::string fitrangeegrp(
"Peak Range Setup");
346 "List of peak parameters' names");
348 "List of peak parameters' value");
351 "Name of the an optional workspace, whose each column "
352 "corresponds to given peak parameter names, "
353 "and each row corresponds to a subset of spectra.");
355 const std::string startvaluegrp(
"Starting Parameters Setup");
362 "Flag for the order to fit peaks. If true, peaks are fitted "
364 "Otherwise peaks are fitted from leftmost.");
366 const std::vector<std::string> minimizerOptions = API::FuncMinimizerFactory::Instance().getKeys();
369 "Minimizer to use for fitting.");
371 const std::array<string, 2> costFuncOptions = {{
"Least squares",
"Rwp"}};
375 auto min_max_iter = std::make_shared<BoundedValidator<int>>();
376 min_max_iter->setLower(49);
377 declareProperty(PropertyNames::MAX_FIT_ITER, 50, min_max_iter,
"Maximum number of function fitting iterations.");
379 const std::string optimizergrp(
"Optimization Setup");
384 std::ostringstream os;
385 os <<
"Deprecated property. Use " << PropertyNames::PEAK_MIN_SIGNAL_TO_NOISE_RATIO <<
" instead.";
389 "Flag whether the input data has high background compared to peak heights.");
392 "List of tolerance on fitted peak positions against given peak positions."
393 "If there is only one value given, then ");
396 "Used for validating peaks before and after fitting. If a peak's observed/estimated or "
397 "fitted height is under this value, the peak will be marked as error.");
400 "If true peak position will be constrained by estimated positions "
401 "(highest Y value position) and "
402 "the peak width either estimted by observation or calculate.");
405 "If true, initial peak parameters (with the exception of peak centre) "
406 "may be copied from the last successfully fit peak in the spectra.");
411 "Name of the output matrix workspace with fitted peak. "
412 "This output workspace has the same dimension as the input workspace."
413 "The Y values belonged to peaks to fit are replaced by fitted value. "
414 "Values of estimated background are used if peak fails to be fit.");
418 "Name of table workspace containing all fitted peak parameters.");
424 "Name of workspace containing all fitted peak parameters' fitting error."
425 "It must be used along with FittedPeaksWorkspace and RawPeakParameters "
429 "false generates table with effective centre/width/height "
430 "parameters. true generates a table with peak function "
434 PropertyNames::PEAK_MIN_SIGNAL_TO_NOISE_RATIO, 0.,
435 "Used for validating peaks before fitting. If the signal-to-noise ratio is under this value, "
436 "the peak will be marked as error. This does not apply to peaks for which the noise cannot be estimated.");
439 "Used for validating peaks before fitting. If the total peak window Y-value count "
440 "is under this value, the peak will be excluded from fitting and calibration.");
443 "Used for validating peaks after fitting. If the signal-to-sigma ratio is under this value, "
444 "the peak will be excluded from fitting and calibration.");
446 const std::string addoutgrp(
"Analysis");
457 map<std::string, std::string> issues;
460 if (!(
isDefault(PropertyNames::START_WKSP_INDEX) &&
isDefault(PropertyNames::STOP_WKSP_INDEX))) {
461 const int startIndex =
getProperty(PropertyNames::START_WKSP_INDEX);
462 const int stopIndex =
getProperty(PropertyNames::STOP_WKSP_INDEX);
463 if (startIndex > stopIndex) {
464 const std::string msg =
465 PropertyNames::START_WKSP_INDEX +
" must be less than or equal to " + PropertyNames::STOP_WKSP_INDEX;
466 issues[PropertyNames::START_WKSP_INDEX] = msg;
467 issues[PropertyNames::STOP_WKSP_INDEX] = msg;
472 bool haveCommonPeakParameters(
false);
473 std::vector<string> suppliedParameterNames =
getProperty(PropertyNames::PEAK_PARAM_NAMES);
474 std::vector<double> peakParamValues =
getProperty(PropertyNames::PEAK_PARAM_VALUES);
475 if ((!suppliedParameterNames.empty()) || (!peakParamValues.empty())) {
476 haveCommonPeakParameters =
true;
477 if (suppliedParameterNames.size() != peakParamValues.size()) {
478 issues[PropertyNames::PEAK_PARAM_NAMES] =
"must have same number of values as PeakParameterValues";
479 issues[PropertyNames::PEAK_PARAM_VALUES] =
"must have same number of values as PeakParameterNames";
484 std::string partablename =
getPropertyValue(PropertyNames::PEAK_PARAM_TABLE);
485 if (!partablename.empty()) {
486 if (haveCommonPeakParameters) {
487 const std::string msg =
"Parameter value table and initial parameter "
488 "name/value vectors cannot be given "
490 issues[PropertyNames::PEAK_PARAM_TABLE] = msg;
491 issues[PropertyNames::PEAK_PARAM_NAMES] = msg;
492 issues[PropertyNames::PEAK_PARAM_VALUES] = msg;
500 if (!suppliedParameterNames.empty()) {
503 std::dynamic_pointer_cast<IPeakFunction>(API::FunctionFactory::Instance().createFunction(peakfunctiontype));
506 std::vector<string> functionParameterNames;
508 functionParameterNames.emplace_back(
m_peakFunction->parameterName(i));
511 const bool failed = std::any_of(suppliedParameterNames.cbegin(), suppliedParameterNames.cend(),
512 [&functionParameterNames](
const auto &parName) {
513 return std::find(functionParameterNames.begin(), functionParameterNames.end(),
514 parName) == functionParameterNames.end();
517 std::string msg =
"Specified invalid parameter for peak function";
518 if (haveCommonPeakParameters)
519 issues[PropertyNames::PEAK_PARAM_NAMES] = msg;
521 issues[PropertyNames::PEAK_PARAM_TABLE] = msg;
526 const std::string error_table_name =
getPropertyValue(PropertyNames::OUTPUT_WKSP_PARAM_ERRS);
527 if (!error_table_name.empty()) {
528 const bool use_raw_params =
getProperty(PropertyNames::RAW_PARAMS);
529 if (!use_raw_params) {
530 issues[PropertyNames::OUTPUT_WKSP_PARAM_ERRS] =
"Cannot be used with " + PropertyNames::RAW_PARAMS +
"=False";
531 issues[PropertyNames::RAW_PARAMS] =
532 "Cannot be False with " + PropertyNames::OUTPUT_WKSP_PARAM_ERRS +
" specified";
571 int start_wi =
getProperty(PropertyNames::START_WKSP_INDEX);
575 int stop_wi =
getProperty(PropertyNames::STOP_WKSP_INDEX);
599 throw std::runtime_error(
"number of peaks to fit is zero.");
605 throw std::runtime_error(
"PeakWidthPercent must be less than 1");
610 double temp =
getProperty(PropertyNames::BACKGROUND_Z_SCORE);
612 std::ostringstream os;
613 os <<
"FitPeaks property \"" << PropertyNames::BACKGROUND_Z_SCORE <<
"\" is deprecated and will be ignored."
645 std::dynamic_pointer_cast<IPeakFunction>(API::FunctionFactory::Instance().createFunction(peakfunctiontype));
649 std::string bkgdname;
650 if (bkgdfunctiontype ==
"Linear")
651 bkgdname =
"LinearBackground";
652 else if (bkgdfunctiontype ==
"Flat") {
653 g_log.
warning(
"There may be problems with Flat background");
654 bkgdname =
"FlatBackground";
656 bkgdname = bkgdfunctiontype;
658 std::dynamic_pointer_cast<IBackgroundFunction>(API::FunctionFactory::Instance().createFunction(bkgdname));
661 API::FunctionFactory::Instance().createFunction(
"LinearBackground"));
667 std::string partablename =
getPropertyValue(PropertyNames::PEAK_PARAM_TABLE);
680 }
else if (peakfunctiontype !=
"Gaussian") {
682 g_log.
warning(
"Neither parameter value table nor initial "
683 "parameter name/value vectors is specified. Fitting might "
684 "not be reliable for peak profile other than Gaussian");
696 std::vector<double> peakwindow =
getProperty(PropertyNames::FIT_WINDOW_LIST);
697 std::string peakwindowname =
getPropertyValue(PropertyNames::FIT_WINDOW_WKSP);
702 if ((!peakwindow.empty()) && peakwindowname.empty()) {
707 throw std::invalid_argument(
708 "Specifying peak windows with a list requires also specifying peak positions with a list.");
711 throw std::invalid_argument(
"Peak window vector must be twice as large as number of peaks.");
716 std::vector<double> peakranges(2);
717 peakranges[0] = peakwindow[i * 2];
718 peakranges[1] = peakwindow[i * 2 + 1];
725 std::stringstream errss;
726 errss <<
"Peak " << i <<
": user specifies an invalid range and peak center against " << peakranges[0] <<
" < "
728 throw std::invalid_argument(errss.str());
731 m_getPeakFitWindow = [
this](std::size_t wi, std::size_t ipeak) -> std::pair<double, double> {
740 }
else if (peakwindow.empty() && peakwindowws !=
nullptr) {
748 if (peakWindowX.empty()) {
749 std::stringstream errss;
750 errss <<
"Peak window required at workspace index " << wi <<
" "
751 <<
"which is undefined in the peak window workspace. "
752 <<
"Ensure workspace indices correspond in peak window workspace and input workspace "
753 <<
"when using start and stop indices.";
754 throw std::invalid_argument(errss.str());
757 if (peakWindowX.size() % 2 != 0) {
758 throw std::invalid_argument(
"The peak window vector must be even, with two edges for each peak center.");
760 if (peakWindowX.size() != peakCenterX.size() * 2) {
761 std::stringstream errss;
762 errss <<
"Peak window workspace index " << wi <<
" has incompatible number of fit windows "
763 << peakWindowX.size() / 2 <<
" with the number of peaks " << peakCenterX.size() <<
" to fit.";
764 throw std::invalid_argument(errss.str());
767 for (
size_t ipeak = 0; ipeak < peakCenterX.size(); ++ipeak) {
768 double left_w_bound = peakWindowX[ipeak * 2];
769 double right_w_bound = peakWindowX[ipeak * 2 + 1];
770 double center = peakCenterX[ipeak];
772 if (!(left_w_bound < center && center < right_w_bound)) {
773 std::stringstream errss;
774 errss <<
"Workspace index " << wi <<
" has incompatible peak window "
775 <<
"(" << left_w_bound <<
", " << right_w_bound <<
") "
776 <<
"with " << ipeak <<
"-th expected peak's center " << center;
777 throw std::runtime_error(errss.str());
781 m_getPeakFitWindow = [
this](std::size_t wi, std::size_t ipeak) -> std::pair<double, double> {
790 }
else if (peakwindow.empty()) {
796 m_getPeakFitWindow = [
this](std::size_t wi, std::size_t ipeak) -> std::pair<double, double> {
805 double left = peak_pos - estimate_peak_width * THREE;
806 double right = peak_pos + estimate_peak_width * THREE;
811 throw std::invalid_argument(
"Without definition of peak window, the "
812 "input workspace must be in unit of dSpacing "
813 "and Delta(D)/D must be given!");
817 throw std::invalid_argument(
"One and only one of peak window array and "
818 "peak window workspace can be specified.");
838 g_log.
notice(
"Peak centers are not specified by peak center workspace");
840 std::string peakpswsname =
getPropertyValue(PropertyNames::PEAK_CENTERS_WKSP);
850 }
else if (
m_peakCenters.empty() && peakcenterws !=
nullptr) {
860 std::stringstream errss;
863 <<
"However, the peak center workspace does not have values defined "
864 <<
"at workspace index " << wi <<
". "
865 <<
"Make sure the workspace indices between input and peak center workspaces correspond.";
867 throw std::invalid_argument(errss.str());
877 std::stringstream errss;
878 errss <<
"One and only one in 'PeakCenters' (vector) and "
879 "'PeakCentersWorkspace' shall be given. "
880 <<
"'PeakCenters' has size " <<
m_peakCenters.size() <<
", and name of peak center workspace "
881 <<
"is " << peakpswsname;
882 throw std::invalid_argument(errss.str());
897 throw std::runtime_error(
"ProcessInputPeakTolerance() must be called after "
898 "ProcessInputPeakCenters()");
915 throw std::runtime_error(
"Number of peak position tolerances and number of "
916 "peaks to fit are inconsistent.");
954 std::map<std::string, size_t> parname_index_map;
955 for (
size_t iparam = 0; iparam <
m_peakFunction->nParams(); ++iparam)
956 parname_index_map.insert(std::make_pair(
m_peakFunction->parameterName(iparam), iparam));
964 auto locator = parname_index_map.find(paramName);
965 if (locator != parname_index_map.end()) {
971 <<
" is not an allowed parameter of peak "
989 std::vector<std::shared_ptr<FitPeaksAlgorithm::PeakFitResult>> fit_result_vector(
m_numSpectraToFit);
991 const int nThreads = FrameworkManager::Instance().getNumOMPThreads();
994 std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> pre_check_result =
995 std::make_shared<FitPeaksAlgorithm::PeakFitPreCheckResult>();
997 PRAGMA_OMP(parallel
for schedule(dynamic, 1) )
998 for (
int ithread = 0; ithread < nThreads; ithread++) {
1004 std::vector<std::vector<double>> lastGoodPeakParameters(
m_numPeaksToFit,
1009 for (
auto wi = iws_begin; wi < iws_end; ++wi) {
1015 std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> fit_result =
1016 std::make_shared<FitPeaksAlgorithm::PeakFitResult>(
m_numPeaksToFit, numfuncparams);
1018 std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> spectrum_pre_check_result =
1019 std::make_shared<FitPeaksAlgorithm::PeakFitPreCheckResult>();
1021 fitSpectrumPeaks(
static_cast<size_t>(wi), expected_peak_centers, fit_result, lastGoodPeakParameters,
1022 lastGoodPeakSpectra, spectrum_pre_check_result);
1025 writeFitResult(
static_cast<size_t>(wi), expected_peak_centers, fit_result);
1027 *pre_check_result += *spectrum_pre_check_result;
1035 return fit_result_vector;
1040bool estimateBackgroundParameters(
const Histogram &histogram,
const std::pair<size_t, size_t> &peak_window,
1043 std::vector<double> &vec_y);
1044template <
typename vector_like>
1045void rangeToIndexBounds(
const vector_like &vecx,
const double range_left,
const double range_right,
size_t &left_index,
1046 size_t &right_index);
1049std::vector<std::string> supported_peak_profiles{
"Gaussian",
"Lorentzian",
"PseudoVoigt",
"Voigt",
1050 "BackToBackExponential"};
1057double estimateBackgroundNoise(
const std::vector<double> &vec_y) {
1059 size_t half_number_of_bkg_datapoints{5};
1060 if (vec_y.size() < 2 * half_number_of_bkg_datapoints + 3 )
1065 std::vector<double> vec_bkg;
1066 vec_bkg.resize(2 * half_number_of_bkg_datapoints);
1067 std::copy(vec_y.begin(), vec_y.begin() + half_number_of_bkg_datapoints, vec_bkg.begin());
1068 std::copy(vec_y.end() - half_number_of_bkg_datapoints, vec_y.end(), vec_bkg.begin() + half_number_of_bkg_datapoints);
1072 std::vector<double> vec_bkg_no_outliers;
1073 vec_bkg_no_outliers.resize(vec_bkg.size());
1074 double zscore_crit = 3.;
1075 for (
size_t ii = 0; ii < vec_bkg.size(); ii++) {
1076 if (zscore_vec[ii] <= zscore_crit)
1077 vec_bkg_no_outliers.push_back(vec_bkg[ii]);
1080 if (vec_bkg_no_outliers.size() < half_number_of_bkg_datapoints)
1084 return intensityStatistics.standard_deviation;
1095template <
typename vector_like>
1096void rangeToIndexBounds(
const vector_like &elems,
const double range_left,
const double range_right,
size_t &left_index,
1097 size_t &right_index) {
1098 const auto left_iter = std::lower_bound(elems.cbegin(), elems.cend(), range_left);
1099 const auto right_iter = std::upper_bound(elems.cbegin(), elems.cend(), range_right);
1101 left_index = std::distance(elems.cbegin(), left_iter);
1102 right_index = std::distance(elems.cbegin(), right_iter);
1103 right_index = std::min(right_index, elems.size() - 1);
1113 std::vector<double> &vec_y) {
1117 bkgd_func->function(vectorx, vector_bkgd);
1120 for (
size_t i = 0; i < vec_y.size(); ++i) {
1121 (vec_y)[i] -= vector_bkgd[i];
1130class LoggingOffsetSentry {
1148 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result,
1149 std::vector<std::vector<double>> &lastGoodPeakParameters,
1150 std::vector<size_t> &lastGoodPeakSpectra,
1151 const std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> &pre_check_result) {
1157 fit_result->setBadRecord(i, -1.);
1158 pre_check_result->setNumberOfSpectrumPeaksWithLowCount(
m_numPeaksToFit);
1167 std::stringstream errss;
1168 errss <<
"The FitPeak algorithm requires the CurveFitting library";
1170 throw std::runtime_error(errss.str());
1177 peak_fitter->setProperty(
"Minimizer",
m_minimizer);
1179 peak_fitter->setProperty(
"CalcErrors",
true);
1186 bool neighborPeakSameSpectrum =
false;
1187 size_t number_of_out_of_range_peaks{0};
1188 for (
size_t fit_index = 0; fit_index <
m_numPeaksToFit; ++fit_index) {
1190 size_t peak_index(fit_index);
1195 for (
size_t i = 0; i < bkgdfunction->nParams(); ++i)
1196 bkgdfunction->setParameter(i, 0.);
1198 double expected_peak_pos = expected_peak_centers[peak_index];
1203 auto peakfunction = std::dynamic_pointer_cast<API::IPeakFunction>(
m_peakFunction->clone());
1204 peakfunction->setCentre(expected_peak_pos);
1205 peakfunction->setMatrixWorkspace(
m_inputMatrixWS, wi, peak_window_i.first, peak_window_i.second);
1207 std::map<size_t, double> keep_values;
1208 for (
size_t ipar = 0; ipar < peakfunction->nParams(); ++ipar) {
1209 if (peakfunction->isFixed(ipar)) {
1213 keep_values[ipar] = peakfunction->getParameter(ipar);
1216 peakfunction->unfix(ipar);
1222 bool samePeakCrossSpectrum = (lastGoodPeakParameters[peak_index].size() >
1223 static_cast<size_t>(std::count_if(lastGoodPeakParameters[peak_index].begin(),
1224 lastGoodPeakParameters[peak_index].end(),
1225 [&](
auto const &val) {
return val <= 1e-10; })));
1230 if (wi > 0 && samePeakCrossSpectrum) {
1231 size_t lastGoodWi = lastGoodPeakSpectra[peak_index];
1232 std::shared_ptr<const Geometry::Detector> pdetector =
1233 std::dynamic_pointer_cast<const Geometry::Detector>(
m_inputMatrixWS->getDetector(lastGoodWi));
1234 std::shared_ptr<const Geometry::Detector> cdetector =
1235 std::dynamic_pointer_cast<const Geometry::Detector>(
m_inputMatrixWS->getDetector(wi));
1238 if (pdetector && cdetector) {
1239 auto prev_id = pdetector->getID();
1240 auto curr_id = cdetector->getID();
1241 if (prev_id + 1 != curr_id)
1242 samePeakCrossSpectrum =
false;
1244 samePeakCrossSpectrum =
false;
1250 samePeakCrossSpectrum =
false;
1252 }
catch (
const std::runtime_error &) {
1255 samePeakCrossSpectrum =
false;
1259 if (samePeakCrossSpectrum) {
1261 for (
size_t i = 0; i < peakfunction->nParams(); ++i) {
1262 peakfunction->setParameter(i, lastGoodPeakParameters[peak_index][i]);
1266 for (
size_t i = 0; i < peakfunction->nParams(); ++i) {
1267 peakfunction->setParameter(i, lastGoodPeakParameters[prev_peak_index][i]);
1272 peakfunction->setCentre(expected_peak_pos);
1275 for (
const auto &[ipar,
value] : keep_values) {
1276 peakfunction->setParameter(ipar,
value);
1279 double cost(DBL_MAX);
1280 if (expected_peak_pos <= x0 || expected_peak_pos >= xf) {
1282 peakfunction->setIntensity(0);
1283 number_of_out_of_range_peaks++;
1293 auto useUserSpecifedIfGiven =
1298 g_log.
warning(
"Peak width can be estimated as ZERO. The result can be wrong");
1302 std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> peak_pre_check_result =
1303 std::make_shared<FitPeaksAlgorithm::PeakFitPreCheckResult>();
1304 cost =
fitIndividualPeak(wi, peak_fitter, expected_peak_pos, peak_window_i, observe_peak_width, peakfunction,
1305 bkgdfunction, peak_pre_check_result);
1306 if (peak_pre_check_result->isIndividualPeakRejected())
1307 fit_result->setBadRecord(peak_index, -1.);
1311 fit_result->setBadRecord(peak_index, -1.);
1316 *pre_check_result += *peak_pre_check_result;
1318 pre_check_result->setNumberOfOutOfRangePeaks(number_of_out_of_range_peaks);
1330 neighborPeakSameSpectrum =
true;
1331 prev_peak_index = peak_index;
1333 for (
size_t i = 0; i < lastGoodPeakParameters[peak_index].size(); ++i) {
1334 lastGoodPeakParameters[peak_index][i] = peakfunction->getParameter(i);
1336 lastGoodPeakSpectra[peak_index] = wi;
1360 if (firstPeakInSpectrum) {
1366 peak_function->setParameter(param_index, param_value);
1374 observe_peak_shape =
true;
1377 return observe_peak_shape;
1393 const std::vector<double> &expected_peak_positions,
1395 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result) {
1397 double postol(DBL_MAX);
1411 throw std::runtime_error(
"Peak tolerance out of index");
1419 bool good_fit(
false);
1420 if ((cost < 0) || (cost >= DBL_MAX - 1.) || std::isnan(cost)) {
1426 }
else if (case23) {
1429 if (fitwindow.first < fitwindow.second) {
1432 if (peak_pos < fitwindow.first || peak_pos > fitwindow.second) {
1435 g_log.
debug() <<
"Peak position " << peak_pos <<
" is out of fit "
1436 <<
"window boundary " << fitwindow.first <<
", " << fitwindow.second <<
"\n";
1437 }
else if (peak_fwhm > (fitwindow.second - fitwindow.first)) {
1440 g_log.
debug() <<
"Peak position " << peak_pos <<
" has fwhm "
1441 <<
"wider than the fit window " << fitwindow.second - fitwindow.first <<
"\n";
1447 double left_bound(-1);
1449 left_bound = 0.5 * (expected_peak_positions[peakindex] - expected_peak_positions[peakindex - 1]);
1450 double right_bound(-1);
1452 right_bound = 0.5 * (expected_peak_positions[peakindex + 1] - expected_peak_positions[peakindex]);
1454 left_bound = right_bound;
1455 if (right_bound < left_bound)
1456 right_bound = left_bound;
1457 if (left_bound < 0 || right_bound < 0)
1458 throw std::runtime_error(
"Code logic error such that left or right "
1459 "boundary of peak position is negative.");
1460 if (peak_pos < left_bound || peak_pos > right_bound) {
1462 }
else if (peak_fwhm > (right_bound - left_bound)) {
1465 g_log.
debug() <<
"Peak position " << peak_pos <<
" has fwhm "
1466 <<
"wider than the fit window " << right_bound - left_bound <<
"\n";
1471 }
else if (
fabs(fitfunction.
peakfunction->centre() - expected_peak_positions[peakindex]) > postol) {
1475 <<
fabs(fitfunction.
peakfunction->centre() - expected_peak_positions[peakindex])
1476 <<
" is out of range of tolerance: " << postol <<
"\n";
1483 double adjust_cost(cost);
1486 adjust_cost = DBL_MAX;
1490 if (adjust_cost > DBL_MAX - 1) {
1495 fit_result->setRecord(peakindex, adjust_cost, peak_pos, fitfunction);
1509 throw std::runtime_error(
"No parameters");
1517 std::size_t iws =
static_cast<std::size_t
>(iiws);
1521 std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> fit_result_i = fit_results[iws -
m_startWorkspaceIndex];
1524 throw std::runtime_error(
"There is something wroing with PeakFitResult vector!");
1528 const double chi2 = fit_result_i->getCost(ipeak);
1532 for (
size_t iparam = 0; iparam < num_peakfunc_params; ++iparam)
1533 peak_function->setParameter(iparam, fit_result_i->getParameterValue(ipeak, iparam));
1534 for (
size_t iparam = 0; iparam < num_bkgdfunc_params; ++iparam)
1535 bkgd_function->setParameter(iparam, fit_result_i->getParameterValue(ipeak, num_peakfunc_params + iparam));
1540 auto start_x_iter = std::lower_bound(vec_x.begin(), vec_x.end(), peakwindow.first);
1541 auto stop_x_iter = std::lower_bound(vec_x.begin(), vec_x.end(), peakwindow.second);
1543 if (start_x_iter == stop_x_iter)
1544 throw std::runtime_error(
"Range size is zero in calculateFittedPeaks");
1549 comp_func->addFunction(peak_function);
1550 comp_func->addFunction(bkgd_function);
1551 comp_func->function(domain, values);
1554 std::size_t istart =
static_cast<size_t>(start_x_iter - vec_x.begin());
1555 std::size_t istop =
static_cast<size_t>(stop_x_iter - vec_x.begin());
1556 for (std::size_t yindex = istart; yindex < istop; ++yindex) {
1570 auto startX = std::lower_bound(vecX.begin(), vecX.end(), peakWindow.first);
1571 auto stopX = std::lower_bound(vecX.begin(), vecX.end(), peakWindow.second);
1576 peakFunction->function(domain, values);
1577 auto peakValues = values.
toVector();
1580 auto startE = errors.begin() + (startX - vecX.begin());
1581 auto stopE = errors.begin() + (stopX - vecX.begin());
1582 std::vector<double> peakErrors(startE, stopE);
1584 double peakSum = std::accumulate(peakValues.cbegin(), peakValues.cend(), 0.0);
1591bool estimateBackgroundParameters(
const Histogram &histogram,
const std::pair<size_t, size_t> &peak_window,
1595 const auto POLYNOMIAL_ORDER = std::min<size_t>(1, bkgd_function->nParams());
1597 if (peak_window.first >= peak_window.second)
1598 throw std::runtime_error(
"Invalid peak window");
1601 const auto nParams = bkgd_function->nParams();
1602 for (
size_t i = 0; i < nParams; ++i)
1603 bkgd_function->setParameter(i, 0.);
1606 const size_t iback_start = peak_window.first + 10;
1607 const size_t iback_stop = peak_window.second - 10;
1612 if (iback_start < iback_stop) {
1616 double chisq{DBL_MAX};
1617 HistogramData::estimateBackground(POLYNOMIAL_ORDER, histogram, peak_window.first, peak_window.second, iback_start,
1618 iback_stop, bkgd_a0, bkgd_a1, bkgd_a2, chisq);
1620 bkgd_function->setParameter(0, bkgd_a0);
1622 bkgd_function->setParameter(1, bkgd_a1);
1640 return (std::find(supported_peak_profiles.begin(), supported_peak_profiles.end(), peakprofile) !=
1641 supported_peak_profiles.end());
1649 constexpr size_t MIN_POINTS{10};
1653 size_t start_index =
findXIndex(points.rawData(), fit_window.first);
1654 size_t expected_peak_index =
findXIndex(points.rawData(), expected_peak_pos, start_index);
1655 size_t stop_index =
findXIndex(points.rawData(), fit_window.second, expected_peak_index);
1658 bool good_fit(
false);
1659 if (expected_peak_index - start_index > MIN_POINTS && stop_index - expected_peak_index > MIN_POINTS) {
1662 const std::pair<double, double> vec_min{fit_window.first, points[expected_peak_index + 5]};
1663 const std::pair<double, double> vec_max{points[expected_peak_index - 5], fit_window.second};
1666 for (
size_t n = 0;
n < bkgd_func->nParams(); ++
n)
1667 bkgd_func->setParameter(
n, 0);
1672 if (chi2 < DBL_MAX - 1) {
1680 g_log.
debug() <<
"Don't know what to do with background fitting with single "
1681 <<
"domain function! " << (expected_peak_index - start_index) <<
" points to the left "
1682 << (stop_index - expected_peak_index) <<
" points to the right\n";
1692 const std::pair<double, double> &fitwindow,
const bool estimate_peak_width,
1695 const std::shared_ptr<FitPeaksAlgorithm::PeakFitPreCheckResult> &pre_check_result) {
1696 pre_check_result->setNumberOfSubmittedIndividualPeaks(1);
1697 double cost(DBL_MAX);
1700 size_t min_required_datapoints{peakfunction->nParams() + bkgdfunc->nParams() + 2};
1702 if (number_of_datapoints < min_required_datapoints) {
1703 pre_check_result->setNumberOfPeaksWithNotEnoughDataPoints(1);
1709 pre_check_result->setNumberOfIndividualPeaksWithLowCount(1);
1715 pre_check_result->setNumberOfPeaksWithLowSignalToNoise(1);
1726 estimate_peak_width,
true);
1741 const std::pair<double, double> &peak_range,
const double &expected_peak_center,
1742 bool estimate_peak_width,
bool estimate_background) {
1743 std::stringstream errorid;
1744 errorid <<
"(WorkspaceIndex=" << wsindex <<
" PeakCentre=" << expected_peak_center <<
")";
1747 if (peak_range.first >= peak_range.second) {
1748 std::stringstream msg;
1749 msg <<
"Invalid peak window: xmin>xmax (" << peak_range.first <<
", " << peak_range.second <<
")" << errorid.str();
1750 throw std::runtime_error(msg.str());
1754 const auto &histogram = dataws->histogram(wsindex);
1755 const auto &vector_x = histogram.points();
1756 const auto start_index =
findXIndex(vector_x, peak_range.first);
1757 const auto stop_index =
findXIndex(vector_x, peak_range.second, start_index);
1758 if (start_index == stop_index)
1759 throw std::runtime_error(
"Range size is zero in fitFunctionSD");
1760 std::pair<size_t, size_t> peak_index_window = std::make_pair(start_index, stop_index);
1763 if (estimate_background) {
1764 if (!estimateBackgroundParameters(histogram, peak_index_window, bkgd_function)) {
1770 peak_function->setCentre(expected_peak_center);
1771 int result =
estimatePeakParameters(histogram, peak_index_window, peak_function, bkgd_function, estimate_peak_width,
1774 if (result !=
GOOD) {
1775 peak_function->setCentre(expected_peak_center);
1783 comp_func->addFunction(peak_function);
1784 comp_func->addFunction(bkgd_function);
1785 IFunction_sptr fitfunc = std::dynamic_pointer_cast<IFunction>(comp_func);
1788 fit->setProperty(
"Function", fitfunc);
1789 fit->setProperty(
"InputWorkspace", dataws);
1790 fit->setProperty(
"WorkspaceIndex",
static_cast<int>(wsindex));
1792 fit->setProperty(
"StartX", peak_range.first);
1793 fit->setProperty(
"EndX", peak_range.second);
1794 fit->setProperty(
"IgnoreInvalidData",
true);
1798 double peak_center = peak_function->centre();
1799 double peak_width = peak_function->fwhm();
1800 std::stringstream peak_center_constraint;
1801 peak_center_constraint << (peak_center - 0.5 * peak_width) <<
" < f0." << peak_function->getCentreParameterName()
1802 <<
" < " << (peak_center + 0.5 * peak_width);
1803 fit->setProperty(
"Constraints", peak_center_constraint.str());
1807 g_log.
debug() <<
"[E1201] FitSingleDomain Before fitting, Fit function: " << fit->asString() <<
"\n";
1808 errorid <<
" starting function [" << comp_func->asString() <<
"]";
1811 g_log.
debug() <<
"[E1202] FitSingleDomain After fitting, Fit function: " << fit->asString() <<
"\n";
1813 if (!fit->isExecuted()) {
1814 g_log.
warning() <<
"Fitting peak SD (single domain) failed to execute. " + errorid.str();
1817 }
catch (std::invalid_argument &e) {
1818 errorid <<
": " << e.what();
1824 std::string fitStatus = fit->getProperty(
"OutputStatus");
1825 double chi2{std::numeric_limits<double>::max()};
1826 if (fitStatus ==
"success") {
1827 chi2 = fit->getProperty(
"OutputChi2overDoF");
1835 const size_t wsindex,
const std::pair<double, double> &vec_xmin,
1836 const std::pair<double, double> &vec_xmax) {
1842 std::stringstream errss;
1843 errss <<
"The FitPeak algorithm requires the CurveFitting library";
1844 throw std::runtime_error(errss.str());
1849 fit->setProperty(
"CalcErrors",
true);
1853 std::shared_ptr<MultiDomainFunction> md_function = std::make_shared<MultiDomainFunction>();
1856 md_function->addFunction(std::move(fit_function));
1859 md_function->clearDomainIndices();
1860 md_function->setDomainIndices(0, {0, 1});
1863 fit->setProperty(
"Function", std::dynamic_pointer_cast<IFunction>(md_function));
1864 fit->setProperty(
"InputWorkspace", dataws);
1865 fit->setProperty(
"WorkspaceIndex",
static_cast<int>(wsindex));
1866 fit->setProperty(
"StartX", vec_xmin.first);
1867 fit->setProperty(
"EndX", vec_xmax.first);
1868 fit->setProperty(
"InputWorkspace_1", dataws);
1869 fit->setProperty(
"WorkspaceIndex_1",
static_cast<int>(wsindex));
1870 fit->setProperty(
"StartX_1", vec_xmin.second);
1871 fit->setProperty(
"EndX_1", vec_xmax.second);
1873 fit->setProperty(
"IgnoreInvalidData",
true);
1877 if (!fit->isExecuted()) {
1878 throw runtime_error(
"Fit is not executed on multi-domain function/data. ");
1882 std::string fitStatus = fit->getProperty(
"OutputStatus");
1884 double chi2 = DBL_MAX;
1885 if (fitStatus ==
"success") {
1886 chi2 = fit->getProperty(
"OutputChi2overDoF");
1895 const size_t &ws_index,
const double &expected_peak_center,
1905 fitBackground(ws_index, fit_window, expected_peak_center, high_bkgd_function);
1908 std::vector<double> vec_x, vec_y, vec_e;
1909 getRangeData(ws_index, fit_window, vec_x, vec_y, vec_e);
1912 reduceByBackground(high_bkgd_function, vec_x, vec_y);
1913 for (std::size_t
n = 0;
n < bkgdfunc->nParams(); ++
n)
1914 bkgdfunc->setParameter(
n, 0);
1920 fitFunctionSD(fit, peakfunction, bkgdfunc, reduced_bkgd_ws, 0, {vec_x.front(), vec_x.back()}, expected_peak_center,
1921 observe_peak_shape,
false);
1924 bkgdfunc->setParameter(0, bkgdfunc->getParameter(0) + high_bkgd_function->getParameter(0));
1925 bkgdfunc->setParameter(1, bkgdfunc->getParameter(1) +
1926 high_bkgd_function->getParameter(1));
1929 expected_peak_center,
false,
false);
1937 const std::vector<double> &vec_y,
1938 const std::vector<double> &vec_e) {
1939 std::size_t size = vec_x.size();
1940 std::size_t ysize = vec_y.size();
1942 HistogramBuilder builder;
1944 builder.setY(ysize);
1947 auto &dataX = matrix_ws->mutableX(0);
1948 auto &dataY = matrix_ws->mutableY(0);
1949 auto &dataE = matrix_ws->mutableE(0);
1951 dataX.assign(vec_x.cbegin(), vec_x.cend());
1952 dataY.assign(vec_y.cbegin(), vec_y.cend());
1953 dataE.assign(vec_e.cbegin(), vec_e.cend());
1969 for (std::size_t ipeak = 0; ipeak < expected_position.size(); ++ipeak) {
1985 const std::vector<std::string> ¶m_names,
bool with_chi2) {
1987 table_ws->addColumn(
"int",
"wsindex");
1988 table_ws->addColumn(
"int",
"peakindex");
1989 for (
const auto ¶m_name : param_names)
1990 table_ws->addColumn(
"double", param_name);
1992 table_ws->addColumn(
"double",
"chi2");
1999 newRow << static_cast<int>(iws);
2000 newRow << static_cast<int>(ipeak);
2001 for (
size_t iparam = 0; iparam < numParam; ++iparam)
2022 std::vector<std::string> param_vec;
2026 param_vec.emplace_back(
"centre");
2027 param_vec.emplace_back(
"width");
2028 param_vec.emplace_back(
"height");
2029 param_vec.emplace_back(
"intensity");
2032 for (
size_t iparam = 0; iparam <
m_bkgdFunction->nParams(); ++iparam)
2040 std::string fiterror_table_name =
getPropertyValue(PropertyNames::OUTPUT_WKSP_PARAM_ERRS);
2042 if (fiterror_table_name.empty()) {
2060 std::string fit_ws_name =
getPropertyValue(PropertyNames::OUTPUT_WKSP_MODEL);
2061 if (fit_ws_name.size() == 0) {
2086 g_log.
debug(
"about to calcualte fitted peaks");
2098 const std::vector<double> &vec_y =
m_inputMatrixWS->histogram(iws).y().rawData();
2099 double total = std::accumulate(vec_y.begin(), vec_y.end(), 0.);
2111 std::vector<double> vec_x, vec_y, vec_e;
2114 double total = std::accumulate(vec_y.begin(), vec_y.end(), 0.);
2125 size_t left_index, right_index;
2127 size_t number_dp = right_index - left_index + 1;
2130 assert(number_dp > 0);
2143void
FitPeaks::histRangeToIndexBounds(
size_t iws, const
std::pair<
double,
double> &range,
size_t &left_index,
2144 size_t &right_index) {
2145 const auto &orig_x = m_inputMatrixWS->histogram(iws).x();
2146 rangeToIndexBounds(orig_x, range.first, range.second, left_index, right_index);
2150 if (left_index >= right_index || (m_inputMatrixWS->isHistogramData() && left_index == right_index - 1)) {
2151 std::stringstream err_ss;
2152 err_ss <<
"Unable to get a valid subset of histogram from given fit window. "
2153 <<
"Histogram X: " << orig_x.front() <<
"," << orig_x.back() <<
"; Range: " << range.first <<
","
2155 throw std::runtime_error(err_ss.str());
2168 std::vector<double> &vec_y, std::vector<double> &vec_e) {
2170 size_t left_index, right_index;
2174 size_t num_elements_x = right_index - left_index;
2176 vec_x.resize(num_elements_x);
2178 std::copy(orig_x.begin() + left_index, orig_x.begin() + right_index, vec_x.begin());
2180 size_t num_datapoints =
m_inputMatrixWS->isHistogramData() ? num_elements_x - 1 : num_elements_x;
2182 const std::vector<double> orig_y =
m_inputMatrixWS->histogram(iws).y().rawData();
2183 const std::vector<double> orig_e =
m_inputMatrixWS->histogram(iws).e().rawData();
2184 vec_y.resize(num_datapoints);
2185 vec_e.resize(num_datapoints);
2186 std::copy(orig_y.begin() + left_index, orig_y.begin() + left_index + num_datapoints, vec_y.begin());
2187 std::copy(orig_e.begin() + left_index, orig_e.begin() + left_index + num_datapoints, vec_e.begin());
2199double
FitPeaks::calculateSignalToNoiseRatio(
size_t iws, const
std::pair<
double,
double> &range,
2202 size_t left_index, right_index;
2203 histRangeToIndexBounds(iws, range, left_index, right_index);
2206 if (!estimateBackgroundParameters(m_inputMatrixWS->histogram(iws), std::pair<size_t, size_t>(left_index, right_index),
2211 std::vector<double> vec_x, vec_y, vec_e;
2212 getRangeData(iws, range, vec_x, vec_y, vec_e);
2217 reduceByBackground(bkgd_function, vec_x, vec_y);
2220 auto it_max = std::max_element(vec_y.begin(), vec_y.end());
2221 double signal = vec_y[it_max - vec_y.begin()];
2222 if (signal <= DBL_MIN)
2227 double noise = estimateBackgroundNoise(vec_y);
2228 if (noise <= DBL_MIN)
2232 return signal / noise;
2240 std::stringstream errss;
2241 errss <<
"Workspace index " << wi <<
" is out of range "
2243 throw std::runtime_error(errss.str());
2250 std::stringstream errss;
2251 errss <<
"Peak index " << ipeak <<
" is out of range (" <<
m_numPeaksToFit <<
")";
2252 throw std::runtime_error(errss.str());
2258 std::stringstream errss;
2259 errss <<
"Peak window is inappropriate for workspace index: " <<
left <<
" >= " <<
right;
2260 throw std::runtime_error(errss.str());
2274 const std::shared_ptr<FitPeaksAlgorithm::PeakFitResult> &fit_result) {
2278 g_log.
error() <<
"workspace index " << wi <<
" is out of output peak position workspace "
2282 throw std::runtime_error(
"Out of boundary to set output peak position workspace");
2287 double exp_peak_pos(expected_positions[ipeak]);
2288 double fitted_peak_pos = fit_result->getPeakPosition(ipeak);
2289 double peak_chi2 = fit_result->getCost(ipeak);
2308 <<
" parameters. Parameter table shall have 3 more "
2309 "columns. But not it has "
2311 throw std::runtime_error(
"Peak parameter vector for one peak has different sizes to output "
2318 std::stringstream err_ss;
2319 err_ss <<
"Peak has 4 effective peak parameters and " <<
m_bkgdFunction->nParams() <<
" background parameters "
2320 <<
". Parameter table shall have 3 more columns. But not it has " <<
m_fittedParamTable->columnCount()
2322 throw std::runtime_error(err_ss.str());
2329 size_t num_peakfunc_params = peak_function->nParams();
2339 for (
size_t iparam = 0; iparam < num_peakfunc_params + num_bkgd_params; ++iparam) {
2340 size_t col_index = iparam + 2;
2342 m_fittedParamTable->cell<
double>(row_index, col_index) = fit_result->getParameterValue(ipeak, iparam);
2345 m_fitErrorTable->cell<
double>(row_index, col_index) = fit_result->getParameterError(ipeak, iparam);
2352 for (
size_t iparam = 0; iparam < num_peakfunc_params; ++iparam)
2353 peak_function->setParameter(iparam, fit_result->getParameterValue(ipeak, iparam));
2362 for (
size_t iparam = 0; iparam < num_bkgd_params; ++iparam)
2364 fit_result->getParameterValue(ipeak, num_peakfunc_params + iparam);
2368 m_fittedParamTable->cell<
double>(row_index, chi2_index) = fit_result->getCost(ipeak);
2376 std::string height_name(
"");
2378 std::vector<std::string> peak_parameters = peak_function->getParameterNames();
2379 for (
const auto &parName : peak_parameters) {
2380 if (parName ==
"Height") {
2381 height_name =
"Height";
2383 }
else if (parName ==
"I") {
2386 }
else if (parName ==
"Intensity") {
2387 height_name =
"Intensity";
2392 if (height_name.empty())
2393 throw std::runtime_error(
"Peak height parameter name cannot be found.");
2401 LoggingOffsetSentry sentry(
this);
#define DECLARE_ALGORITHM(classname)
double value
The value of the point.
#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_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.
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
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.
TableRow represents a row in a TableWorkspace.
A property class for workspaces.
size_t m_submitted_spectrum_peaks
void setNumberOfSpectrumPeaksWithLowCount(const size_t n)
std::string getReport() const
size_t m_low_count_individual
PeakFitPreCheckResult & operator+=(const PeakFitPreCheckResult &another)
size_t m_not_enough_datapoints
void setNumberOfSubmittedIndividualPeaks(const size_t n)
void setNumberOfPeaksWithNotEnoughDataPoints(const size_t n)
size_t m_low_count_spectrum
void setNumberOfOutOfRangePeaks(const size_t n)
void setNumberOfPeaksWithLowSignalToNoise(const size_t n)
size_t m_submitted_individual_peaks
void setNumberOfIndividualPeaksWithLowCount(const size_t n)
void setNumberOfSubmittedSpectrumPeaks(const size_t n)
bool isIndividualPeakRejected() const
size_t m_function_parameters_number
number of function parameters
double getPeakPosition(size_t ipeak) const
size_t getNumberPeaks() const
std::vector< std::vector< double > > m_function_parameters_vector
PeakFitResult(size_t num_peaks, size_t num_params)
Holds all of the fitting information for a single spectrum.
double getParameterValue(size_t ipeak, size_t iparam) const
get the fitted value of a particular parameter
std::vector< double > m_costs
std::vector< std::vector< double > > m_function_errors_vector
fitted peak and background parameters' fitting error
size_t getNumberParameters() const
double getCost(size_t ipeak) const
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
double getParameterError(size_t ipeak, size_t iparam) const
get the fitting error of a particular parameter
void setBadRecord(size_t ipeak, const double peak_position)
The peak postition should be negative and indicates what went wrong.
std::vector< double > m_fitted_peak_positions
Algorithms::PeakParameterHelper::EstimatePeakWidth m_peakWidthEstimateApproach
Flag for observing peak width: there are 3 states (1) no estimation (2) from 'observation' (3) calcul...
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
void generateFittedParametersValueWorkspaces()
Generate output workspaces.
void setupParameterTableWorkspace(const API::ITableWorkspace_sptr &table_ws, const std::vector< std::string > ¶m_names, bool with_chi2)
Set up parameter table (parameter value or error)
bool m_copyLastGoodPeakParameters
std::vector< double > m_peakPosTolerances
tolerances for fitting peak positions
API::IPeakFunction_sptr m_peakFunction
Peak profile name.
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
API::MatrixWorkspace_sptr m_fittedPeakWS
matrix workspace contained calcalated peaks+background from fitted result it has same number of spect...
API::ITableWorkspace_const_sptr m_profileStartingValueTable
table workspace for profile parameters' starting value
std::string m_minimizer
Minimzer.
API::IBackgroundFunction_sptr m_linearBackgroundFunction
Linear background function for high background fitting.
bool m_constrainPeaksPosition
bool m_peakPosTolCase234
peak positon tolerance case b, c and d
void fitSpectrumPeaks(size_t wi, const std::vector< double > &expected_peak_centers, const std::shared_ptr< FitPeaksAlgorithm::PeakFitResult > &fit_result, std::vector< std::vector< double > > &lastGoodPeakParameters, std::vector< size_t > &lastGoodPeakSpectra, const std::shared_ptr< FitPeaksAlgorithm::PeakFitPreCheckResult > &pre_check_result)
fit peaks in a same spectrum
std::map< std::string, std::string > validateInputs() override
Validate inputs.
double m_peakWidthPercentage
flag to estimate peak width from
API::IBackgroundFunction_sptr m_bkgdFunction
Background function.
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
std::vector< std::shared_ptr< FitPeaksAlgorithm::PeakFitResult > > fitPeaks()
suites of method to fit peaks
double m_minPeakHeight
minimum peak height without background and it also serves as the criteria for observed peak parameter
std::string m_costFunction
Cost function.
void processInputFunctions()
process inputs for peak and background functions
void processInputPeakTolerance()
process inputs about fitted peak positions' tolerance
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.
void init() override
Init.
void logNoOffset(const size_t &priority, const std::string &msg)
std::size_t m_numSpectraToFit
total number of spectra to be fit
std::vector< std::string > m_peakParamNames
input peak parameters' names
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
std::size_t m_numPeaksToFit
the number of peaks to fit in all spectra
std::size_t m_startWorkspaceIndex
start index
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)
double m_minPeakTotalCount
std::string getPeakHeightParameterName(const API::IPeakFunction_const_sptr &peak_function)
Get the parameter name for peak height (I or height or etc)
bool m_uniformPeakPositions
API::ITableWorkspace_sptr m_fittedParamTable
output analysis workspaces table workspace for fitted parameters
double m_minSignalToNoiseRatio
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.
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
std::vector< double > m_initParamValues
input peak parameters' starting values corresponding to above peak parameter names
std::vector< double > m_peakCenters
Designed peak positions and tolerance.
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
std::function< std::pair< double, double >(std::size_t const &, std::size_t const &)> m_getPeakFitWindow
API::MatrixWorkspace_sptr m_inputMatrixWS
mandatory input and output workspaces
std::vector< size_t > m_initParamIndexes
input starting parameters' indexes in peak function
bool m_fitPeaksFromRight
Fit from right or left.
std::function< std::vector< double >(std::size_t const &)> m_getExpectedPeakPositions
bool m_rawPeaksTable
flag to show that the pamarameters in table are raw parameters or effective parameters
void processInputs()
process inputs (main and child algorithms)
void checkWorkspaceIndices(std::size_t const &)
Get the expected peak's position.
void checkPeakWindowEdgeOrder(double const &, double const &)
int m_fitIterations
Fit iterations.
double m_minSignalToSigmaRatio
double numberCounts(size_t iws)
sum up all counts in histogram
API::MatrixWorkspace_const_sptr m_peakWindowWorkspace
bool m_uniformProfileStartingValue
flag for profile startng value being uniform or not
API::ITableWorkspace_sptr m_fitErrorTable
table workspace for fitted parameters' fitting error. This is optional
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)
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.
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.
void notice(const std::string &msg)
Logs at notice level.
void error(const std::string &msg)
Logs at error level.
void warning(const std::string &msg)
Logs at warning level.
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
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.
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.
constexpr double EMPTY_DBL() noexcept
Returns what we consider an "empty" double within a property.
API::IPeakFunction_sptr peakfunction
API::IBackgroundFunction_sptr bkgdfunction
@ Input
An input workspace.
@ Output
An output workspace.
Functor to accumulate a sum of squares.