Mantid
Loading...
Searching...
No Matches
FitPowderDiffPeaks.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 +
8
12
13#include "MantidAPI/Column.h"
18#include "MantidAPI/IFunction.h"
21#include "MantidAPI/TableRow.h"
22#include "MantidAPI/TextAxis.h"
24
35
36#include <cmath>
37#include <utility>
38
39#include <gsl/gsl_sf_erf.h>
40
42#define PEAKFITRANGEFACTOR 5.0
43
45#define PEAKBOUNDARYFACTOR 2.5
46
48#define EXCLUDEPEAKRANGEFACTOR 8.0
50#define WINDOWSIZE 3.0
51
52using namespace Mantid;
53using namespace Mantid::API;
54using namespace Mantid::Kernel;
55using namespace Mantid::DataObjects;
58
59using namespace std;
60
62
64
65//----------------------------------------------------------------------------------------------
69 : m_wsIndex(-1), m_tofMin(0.), m_tofMax(0.), m_useGivenTOFh(false), m_confidentInInstrumentParameters(false),
70 m_minimumHKL(), m_numPeaksLowerToMin(-1), m_indexGoodFitPeaks(), m_chi2GoodFitPeaks(), m_fitMode(ROBUSTFIT),
71 m_genPeakStartingValue(HKLCALCULATION), m_rightmostPeakHKL(), m_rightmostPeakLeftBound(0.),
72 m_rightmostPeakRightBound(0.), m_minPeakHeight(0.), m_unitCell(), m_fitPeakBackgroundComposite(false) {}
73
74//----------------------------------------------------------------------------------------------
78 // Input data workspace
79 declareProperty(std::make_unique<WorkspaceProperty<MatrixWorkspace>>("InputWorkspace", "Anonymous", Direction::Input),
80 "Input workspace for data (diffraction pattern). ");
81
82 // Output workspace
83 declareProperty(std::make_unique<WorkspaceProperty<Workspace2D>>("OutputWorkspace", "Anonymous2", Direction::Output),
84 "Output Workspace2D for the fitted peaks. ");
85
86 // Input/output peaks table workspace
87 declareProperty(std::make_unique<WorkspaceProperty<TableWorkspace>>("BraggPeakParameterWorkspace", "AnonymousPeak",
89 "TableWorkspace containg all peaks' parameters.");
90
91 // Input and output instrument parameters table workspace
92 declareProperty(std::make_unique<WorkspaceProperty<TableWorkspace>>("InstrumentParameterWorkspace",
93 "AnonymousInstrument", Direction::InOut),
94 "TableWorkspace containg instrument's parameters.");
95
96 // Workspace to output fitted peak parameters
97 declareProperty(std::make_unique<WorkspaceProperty<TableWorkspace>>("OutputBraggPeakParameterWorkspace",
98 "AnonymousOut2", Direction::Output),
99 "Output TableWorkspace containing the fitted peak parameters for each "
100 "peak.");
101
102 // Data workspace containing fitted peak parameters
103 declareProperty(std::make_unique<WorkspaceProperty<Workspace2D>>("OutputBraggPeakParameterDataWorkspace",
104 "ParameterData", Direction::Output),
105 "Output Workspace2D containing fitted peak parameters for "
106 "further refinement.");
107
108 // Zscore table workspace
110 std::make_unique<WorkspaceProperty<TableWorkspace>>("OutputZscoreWorkspace", "ZscoreTable", Direction::Output),
111 "Output TableWorkspace containing the Zscore of the fitted "
112 "peak parameters. ");
113
114 // Workspace index of the
115 declareProperty("WorkspaceIndex", 0, "Worskpace index for the data to refine against.");
116
117 // Range of the peaks to fit
118 declareProperty("MinTOF", EMPTY_DBL(), "Minimum TOF to fit peaks. ");
119 declareProperty("MaxTOF", EMPTY_DBL(), "Maximum TOF to fit peaks. ");
120
121 vector<string> fitmodes(2);
122 fitmodes[0] = "Robust";
123 fitmodes[1] = "Confident";
124 auto fitvalidator = std::make_shared<StringListValidator>(fitmodes);
125 declareProperty("FittingMode", "Robust", fitvalidator,
126 "Fitting mode such that user can determine"
127 "whether the input parameters are trustful or not.");
128
129 // Option to calculate peak position from (HKL) and d-spacing data
130 declareProperty("UseGivenPeakCentreTOF", true,
131 "Use each Bragg peak's centre in TOF given in "
132 "BraggPeakParameterWorkspace."
133 "Otherwise, calculate each peak's centre from d-spacing.");
134
135 vector<string> genpeakoptions{"(HKL) & Calculation", "From Bragg Peak Table"};
136 auto propvalidator = std::make_shared<StringListValidator>(genpeakoptions);
137 declareProperty("PeakParametersStartingValueFrom", "(HKL) & Calculation", propvalidator,
138 "Choice of how to generate starting values of "
139 "Bragg peak profile parmeters.");
140
141 declareProperty("MinimumPeakHeight", 0.20,
142 "Minimum peak height (with background removed) "
143 "Any peak whose maximum height under this value will be "
144 "treated as zero intensity. ");
145
146 // Flag to calculate and trust peak parameters from instrument
147 // declareProperty("ConfidentInInstrumentParameters", false, "Option to
148 // calculate peak parameters from "
149 // "instrument parameters.");
150
151 // Option to denote that peaks are related
152 declareProperty("PeaksCorrelated", false,
153 "Flag for fact that all peaks' corresponding profile parameters "
154 "are correlated by an analytical function");
155
156 // Option for peak's HKL for minimum d-spacing
157 auto arrayprop = std::make_unique<ArrayProperty<int>>("MinimumHKL", "");
158 declareProperty(std::move(arrayprop), "Miller index of the left most peak (peak with "
159 "minimum d-spacing) to be fitted. ");
160
161 // Number of the peaks to fit left to peak with minimum HKL
162 declareProperty("NumberPeaksToFitBelowLowLimit", 0,
163 "Number of peaks to fit with d-spacing value "
164 "less than specified minimum. ");
165
166 // Right most peak property
167 auto righthklprop = std::make_unique<ArrayProperty<int>>("RightMostPeakHKL", "");
168 declareProperty(std::move(righthklprop), "Miller index of the right most peak. "
169 "It is only required and used in RobustFit mode.");
170
171 declareProperty("RightMostPeakLeftBound", EMPTY_DBL(),
172 "Left bound of the right most peak. "
173 "Used in RobustFit mode.");
174
175 declareProperty("RightMostPeakRightBound", EMPTY_DBL(),
176 "Right bound of the right most peak. "
177 "Used in RobustFit mode.");
178
179 // Fit option
180 declareProperty("FitCompositePeakBackground", true,
181 "Flag to do fit to both peak and background in a composite "
182 "function as last fit step.");
183}
184
185//----------------------------------------------------------------------------------------------
189 // 1. Get input
191
192 // 2. Crop input workspace
194
195 // 3. Parse input table workspace
197
198 // 4. Unit cell
199 double latticesize = getParameter("LatticeConstant");
200 if (latticesize == EMPTY_DBL())
201 throw runtime_error("Input instrument table workspace lacks LatticeConstant. "
202 "Unable to complete processing.");
203 m_unitCell.set(latticesize, latticesize, latticesize, 90.0, 90.0, 90.0);
204
205 // 5. Generate peaks
207
208 // 6. Fit peaks & get peak centers
209 m_indexGoodFitPeaks.clear();
210 m_chi2GoodFitPeaks.clear();
211 size_t numpts = m_dataWS->x(m_wsIndex).size();
212 m_peakData.resize(numpts);
213
214 g_log.information() << "[FitPeaks] Total Number of Peak = " << m_vecPeakFunctions.size() << '\n';
215 m_peakFitChi2.resize(m_vecPeakFunctions.size(), -1.0 * DBL_MIN);
216 m_goodFit.resize(m_vecPeakFunctions.size(), false);
217
218 if (m_fitMode == ROBUSTFIT) {
219 g_log.information() << "Fit (Single) Peaks In Robust Mode.\n";
221 } else if (m_fitMode == TRUSTINPUTFIT) {
222 g_log.information() << "Fit Peaks In Trust Mode. Possible to fit "
223 "overlapped "
224 "peaks.\n";
226 } else {
227 g_log.error("Unsupported fit mode!");
228 throw runtime_error("Unsupported fit mode.");
229 }
230
231 // 5. Create Output
232 // a) Create a Table workspace for peak profile
233 pair<TableWorkspace_sptr, TableWorkspace_sptr> tables = genPeakParametersWorkspace();
234 TableWorkspace_sptr outputpeaksws = tables.first;
235 TableWorkspace_sptr ztablews = tables.second;
236 setProperty("OutputBraggPeakParameterWorkspace", outputpeaksws);
237 setProperty("OutputZscoreWorkspace", ztablews);
238
239 // b) Create output data workspace (as a middle stage product)
241 setProperty("OutputWorkspace", outdataws);
242
243 // c) Create data workspace for X0, A, B and S of peak with good fit
245 setProperty("OutputBraggPeakParameterDataWorkspace", peakparamvaluews);
246}
247
248//----------------------------------------------------------------------------------------------
252 // data workspace
253 m_dataWS = this->getProperty("InputWorkspace");
254 m_wsIndex = this->getProperty("WorkspaceIndex");
255 if (m_wsIndex < 0 || m_wsIndex > static_cast<int>(m_dataWS->getNumberHistograms())) {
256 stringstream errss;
257 errss << "Input workspace = " << m_wsIndex << " is out of range [0, " << m_dataWS->getNumberHistograms();
258 g_log.error(errss.str());
259 throw std::invalid_argument(errss.str());
260 }
261
262 // table workspaces for parameters
263 m_peakParamTable = this->getProperty("BraggPeakParameterWorkspace");
264 m_profileTable = this->getProperty("InstrumentParameterWorkspace");
265
266 // fitting range
267 m_tofMin = getProperty("MinTOF");
268 m_tofMax = getProperty("MaxTOF");
269 if (m_tofMin == EMPTY_DBL())
270 m_tofMin = m_dataWS->x(m_wsIndex)[0];
271 if (m_tofMax == EMPTY_DBL())
272 m_tofMax = m_dataWS->x(m_wsIndex).back();
273
274 m_minimumHKL = getProperty("MinimumHKL");
275 m_numPeaksLowerToMin = getProperty("NumberPeaksToFitBelowLowLimit");
276
277 // fitting algorithm option
278 string fitmode = getProperty("FittingMode");
279 if (fitmode == "Robust") {
281 } else if (fitmode == "Confident") {
283 } else {
284 throw runtime_error("Input fit mode can only accept either Robust or Confident. ");
285 }
286
287 m_useGivenTOFh = getProperty("UseGivenPeakCentreTOF");
288 // m_confidentInInstrumentParameters =
289 // getProperty("ConfidentInInstrumentParameters");
290
291 // peak parameter generation option
292 string genpeakparamalg = getProperty("PeakParametersStartingValueFrom");
293 if (genpeakparamalg == "(HKL) & Calculation") {
295 } else if (genpeakparamalg == "From Bragg Peak Table") {
297 } else {
298 throw runtime_error("Input option from PeakParametersStaringValueFrom is not supported.");
299 }
300
301 // Right most peak information
302 m_rightmostPeakHKL = getProperty("RightMostPeakHKL");
303 m_rightmostPeakLeftBound = getProperty("RightMostPeakLeftBound");
304 m_rightmostPeakRightBound = getProperty("RightMostPeakRightBound");
305
306 if (m_fitMode == ROBUSTFIT) {
309 stringstream errss;
310 errss << "If fit mode is 'RobustFit', then user must specify all 3 "
311 "properties of right most peak "
312 << "(1) Miller Index (given size = " << m_rightmostPeakHKL.size() << "), "
313 << "(2) Left boundary (given value = " << m_rightmostPeakLeftBound << "), "
314 << "(3) Right boundary (given value = " << m_rightmostPeakRightBound << "). ";
315 g_log.error(errss.str());
316 throw runtime_error(errss.str());
317 }
318 }
319
320 m_minPeakHeight = getProperty("MinimumPeakHeight");
321
322 m_fitPeakBackgroundComposite = getProperty("FitCompositePeakBackground");
323}
324
325//================================= Fit Peaks In Robust Mode
326//==================================
327
328//----------------------------------------------------------------------------------------------
338 // I. Prepare
340 bool isrightmost = true;
341 size_t numpeaks = m_vecPeakFunctions.size();
342 if (numpeaks == 0)
343 throw runtime_error("There is no peak to fit!");
344
345 vector<string> peakparnames = m_vecPeakFunctions[0].second.second->getParameterNames();
346
347 // II. Create local background function.
348 Polynomial_sptr backgroundfunction = std::make_shared<Polynomial>();
349 backgroundfunction->setAttributeValue("n", 1);
350 backgroundfunction->initialize();
351
352 // III. Fit peaks
353 double chi2;
354 double refpeakshift = 0.0;
355
356 for (int peakindex = static_cast<int>(numpeaks) - 1; peakindex >= 0; --peakindex) {
357 vector<int> peakhkl = m_vecPeakFunctions[peakindex].second.first;
358 BackToBackExponential_sptr thispeak = m_vecPeakFunctions[peakindex].second.second;
359
360 stringstream infoss;
361
362 bool goodfit = false;
363
364 if (isrightmost && peakhkl == m_rightmostPeakHKL) {
365 // It is the specified right most peak. Estimate background, peak height,
366 // fwhm, ...
367 // 1. Determine the starting value of the peak
368 double peakleftbound, peakrightbound;
369 peakleftbound = m_rightmostPeakLeftBound;
370 peakrightbound = m_rightmostPeakRightBound;
371
372 double predictpeakcentre = thispeak->centre();
373
374 infoss << "[DBx102] The " << numpeaks - 1 - peakindex << "-th rightmost peak's miller index = " << peakhkl[0]
375 << ", " << peakhkl[1] << ", " << peakhkl[2] << ", predicted at TOF = " << thispeak->centre()
376 << "; User specify boundary = [" << peakleftbound << ", " << peakrightbound << "].";
377 g_log.information() << infoss.str() << '\n';
378
379 map<string, double> rightpeakparameters;
380 goodfit = fitSinglePeakRobust(thispeak, std::dynamic_pointer_cast<BackgroundFunction>(backgroundfunction),
381 peakleftbound, peakrightbound, rightpeakparameters, chi2);
382
383 m_peakFitChi2[peakindex] = chi2;
384
385 if (!goodfit)
386 throw runtime_error("Failed to fit the right most peak. Unable to process. ");
387
388 stringstream robmsgss;
389 for (auto &parname : peakparnames) {
390 robmsgss << parname << " = " << thispeak->getParameter(parname) << '\n';
391 }
392 g_log.information() << "[DB1151] Robust Fit Result: Chi^2 = " << chi2 << '\n' << robmsgss.str();
393
394 rightpeak = thispeak;
395 isrightmost = false;
396
397 // iii. Reference peak shift
398 refpeakshift = thispeak->centre() - predictpeakcentre;
399
400 g_log.notice() << "[DBx332] Peak -" << static_cast<int>(numpeaks) - peakindex - 1
401 << ": shifted = " << refpeakshift << '\n';
402 } else if (!isrightmost) {
403 // All peaks but not the right most peak
404 // 1. Validate inputs
405 if (peakindex == static_cast<int>(numpeaks) - 1)
406 throw runtime_error("Impossible to have peak index as the right most peak here!");
407
408 double predictcentre = thispeak->centre();
409
410 // 2. Determine the peak range by observation
411 double peakleftbound, peakrightbound;
412 observePeakRange(thispeak, rightpeak, refpeakshift, peakleftbound, peakrightbound);
413
414 stringstream dbxss;
415 dbxss << '\n';
416 for (int i = 0; i < 10; ++i)
417 dbxss << "==\n[DBx323] Peak (" << peakhkl[0] << ", " << peakhkl[1] << "," << peakhkl[2]
418 << "). Centre predicted @ TOF = " << predictcentre << ". Observed range = " << peakleftbound << ", "
419 << peakrightbound;
420 g_log.notice(dbxss.str());
421
422 // 3. Fit peak
423 /* Disabled and replaced by fit1PeakRobust
424 goodfit = fitSinglePeakRefRight(thispeak, backgroundfunction, rightpeak,
425 peakleftbound,
426 peakrightbound, chi2);
427 */
428
429 map<string, double> rightpeakparameters;
430 storeFunctionParameters(rightpeak, rightpeakparameters);
431 goodfit = fitSinglePeakRobust(thispeak, std::dynamic_pointer_cast<BackgroundFunction>(backgroundfunction),
432 peakleftbound, peakrightbound, rightpeakparameters, chi2);
433
434 if (goodfit) {
435 // Fit successful
436 m_peakFitChi2[peakindex] = chi2;
437 // Update right peak and reference peak shift if peak is not trivial
438 if (thispeak->height() >= m_minPeakHeight) {
439 rightpeak = thispeak;
440 refpeakshift = thispeak->centre() - predictcentre;
441 }
442
443 } else {
444 // Bad fit
445 m_peakFitChi2[peakindex] = -1.0;
446 g_log.warning() << "Fitting peak @ " << thispeak->centre() << " failed. \n";
447 }
448
449 } else {
450 // It is right to the specified right most peak. Skip to next peak
451 double peakleftbound, peakrightbound;
452 peakleftbound = m_rightmostPeakLeftBound;
453 peakrightbound = m_rightmostPeakRightBound;
454
455 infoss << "[DBx102] The " << numpeaks - 1 - peakindex << "-th rightmost peak's miller index = " << peakhkl[0]
456 << ", " << peakhkl[1] << ", " << peakhkl[2] << ", predicted at TOF = " << thispeak->centre() << "; "
457 << "User specify right most peak's miller index = " << m_rightmostPeakHKL[0] << ", "
458 << m_rightmostPeakHKL[1] << ", " << m_rightmostPeakHKL[2] << " User specify boundary = [" << peakleftbound
459 << ", " << peakrightbound << "].";
460 g_log.information() << infoss.str() << '\n';
461 continue;
462 }
463
464 } // ENDFOR Peaks
465}
466
467//----------------------------------------------------------------------------------------------
472 const BackToBackExponential_sptr &rightpeak, double refpeakshift,
473 double &peakleftbound, double &peakrightbound) {
474 double predictcentre = thispeak->centre();
475 double rightfwhm = rightpeak->fwhm();
476
477 // 1. Roughly determine the peak range from this peak's starting values and
478 // right peak' fitted
479 // parameters values
480 if (refpeakshift > 0) {
481 // tend to shift to right
482 peakleftbound = predictcentre - rightfwhm;
483 peakrightbound = predictcentre + rightfwhm + refpeakshift;
484 } else {
485 // tendency to shift to left
486 peakleftbound = predictcentre - rightfwhm + refpeakshift;
487 peakrightbound = predictcentre + rightfwhm;
488 }
489 if (peakrightbound > rightpeak->centre() - 3 * rightpeak->fwhm()) {
490 // the search of peak's right end shouldn't exceed the left tail of its real
491 // right peak!
492 // Remember this is robust mode. Any 2 adjacent peaks should be faw enough.
493 peakrightbound = rightpeak->centre() - 3 * rightpeak->fwhm();
494 }
495
496 // 2. Search for maximum
497 const auto &vecX = m_dataWS->x(m_wsIndex);
498
499 size_t icentre = findMaxValue(m_dataWS, m_wsIndex, peakleftbound, peakrightbound);
500
501 // 3. Narrow now the peak range
502 peakleftbound = vecX[icentre] - 4.0 * rightfwhm;
503 peakrightbound = vecX[icentre] + 4.0 * rightfwhm;
504
505 double rightpeakleftbound = rightpeak->centre() - 3 * rightfwhm;
506 if (peakrightbound > rightpeakleftbound) {
507 peakrightbound = rightpeakleftbound;
508 double peakcentre = vecX[icentre];
509 if (peakrightbound < 2.0 * rightfwhm + peakcentre)
510 g_log.warning() << "Peak @ " << peakcentre << "'s right boundary is too close to its right peak!\n";
511 }
512}
513
514//----------------------------------------------------------------------------------------------
554 const BackgroundFunction_sptr &backgroundfunction, double peakleftbound,
555 double peakrightbound, const map<string, double> &rightpeakparammap,
556 double &finalchi2) {
557 // 1. Build partial workspace
558 Workspace2D_sptr peakws = buildPartialWorkspace(m_dataWS, m_wsIndex, peakleftbound, peakrightbound);
559 g_log.debug() << "[DB1208] Build partial workspace for peak @ " << peak->centre() << " (predicted).\n";
560
561 // 2. Estimate and remove background
562 size_t rawdata_wsindex = 0;
563 size_t estbkgd_wsindex = 2;
564 size_t peak_wsindex = 1;
565 estimateBackgroundCoarse(peakws, backgroundfunction, rawdata_wsindex, estbkgd_wsindex, peak_wsindex);
566
567 stringstream dbss;
568 dbss << "[DBx203] Removed background peak data: \n";
569 for (size_t i = 0; i < peakws->x(peak_wsindex).size(); ++i)
570 dbss << peakws->x(peak_wsindex)[i] << "\t\t" << peakws->y(peak_wsindex)[i] << "\t\t" << peakws->e(peak_wsindex)[i]
571 << '\n';
572 g_log.debug(dbss.str());
573
574 // 3. Estimate FWHM, peak centre, and height
575 double centre, fwhm, height;
576 string errmsg;
577 bool pass = observePeakParameters(peakws, 1, centre, height, fwhm, errmsg);
578 if (!pass) {
579// If estiamtion fails, quit b/c first/rightmost peak must be fitted.
580#if 0
581 g_log.error(errmsg);
582 throw runtime_error(errmsg);
583#else
584 g_log.notice("Unable to observe peak parameters. Proceed to next peak.");
585 return false;
586#endif
587 } else if (height < m_minPeakHeight) {
588 g_log.notice() << "[FLAGx409] Peak proposed @ TOF = " << peak->centre() << " has a trivial "
589 << "peak height = " << height << " by observation. Skipped.\n";
590 return false;
591 } else {
592 g_log.information() << "[DBx201] Peak Predicted @ TOF = " << peak->centre()
593 << ", Estimated (observation) Centre = " << centre << ", FWHM = " << fwhm
594 << " Height = " << height << '\n';
595 }
596
597 // 4. Fit by Gaussian to get some starting value
598 double tof_h, sigma;
599 doFitGaussianPeak(peakws, peak_wsindex, centre, fwhm, fwhm, tof_h, sigma, height);
600
601 // 5. Fit by various methods
602 // Set all parameters for fit
603 vector<string> peakparnames = peak->getParameterNames();
604 for (size_t i = 0; i < peakparnames.size(); ++i)
605 peak->unfix(i);
606
607 // Set up the universal starting parameter
608 peak->setParameter("I", height * fwhm);
609 peak->setParameter("X0", tof_h);
610
611 size_t numsteps = 2;
612 vector<string> minimizers(numsteps);
613 minimizers[0] = "Simplex";
614 minimizers[1] = "Levenberg-Marquardt";
615 vector<size_t> maxiterations(numsteps, 10000);
616 vector<double> dampfactors(numsteps, 0.0);
617
618 // Record the start value
619 map<string, double> origparammap;
620 storeFunctionParameters(peak, origparammap);
621
622 vector<double> chi2s;
623 vector<bool> goodfits;
624 vector<map<string, double>> solutions;
625
626 // a) Fit by using input peak parameters
627 string peakinfoa0 = getFunctionInfo(std::dynamic_pointer_cast<IFunction>(peak));
628 g_log.notice() << "[DBx533A] Approach A: Starting Peak Function Information: \n" << peakinfoa0 << '\n';
629
630 double chi2a;
631 bool fitgooda = doFit1PeakSequential(peakws, peak_wsindex, peak, minimizers, maxiterations, dampfactors, chi2a);
632 map<string, double> solutiona;
633 storeFunctionParameters(peak, solutiona);
634
635 chi2s.emplace_back(chi2a);
636 goodfits.emplace_back(fitgooda);
637 solutions.emplace_back(solutiona);
638
639 string peakinfoa1 = getFunctionInfo(std::dynamic_pointer_cast<IFunction>(peak));
640 g_log.notice() << "[DBx533A] Approach A: Fit Successful = " << fitgooda << ", Chi2 = " << chi2a
641 << ", Peak Function Information: \n"
642 << peakinfoa1 << '\n';
643
644 // b) Fit by using Gaussian result (Sigma)
645 restoreFunctionParameters(peak, origparammap);
646 peak->setParameter("S", sigma);
647
648 string peakinfob0 = getFunctionInfo(std::dynamic_pointer_cast<IFunction>(peak));
649 g_log.notice() << "[DBx533B] Approach B: Starting Peak Function Information: \n" << peakinfob0 << '\n';
650
651 double chi2b;
652 bool fitgoodb = doFit1PeakSequential(peakws, peak_wsindex, peak, minimizers, maxiterations, dampfactors, chi2b);
653
654 map<string, double> solutionb;
655 storeFunctionParameters(peak, solutionb);
656
657 chi2s.emplace_back(chi2b);
658 goodfits.emplace_back(fitgoodb);
659 solutions.emplace_back(solutionb);
660
661 string peakinfob1 = getFunctionInfo(std::dynamic_pointer_cast<IFunction>(peak));
662 g_log.notice() << "[DBx533B] Approach 2: Fit Successful = " << fitgoodb << ", Chi2 = " << chi2b
663 << ", Peak Function Information: \n"
664 << peakinfob1 << '\n';
665
666 // c) Fit peak parameters by the value from right peak
667 if (!rightpeakparammap.empty()) {
668 restoreFunctionParameters(peak, rightpeakparammap);
669 peak->setParameter("X0", tof_h);
670 peak->setParameter("I", height * fwhm);
671
672 string peakinfoc0 = getFunctionInfo(std::dynamic_pointer_cast<IFunction>(peak));
673 g_log.notice() << "[DBx533C] Approach C: Starting Peak Function Information: \n" << peakinfoc0 << '\n';
674
675 double chi2c;
676 bool fitgoodc = doFit1PeakSequential(peakws, peak_wsindex, peak, minimizers, maxiterations, dampfactors, chi2c);
677 map<string, double> solutionc;
678 storeFunctionParameters(peak, solutionc);
679
680 chi2s.emplace_back(chi2c);
681 goodfits.emplace_back(fitgoodc);
682 solutions.emplace_back(solutionc);
683
684 string peakinfoc1 = getFunctionInfo(std::dynamic_pointer_cast<IFunction>(peak));
685 g_log.notice() << "[DBx533C] Approach C: Fit Successful = " << fitgoodc << ", Chi2 = " << chi2c
686 << ", Peak Function Information: \n"
687 << peakinfoc1 << '\n';
688 } else {
689 // No right peak information: set a error entry
690 chi2s.emplace_back(DBL_MAX);
691 goodfits.emplace_back(false);
692 solutions.emplace_back(rightpeakparammap);
693 }
694
695 // 6. Summarize the above 3 approach
696 size_t bestapproach = goodfits.size() + 1;
697 double bestchi2 = DBL_MAX;
698 for (size_t i = 0; i < goodfits.size(); ++i) {
699 if (goodfits[i] && chi2s[i] < bestchi2) {
700 bestapproach = i;
701 bestchi2 = chi2s[i];
702 }
703 }
704
705 stringstream fitsumss;
706 fitsumss << "Best fit result is obtained by approach " << bestapproach << " of total " << goodfits.size()
707 << " approaches. Best Chi^2 = " << bestchi2 << ", Peak Height = " << peak->height();
708 g_log.notice() << "[DB1127] " << fitsumss.str() << '\n';
709
710 bool fitgood = true;
711 if (bestapproach < goodfits.size()) {
712 restoreFunctionParameters(peak, solutions[bestapproach]);
713 } else {
714 fitgood = false;
715 }
716
717 // 7. Fit by Monte Carlo if previous failed
718 if (!fitgood) {
719 peak->setParameter("S", sigma);
720 peak->setParameter("I", height * fwhm);
721 peak->setParameter("X0", tof_h);
722
723 vector<string> paramsinmc{"A", "B"};
724
725 fitSinglePeakSimulatedAnnealing(peak, paramsinmc);
726 }
727
728 // 8. Fit with background
730 // Fit peak + background
731 double chi2compf;
732 bool fitcompfunsuccess = doFit1PeakBackground(peakws, rawdata_wsindex, peak, backgroundfunction, chi2compf);
733 if (fitcompfunsuccess) {
734 finalchi2 = chi2compf;
735 } else {
736 finalchi2 = bestchi2;
737 g_log.warning("Fit peak-background composite function failed! Need to "
738 "find out how this case peak value is changed from best "
739 "fit.");
740 }
741 } else {
742 // Flag is turned off
743 finalchi2 = bestchi2;
744 }
745
746 // 9. Plot function
747 FunctionDomain1DVector domain(peakws->x(0).rawData());
748 plotFunction(peak, backgroundfunction, domain);
749
750 return fitgood;
751}
752
753//----------------------------------------------------------------------------------------------
758 const BackToBackExponential_sptr &peak,
759 const BackgroundFunction_sptr &backgroundfunction, double &chi2) {
760 // 0. Set fit parameters
761 string minimzername("Levenberg-MarquardtMD");
762 double startx = peak->centre() - 3.0 * peak->fwhm();
763 double endx = peak->centre() + 3.0 * peak->fwhm();
764
765 // 1. Create composite function
767 compfunc->addFunction(peak);
768 compfunc->addFunction(backgroundfunction);
769
770 // 2. Unfix all parameters
771 vector<string> comparnames = compfunc->getParameterNames();
772 for (size_t ipar = 0; ipar < comparnames.size(); ++ipar)
773 compfunc->unfix(ipar);
774
775 // 3. Fit
776 string cominfoa = getFunctionInfo(compfunc);
777 g_log.notice() << "[DBx533X-0] Fit All: Starting Peak Function Information: \n"
778 << cominfoa << "Fit range = " << startx << ", " << endx << '\n';
779
780 // 3. Set
781 auto fitalg = createChildAlgorithm("Fit", -1, -1, true);
782 fitalg->initialize();
783
784 fitalg->setProperty("Function", std::dynamic_pointer_cast<API::IFunction>(compfunc));
785 fitalg->setProperty("InputWorkspace", dataws);
786 fitalg->setProperty("WorkspaceIndex", static_cast<int>(wsindex));
787 fitalg->setProperty("Minimizer", minimzername);
788 fitalg->setProperty("CostFunction", "Least squares");
789 fitalg->setProperty("MaxIterations", 1000);
790 fitalg->setProperty("Output", "FitPeakBackground");
791 fitalg->setProperty("StartX", startx);
792 fitalg->setProperty("EndX", endx);
793
794 // 3. Execute and parse the result
795 bool isexecute = fitalg->execute();
796 bool fitsuccess;
797 chi2 = DBL_MAX;
798
799 if (isexecute) {
800 std::string fitresult = parseFitResult(fitalg, chi2, fitsuccess);
801
802 // Figure out result
803 stringstream cominfob;
804 cominfob << "[DBx533X] Fit All: Fit Successful = " << fitsuccess << ", Chi^2 = " << chi2 << '\n';
805 cominfob << "Detailed info = " << fitresult << '\n';
806 string fitinfo = getFunctionInfo(compfunc);
807 cominfob << fitinfo;
808
809 g_log.notice(cominfob.str());
810 } else {
811 g_log.notice() << "[DB1203B] Failed To Fit Peak+Background @ " << peak->centre() << '\n';
812 fitsuccess = false;
813 }
814
815 return fitsuccess;
816}
817
818//----------------------------------------------------------------------------------------------
822 const vector<string> &paramtodomc) {
823 UNUSED_ARG(peak);
824 UNUSED_ARG(paramtodomc);
825 throw runtime_error("To Be Implemented Soon!");
826
827 /*
828 // 6. Fit peak by the result from Gaussian
829 pair<bool, double> fitresult;
830
831 double varA = 1.0;
832 double varB = 1.0;
833
834
835 // a) Get initial chi2
836 double startchi2 = DBL_MAX;
837 doFit1PeakSimple(peakws, 1, peak, "Levenberg-MarquardtMD", 0, newchi2);
838 g_log.debug() << "[DBx401] Peak @ TOF = " << peak->centre() << ", Starting
839 Chi^2 = "
840 << newchi2 << '\n';
841
842 bool continuetofit = true;
843 bool goodfit;
844
845 int count = 0;
846 size_t numfittohave = 5;
847 int maxnumloops = 100;
848 vector<pair<double, map<string, double> > > fitparammaps; // <chi2,
849 <parameter, value> >
850
851 map<string, double> curparammap;
852 double curchi2 = startchi2;
853
854 srand(0);
855 size_t reject = 0;
856 size_t count10 = 0;
857 double stepratio = 1.0;
858 while(continuetofit)
859 {
860 // a) Store
861 storeFunctionParameters(peak, curparammap);
862
863 // b) Monte Carlo
864 double aorb = static_cast<double>(rand())/static_cast<double>(RAND_MAX)-0.5;
865 double change =
866 static_cast<double>(rand())/static_cast<double>(RAND_MAX)-0.5;
867 if (aorb > 0)
868 {
869 // A
870 varA = varA*(1+stepratio*change);
871 if (varA <= DBL_MIN)
872 varA = fabs(varA) + DBL_MIN;
873 }
874 else
875 {
876 // B
877 varB = varB*(1+stepratio*change);
878 if (varB <= DBL_MIN)
879 varB = fabs(varB) + DBL_MIN;
880 }
881
882 // b) Set A and B
883 peak->setParameter("A", varA);
884 peak->setParameter("B", varB);
885
886 // b) Fit
887 g_log.debug() << "DBx329: Proposed A = " << varA << ", B = " << varB <<
888 '\n';
889 fitresult = doFitPeak(peakws, peak, fwhm);
890
891 // c) Get fit result
892 goodfit = fitresult.first;
893 newchi2 = fitresult.second;
894
895 // d) Take it or not!
896 ++ count10;
897 if (!goodfit)
898 {
899 // Bad fit, reverse
900 restoreFunctionParameters(peak, curparammap);
901 ++ reject;
902 }
903 else if (newchi2 < curchi2)
904 {
905 // Good. Take it
906 curchi2 = newchi2;
907 }
908 else
909 {
910 // Toss the dice again
911 double bar = exp(-(newchi2-curchi2)/curchi2);
912 double dice = static_cast<double>(rand())/static_cast<double>(RAND_MAX);
913 if (dice < bar)
914 {
915 // Take it
916 curchi2 = newchi2;
917 }
918 else
919 {
920 // Reject it
921 restoreFunctionParameters(peak, curparammap);
922 ++ reject;
923 }
924 }
925
926 // Adjust step size
927 if (count10 >=10 && reject >= 8)
928 {
929 reject = 10;
930 count10 = 0;
931 stepratio *= 5;
932 }
933
934 // e) Record
935 if (goodfit && newchi2 < startchi2)
936 {
937 // i. Store parameters;
938 map<string,double> parammap;
939 for (size_t i = 0; i < peakparnames.size(); ++i)
940 parammap.emplace(peakparnames[i],
941 peak->getParameter(peakparnames[i]));
942 fitparammaps.emplace_back(newchi2, parammap);
943
944 // ii. sort
945 sort(fitparammaps.begin(), fitparammaps.end());
946
947 // iii. keep size
948 if (fitparammaps.size() > numfittohave)
949 fitparammaps.pop_back();
950 }
951
952 // f) Loop control & MC
953 count += 1;
954 if (count >= maxnumloops)
955 {
956 continuetofit = false;
957 }
958 } // ENDWHILE For Fit
959
960 // 7. Process MC fit results
961 if (fitparammaps.size() > 0)
962 {
963 // Good results
964 goodfit = true;
965 sort(fitparammaps.begin(), fitparammaps.end());
966 newchi2 = fitparammaps[0].first;
967
968 // Set the best result to peak
969 for (size_t i = 0; i < peakparnames.size(); ++i)
970 {
971 string& parname = peakparnames[i];
972 map<string, double>& thisparams = fitparammaps[0].second;
973 map<string, double>::iterator liter = thisparams.find(parname);
974 if (liter == thisparams.end())
975 {
976 stringstream errmsg;
977 errmsg << "In parameter value map, there is no key named " << parname <<
978 '\n';
979 throw runtime_error(errmsg.str());
980 }
981 peak->setParameter(parname, liter->second);
982 }
983
984 // Plot the peak
985 FunctionDomain1DVector domain(peakws->x(1));
986 plotFunction(peak, backgroundfunction, domain);
987
988 // Debug print the best solutions
989 for (size_t i = 0; i < fitparammaps.size(); ++i)
990 {
991 double thischi2 = fitparammaps[i].first;
992 map<string, double>& thisparams = fitparammaps[i].second;
993
994 g_log.information() << "[DB1055] Chi2 = " << thischi2 << ", A = " <<
995 thisparams["A"]
996 << ", B = " << thisparams["B"] << '\n';
997 }
998 }
999 else
1000 {
1001 // Bad result
1002 goodfit = false;
1003 stringstream errmsg;
1004 errmsg << "Failed to fit peak @ " << tof_h << " in robust mode!\n";
1005 g_log.error(errmsg.str());
1006 }
1007
1008 return goodfit;
1009 */
1010
1011 return false;
1012}
1013
1014//============================== Fit Peaks With Good Starting Values
1015//==========================
1016
1017//----------------------------------------------------------------------------------------------
1026 // 1. Initialize (local) background function
1027 Polynomial_sptr backgroundfunction = std::make_shared<Polynomial>();
1028 backgroundfunction->setAttributeValue("n", 1);
1029 backgroundfunction->initialize();
1030
1031 // 2. Fit peak / peaks group
1032 int ipeak = static_cast<int>(m_vecPeakFunctions.size()) - 1;
1033 double chi2;
1034
1035 // BackToBackExponential_sptr rightpeak = m_peaks[ipeak].second.second;
1036
1037 while (ipeak >= 0) {
1038 // 1. Make a peak group
1039 vector<size_t> indexpeakgroup;
1040
1041 bool makegroup = true;
1042 // Loop over 2nd level: peak groups: 1 peak or overlapped peaks
1043 while (makegroup) {
1044 // There is no need to worry about its right neighbor, b/c
1045 // this situation is already considered as its right neighbor is
1046 // treated;
1047
1048 // a) Add this peak,
1049 BackToBackExponential_sptr thispeak = m_vecPeakFunctions[ipeak].second.second;
1050 indexpeakgroup.emplace_back(ipeak);
1051
1052 // b) update the peak index
1053 --ipeak;
1054
1055 if (ipeak < 0) {
1056 // this is last peak. next peak does not exist
1057 makegroup = false;
1058 } else {
1059 // this is not the last peak. search the left one.
1060 double thispeakleftbound = thispeak->centre() - thispeak->fwhm() * 2.5;
1061 BackToBackExponential_sptr leftpeak = m_vecPeakFunctions[ipeak].second.second;
1062 double leftpeakrightbound = leftpeak->centre() + leftpeak->fwhm() * 2.5;
1063 if (thispeakleftbound > leftpeakrightbound) {
1064 // This peak and next peak is faw enough!
1065 makegroup = false;
1066 }
1067 }
1068 } // while
1069
1070 if (indexpeakgroup.size() == 1) {
1071 // Fit a single peak
1072 const size_t ifit = indexpeakgroup[0];
1073 double peakfitleftbound, peakfitrightbound;
1074 calculatePeakFitBoundary(ifit, ifit, peakfitleftbound, peakfitrightbound);
1075
1076 g_log.information() << "\n[T] Fit Peak Indexed " << ifit << " (" << m_vecPeakFunctions.size() - 1 - ifit
1077 << ")\t----------------------------------\n";
1078
1079 BackToBackExponential_sptr thispeak = m_vecPeakFunctions[ifit].second.second;
1080 bool annihilatedpeak;
1081 m_goodFit[ifit] = fitSinglePeakConfident(thispeak, backgroundfunction, peakfitleftbound, peakfitrightbound, chi2,
1082 annihilatedpeak);
1083 m_peakFitChi2[ifit] = chi2;
1084 if (annihilatedpeak)
1085 thispeak->setHeight(0.0);
1086
1087 // Debug output
1088 const vector<int> &hkl = m_vecPeakFunctions[ifit].second.first;
1089 stringstream dbss;
1090 dbss << "Peak [" << hkl[0] << ", " << hkl[1] << ", " << hkl[2] << "] expected @ TOF = " << thispeak->centre()
1091 << ": \t";
1092 if (annihilatedpeak)
1093 dbss << "Annihilated!";
1094 else
1095 dbss << "Fit Status = " << m_goodFit[ifit] << ", Chi2 = " << chi2;
1096 g_log.information() << "[DB531] " << dbss.str() << '\n';
1097 } else {
1098 // Fit overlapped peaks
1099 vector<BackToBackExponential_sptr> peaksgroup;
1100 for (auto ipk : indexpeakgroup) {
1101 BackToBackExponential_sptr temppeak = m_vecPeakFunctions[ipk].second.second;
1102 peaksgroup.emplace_back(temppeak);
1103 }
1104
1105 fitOverlappedPeaks(peaksgroup, backgroundfunction, -1.0);
1106 // rightpeak = indexpeaksgroup[0];
1107 }
1108 } // ENDWHILE
1109
1110 // 2. Output
1111
1112 g_log.information() << "DBx415: Number of good fit peaks = " << m_indexGoodFitPeaks.size() << '\n';
1113
1114 // 3. Clean up
1115 g_log.information() << "[FitPeaks] Number of peak of good chi2 = " << m_chi2GoodFitPeaks.size() << '\n';
1116}
1117
1118//----------------------------------------------------------------------------------------------
1128 const BackgroundFunction_sptr &backgroundfunction, double leftbound,
1129 double rightbound, double &chi2, bool &annhilatedpeak) {
1130 // 1. Build the partial workspace
1131 // a) Determine boundary if necessary
1132 if (leftbound < 0 || leftbound >= peak->centre())
1133 leftbound = peak->centre() - 5 * peak->fwhm();
1134 if (rightbound < 0 || rightbound <= peak->centre())
1135 rightbound = peak->centre() + 5 * peak->fwhm();
1136
1137 // b) Build partial
1138 Workspace2D_sptr peakdataws = buildPartialWorkspace(m_dataWS, m_wsIndex, leftbound, rightbound);
1139
1140 // 2. Remove background
1141 estimateBackgroundCoarse(peakdataws, backgroundfunction, 0, 2, 1);
1142
1143 stringstream dbss2;
1144 dbss2 << "[T] Partial workspace No Background: \n";
1145 for (size_t i = 0; i < peakdataws->x(1).size(); ++i)
1146 dbss2 << peakdataws->x(1)[i] << "\t\t" << peakdataws->y(1)[i] << "\t\t" << peakdataws->e(1)[i] << "\t\t"
1147 << peakdataws->y(0)[i] << '\n';
1148 g_log.notice(dbss2.str());
1149
1150 // 3. Estimate peak heights
1151 size_t imaxheight = findMaxValue(peakdataws->y(1).rawData());
1152 double maxheight = peakdataws->y(1)[imaxheight];
1153 if (maxheight <= m_minPeakHeight) {
1154 // Max height / peak height is smaller than user defined minimum height. No
1155 // fit, Zero
1156 annhilatedpeak = true;
1157 return false;
1158 } else {
1159 // Max hieght is larger than user defined minimum. Fit it
1160 annhilatedpeak = false;
1161 }
1162
1163 // 4. Set the constraint and height
1164 // a) Peak centre
1165 double peakcentreleftbound = peak->centre() - peak->fwhm();
1166 double peakcentrerightbound = peak->centre() + peak->fwhm();
1167 auto x0bc = std::make_unique<BoundaryConstraint>(peak.get(), "X0", peakcentreleftbound, peakcentrerightbound);
1168 peak->addConstraint(std::move(x0bc));
1169
1170 // b) A
1171 auto abc = std::make_unique<BoundaryConstraint>(peak.get(), "A", 1.0E-10, false);
1172 peak->addConstraint(std::move(abc));
1173
1174 // c) B
1175 auto bbc = std::make_unique<BoundaryConstraint>(peak.get(), "B", 1.0E-10, false);
1176 peak->addConstraint(std::move(bbc));
1177
1178 // d) Guessed height
1179 peak->setHeight(maxheight);
1180
1181 /* Debug information */
1182 stringstream dbss0;
1183 dbss0 << "[DBx100] Peak @" << peak->centre() << ", FWHM = " << peak->fwhm() << '\n';
1184 vector<string> peakparams = peak->getParameterNames();
1185 for (size_t i = 0; i < peakparams.size(); ++i)
1186 dbss0 << peakparams[i] << " = " << peak->getParameter(i) << '\n';
1187 g_log.notice(dbss0.str());
1188
1189 // 5. Fit peak with simple scheme
1190 vector<string> peakparamnames = peak->getParameterNames();
1191 vector<map<string, double>> fitparamvaluemaps;
1192 vector<pair<double, size_t>> chi2indexvec;
1193
1194 // a) Fit peak height
1195 for (size_t iparam = 0; iparam < peakparamnames.size(); ++iparam) {
1196 const string &parname = peakparams[iparam];
1197 if (parname == "I")
1198 peak->unfix(iparam);
1199 else
1200 peak->fix(iparam);
1201 }
1202 double chi2height;
1203 bool goodfit1 = doFit1PeakSimple(peakdataws, 1, peak, "Levenberg-MarquardtMD", 10000, chi2height);
1204
1205 // Store parameters
1206 map<string, double> step1params;
1207 storeFunctionParameters(peak, step1params);
1208 fitparamvaluemaps.emplace_back(step1params);
1209 if (!goodfit1)
1210 chi2height = 1.0E20;
1211 chi2indexvec.emplace_back(chi2height, 0);
1212
1213 // Fix background
1214 vector<string> bkgdparnames = backgroundfunction->getParameterNames();
1215 for (size_t iname = 0; iname < bkgdparnames.size(); ++iname)
1216 backgroundfunction->fix(iname);
1217
1218 // b) Plan A: fit all parameters
1219 double chi2planA;
1220 for (size_t iparam = 0; iparam < peakparamnames.size(); ++iparam)
1221 peak->unfix(iparam);
1222 bool goodfitA = doFit1PeakSimple(peakdataws, 1, peak, "Simplex", 10000, chi2planA);
1223
1224 // Store A's result
1225 map<string, double> planAparams;
1226 storeFunctionParameters(peak, planAparams);
1227 if (!goodfitA)
1228 chi2planA = 1.0E20;
1229 fitparamvaluemaps.emplace_back(planAparams);
1230 chi2indexvec.emplace_back(chi2planA, 1);
1231
1232 // c) Plan B: fit parameters in two groups in 2 steps
1233 // i. Restore step 1's result
1234 restoreFunctionParameters(peak, step1params);
1235
1236 // ii. Fit peak height and everything else but "A"
1237 double chi2planB;
1238 for (size_t iparam = 0; iparam < peakparamnames.size(); ++iparam) {
1239 string parname = peakparams[iparam];
1240 if (parname == "A")
1241 peak->fix(iparam);
1242 else
1243 peak->unfix(iparam);
1244 }
1245 bool goodfitB1 = doFit1PeakSimple(peakdataws, 1, peak, "Levenberg-MarquardtMD", 10000, chi2planB);
1246
1247 // iii. Fit "A" only
1248 for (size_t iparam = 0; iparam < peakparamnames.size(); ++iparam) {
1249 string parname = peakparams[iparam];
1250 if (parname == "A" || parname == "I")
1251 peak->unfix(iparam);
1252 else
1253 peak->fix(iparam);
1254 }
1255 bool goodfitB2 = doFit1PeakSimple(peakdataws, 1, peak, "Levenberg-MarquardtMD", 10000, chi2planB);
1256
1257 bool goodfitB = goodfitB1 || goodfitB2;
1258 map<string, double> planBparams;
1259 storeFunctionParameters(peak, planBparams);
1260 if (!goodfitB)
1261 chi2planB = 1.0E20;
1262 fitparamvaluemaps.emplace_back(planBparams);
1263 chi2indexvec.emplace_back(chi2planB, 2);
1264
1265 // d) Plan C: fit parameters in two groups in 2 steps in alternate order
1266 // i. Restore step 1's result
1267 restoreFunctionParameters(peak, step1params);
1268
1269 // ii. Fit "A"
1270 double chi2planC;
1271 for (size_t iparam = 0; iparam < peakparamnames.size(); ++iparam) {
1272 string parname = peakparams[iparam];
1273 if (parname != "A")
1274 peak->fix(iparam);
1275 else
1276 peak->unfix(iparam);
1277 }
1278 bool goodfitC1 = doFit1PeakSimple(peakdataws, 1, peak, "Levenberg-MarquardtMD", 10000, chi2planC);
1279
1280 // iii. Fit peak height and everything else but "A"
1281 for (size_t iparam = 0; iparam < peakparamnames.size(); ++iparam) {
1282 string parname = peakparams[iparam];
1283 if (parname == "A")
1284 peak->fix(iparam);
1285 else
1286 peak->unfix(iparam);
1287 }
1288 bool goodfitC2 = doFit1PeakSimple(peakdataws, 1, peak, "Levenberg-MarquardtMD", 10000, chi2planB);
1289
1290 bool goodfitC = goodfitC1 || goodfitC2;
1291 map<string, double> planCparams;
1292 storeFunctionParameters(peak, planCparams);
1293 if (!goodfitC)
1294 chi2planC = 1.0E20;
1295 fitparamvaluemaps.emplace_back(planCparams);
1296 chi2indexvec.emplace_back(chi2planC, 3);
1297
1298 // d) Summarize and compare result
1299 stringstream sumss;
1300 sumss << "[DBx833] Confident fit on peak summary: \n";
1301 for (size_t i = 0; i < 4; ++i)
1302 sumss << "Index " << i << ": chi^2 = " << chi2indexvec[i].first << '\n';
1303 g_log.notice(sumss.str());
1304
1305 sort(chi2indexvec.begin(), chi2indexvec.end());
1306
1307 bool goodfit;
1308 if (chi2indexvec[0].first < 1.0E19) {
1309 // There is good fit.
1310 size_t goodindex = chi2indexvec[0].second;
1311 restoreFunctionParameters(peak, fitparamvaluemaps[goodindex]);
1312 chi2 = chi2indexvec[0].first;
1313 goodfit = true;
1314 } else {
1315 // There is no good fit
1316 chi2 = DBL_MAX;
1317 goodfit = false;
1318 }
1319
1320 // 6. Plot the peak in the output workspace data
1321 if (goodfit) {
1322 FunctionDomain1DVector domain(peakdataws->x(1).rawData());
1323 plotFunction(peak, backgroundfunction, domain);
1324 } else {
1325 // Throw exception if fit peak bad. This is NOT a PERMANANT solution.
1326 stringstream errss;
1327 errss << "Fit Peak @ " << peak->centre() << " Error! Chi^2 (false) = " << chi2
1328 << ". Do Not Know How To Proceed To Next Peak!";
1329 g_log.error(errss.str());
1330 throw runtime_error(errss.str());
1331 }
1332
1333 // 7. Debug output
1334 vector<string> parnames = peak->getParameterNames();
1335 stringstream debugss;
1336 debugss << "DB1251 Single Peak Confident Fit Result: Chi^2 = " << chi2 << '\n';
1337 for (auto &parname : parnames) {
1338 debugss << parname << " = " << peak->getParameter(parname) << '\n';
1339 }
1340 g_log.notice(debugss.str());
1341
1342 return goodfit;
1343}
1344
1345//----------------------------------------------------------------------------------------------
1357void FitPowderDiffPeaks::calculatePeakFitBoundary(size_t ileftpeak, size_t irightpeak, double &peakleftboundary,
1358 double &peakrightboundary) {
1359 BackToBackExponential_sptr leftpeak = m_vecPeakFunctions[ileftpeak].second.second;
1360 BackToBackExponential_sptr rightpeak = m_vecPeakFunctions[irightpeak].second.second;
1361
1362 // 1. Determine its left boundary
1363 peakleftboundary = leftpeak->centre() - PEAKFITRANGEFACTOR * leftpeak->fwhm();
1364
1365 int ileftneighbor = static_cast<int>(ileftpeak) - 1;
1366 if (ileftneighbor < 0) {
1367 // a) No left neighbor, compare to TOF_Min
1368 if (peakleftboundary < m_tofMin)
1369 peakleftboundary = m_tofMin;
1370 } else {
1371 // b) Compare to the right peak boundary of its left neighbor
1372 BackToBackExponential_sptr leftneighbor = m_vecPeakFunctions[ileftneighbor].second.second;
1373 double leftneighborrightbound = leftneighbor->centre() + PEAKBOUNDARYFACTOR * leftneighbor->fwhm();
1374 if (leftneighborrightbound > peakleftboundary)
1375 peakleftboundary = leftneighborrightbound;
1376 }
1377
1378 // 2. Determine its right boundary
1379 peakrightboundary = rightpeak->centre() + PEAKFITRANGEFACTOR * rightpeak->fwhm();
1380
1381 size_t irightneighbor = irightpeak + 1;
1382 if (irightneighbor < m_vecPeakFunctions.size()) {
1383 // a) right peak exists
1384 BackToBackExponential_sptr rightneighbor = m_vecPeakFunctions[irightneighbor].second.second;
1385 double rightneighborleftbound = rightneighbor->centre() - PEAKBOUNDARYFACTOR * rightneighbor->fwhm();
1386 if (rightneighborleftbound < peakrightboundary)
1387 peakrightboundary = rightneighborleftbound;
1388 }
1389}
1390
1391//======================= Fit 1 Set of Overlapped Peaks ======================
1392
1393//----------------------------------------------------------------------------
1405std::pair<bool, double> FitPowderDiffPeaks::doFitPeak(const Workspace2D_sptr &dataws,
1406 const BackToBackExponential_sptr &peakfunction,
1407 double guessedfwhm) {
1408 // 1. Set up boundary
1409 if (guessedfwhm > 0) {
1410 double tof_h = peakfunction->centre();
1411 double centerleftend = tof_h - guessedfwhm * 3.0;
1412 double centerrightend = tof_h + guessedfwhm * 3.0;
1413 auto centerbound =
1414 std::make_unique<BoundaryConstraint>(peakfunction.get(), "X0", centerleftend, centerrightend, false);
1415 peakfunction->addConstraint(std::move(centerbound));
1416
1417 g_log.debug() << "[DoFitPeak] Peak Center Boundary = " << centerleftend << ", " << centerrightend << '\n';
1418 }
1419
1420 // A > 0, B > 0, S > 0
1421 auto abound = std::make_unique<BoundaryConstraint>(peakfunction.get(), "A", 0.0000001, DBL_MAX, false);
1422 peakfunction->addConstraint(std::move(abound));
1423
1424 auto bbound = std::make_unique<BoundaryConstraint>(peakfunction.get(), "B", 0.0000001, DBL_MAX, false);
1425 peakfunction->addConstraint(std::move(bbound));
1426
1427 auto sbound = std::make_unique<BoundaryConstraint>(peakfunction.get(), "S", 0.0001, DBL_MAX, false);
1428 peakfunction->addConstraint(std::move(sbound));
1429
1430 // 2. Unfix all parameters
1431 vector<string> paramnames = peakfunction->getParameterNames();
1432 size_t numparams = paramnames.size();
1433 for (size_t i = 0; i < numparams; ++i)
1434 peakfunction->unfix(i);
1435
1436 // 3. Set up the fitting scheme
1437 vector<vector<string>> vecMinimizers;
1438 vector<vector<size_t>> vecMaxIterations;
1439 vector<vector<double>> vecDampings;
1440
1441 // a) Scheme 1: Levenberg and Levenger
1442 /*
1443 vector<string> minimizers(2);
1444 minimizers[0] = "Levenberg-MarquardtMD";
1445 minimizers[1] = "Levenberg-Marquardt";
1446 vector<size_t> maxiterations(2, 1000);
1447 vector<double> dampings(2, 0.0);
1448 vecMinimizers.emplace_back(minimizers);
1449 vecMaxIterations.emplace_back(maxiterations);
1450 vecDampings.emplace_back(dampings);
1451 */
1452
1453 vector<string> minimizers2(3);
1454 minimizers2[0] = "Simplex";
1455 minimizers2[1] = "Levenberg-MarquardtMD";
1456 minimizers2[2] = "Levenberg-Marquardt";
1457 vector<size_t> maxiterations2(3, 1000);
1458 maxiterations2[0] = 10000;
1459 vector<double> dampings2(3, 0.0);
1460 vecMinimizers.emplace_back(minimizers2);
1461 vecMaxIterations.emplace_back(maxiterations2);
1462 vecDampings.emplace_back(dampings2);
1463
1464 // 4. Fit in different sequential
1465 bool goodfit = false;
1466 size_t numschemes = vecMinimizers.size();
1467
1468 map<string, double> origparams, bestparams;
1469 storeFunctionParameters(peakfunction, origparams);
1470 bestparams = origparams;
1471 double bestchi2 = DBL_MAX;
1472
1473 for (size_t is = 0; is < numschemes; ++is) {
1474 // a) Restore the starting value
1475 restoreFunctionParameters(peakfunction, origparams);
1476
1477 // b) Fit in multiple steps
1478 double thischi2;
1479 bool localgoodfit = doFit1PeakSequential(dataws, 1, peakfunction, vecMinimizers[is], vecMaxIterations[is],
1480 vecDampings[is], thischi2);
1481
1482 // c) Book keep
1483 if (localgoodfit && !goodfit) {
1484 // First local good fit
1485 bestchi2 = thischi2;
1486 storeFunctionParameters(peakfunction, bestparams);
1487 goodfit = true;
1488 } else if (localgoodfit && goodfit) {
1489 // Not the first time to have a good fit
1490 if (thischi2 < bestchi2) {
1491 bestchi2 = thischi2;
1492 storeFunctionParameters(peakfunction, bestparams);
1493 }
1494 }
1495 } // ENDFOR
1496
1497 pair<bool, double> returnvalue = make_pair(goodfit, bestchi2);
1498
1499 return returnvalue;
1500}
1501
1502//----------------------------------------------------------------------------
1505void FitPowderDiffPeaks::storeFunctionParameters(const IFunction_sptr &function, std::map<string, double> &parammaps) {
1506 vector<string> paramnames = function->getParameterNames();
1507 parammaps.clear();
1508 for (auto &paramname : paramnames)
1509 parammaps.emplace(paramname, function->getParameter(paramname));
1510}
1511
1512//----------------------------------------------------------------------------
1515void FitPowderDiffPeaks::restoreFunctionParameters(const IFunction_sptr &function, map<string, double> parammap) {
1516 vector<string> paramnames = function->getParameterNames();
1517 for (auto &parname : paramnames) {
1518 auto miter = parammap.find(parname);
1519 if (miter != parammap.end())
1520 function->setParameter(parname, miter->second);
1521 }
1522}
1523
1524//----------------------------------------------------------------------------
1537bool FitPowderDiffPeaks::doFit1PeakSimple(const Workspace2D_sptr &dataws, size_t workspaceindex,
1538 const BackToBackExponential_sptr &peakfunction, const string &minimzername,
1539 size_t maxiteration, double &chi2) {
1540 stringstream dbss;
1541 dbss << peakfunction->asString() << '\n';
1542 dbss << "Starting Value: ";
1543 vector<string> paramNames = peakfunction->getParameterNames();
1544 for (auto &paramName : paramNames)
1545 dbss << paramName << "= " << peakfunction->getParameter(paramName) << ", \t";
1546 for (size_t i = 0; i < dataws->x(workspaceindex).size(); ++i)
1547 dbss << dataws->x(workspaceindex)[i] << "\t\t" << dataws->y(workspaceindex)[i] << "\t\t"
1548 << dataws->e(workspaceindex)[i] << '\n';
1549 g_log.debug() << "DBx430 " << dbss.str() << '\n';
1550
1551 // 1. Peak height
1552 if (peakfunction->height() < 1.0E-5)
1553 peakfunction->setHeight(4.0);
1554
1555 // 2. Create fit
1556 Algorithm_sptr fitalg = createChildAlgorithm("Fit", -1, -1, true);
1557 fitalg->initialize();
1558
1559 // 3. Set
1560 fitalg->setProperty("Function", std::dynamic_pointer_cast<API::IFunction>(peakfunction));
1561 fitalg->setProperty("InputWorkspace", dataws);
1562 fitalg->setProperty("WorkspaceIndex", static_cast<int>(workspaceindex));
1563 fitalg->setProperty("Minimizer", minimzername);
1564 fitalg->setProperty("CostFunction", "Least squares");
1565 fitalg->setProperty("MaxIterations", static_cast<int>(maxiteration));
1566 fitalg->setProperty("Output", "FitPeak");
1567
1568 // 3. Execute and parse the result
1569 bool isexecute = fitalg->execute();
1570 bool fitsuccess = false;
1571 chi2 = DBL_MAX;
1572
1573 if (isexecute) {
1574 std::string fitresult = parseFitResult(fitalg, chi2, fitsuccess);
1575
1576 // Figure out result
1577 g_log.information() << "[DBx138A] Fit Peak @ " << peakfunction->centre() << " Result:" << fitsuccess << "\n"
1578 << "Detailed info = " << fitresult << ", Chi^2 = " << chi2 << '\n';
1579
1580 // Debug information output
1581 API::ITableWorkspace_sptr paramws = fitalg->getProperty("OutputParameters");
1582 std::string infofit = parseFitParameterWorkspace(paramws);
1583 g_log.information() << "Fitted Parameters: \n" << infofit << '\n';
1584 } else {
1585 g_log.error() << "[DBx128B] Failed to execute fitting peak @ " << peakfunction->centre() << '\n';
1586 }
1587
1588 return fitsuccess;
1589}
1590
1591//----------------------------------------------------------------------------------------------
1601bool FitPowderDiffPeaks::doFit1PeakSequential(const Workspace2D_sptr &dataws, size_t workspaceindex,
1602 const BackToBackExponential_sptr &peakfunction,
1603 vector<string> minimzernames, vector<size_t> maxiterations,
1604 const vector<double> &dampfactors, double &chi2) {
1605 // 1. Check
1606 if (minimzernames.size() != maxiterations.size() && minimzernames.size() != dampfactors.size()) {
1607 throw runtime_error("doFit1PeakSequential should have the input vectors of same size.");
1608 }
1609
1610 // 2. Start Chi2
1611 map<string, double> parambeforefit;
1612 storeFunctionParameters(peakfunction, parambeforefit);
1613
1614 double startchi2;
1615 doFit1PeakSimple(dataws, workspaceindex, peakfunction, "Levenberg-MarquardtMD", 0, startchi2);
1616
1617 restoreFunctionParameters(peakfunction, parambeforefit);
1618
1619 double currchi2 = startchi2;
1620 bool goodfit = false;
1621
1622 // 3. Fit sequentially
1623 for (size_t i = 0; i < minimzernames.size(); ++i) {
1624 string minimizer = minimzernames[i];
1625 size_t maxiteration = maxiterations[i];
1626 g_log.notice() << "DBx315 Start Chi2 = " << startchi2 << ", Minimizer = " << minimizer
1627 << ", Max Iterations = " << maxiteration << ", Workspace Index = " << workspaceindex
1628 << ", Data Range = " << dataws->x(workspaceindex)[0] << ", " << dataws->x(workspaceindex).back()
1629 << '\n';
1630
1631 storeFunctionParameters(peakfunction, parambeforefit);
1632
1633 double newchi2;
1634 bool localgoodfit = doFit1PeakSimple(dataws, workspaceindex, peakfunction, minimizer, maxiteration, newchi2);
1635
1636 if (localgoodfit && newchi2 < currchi2) {
1637 // A better solution
1638 currchi2 = newchi2;
1639 goodfit = true;
1640 } else {
1641 // A same of worse one
1642 restoreFunctionParameters(peakfunction, parambeforefit);
1643 g_log.information() << "DBx315C Fit Failed. Fit Good = " << localgoodfit << '\n';
1644 }
1645 }
1646
1647 // 4. Return
1648 chi2 = currchi2;
1649
1650 return goodfit;
1651}
1652
1653//----------------------------------------------------------------------------
1657 double in_center, double leftfwhm, double rightfwhm, double &center,
1658 double &sigma, double &height) {
1659 // 1. Estimate
1660 const auto &X = dataws->x(workspaceindex);
1661 const auto &Y = dataws->y(workspaceindex);
1662
1663 height = 0;
1664 for (size_t i = 1; i < X.size(); ++i) {
1665 height += (X[i] - X[i - 1]) * Y[i];
1666 }
1667 sigma = (leftfwhm + rightfwhm) * 0.5;
1668
1669 // 2. Use factory to generate Gaussian
1670 auto temppeak = API::FunctionFactory::Instance().createFunction("Gaussian");
1671 auto gaussianpeak = std::dynamic_pointer_cast<API::IPeakFunction>(temppeak);
1672 gaussianpeak->setHeight(height);
1673 gaussianpeak->setCentre(in_center);
1674 gaussianpeak->setFwhm(sigma);
1675
1676 // b) Constraint
1677 double centerleftend = in_center - leftfwhm * 0.5;
1678 double centerrightend = in_center + rightfwhm * 0.5;
1679 auto centerbound =
1680 std::make_unique<BoundaryConstraint>(gaussianpeak.get(), "PeakCentre", centerleftend, centerrightend, false);
1681 gaussianpeak->addConstraint(std::move(centerbound));
1682
1683 // 3. Fit
1684 auto fitalg = createChildAlgorithm("Fit", -1, -1, true);
1685 fitalg->initialize();
1686
1687 fitalg->setProperty("Function", std::dynamic_pointer_cast<API::IFunction>(gaussianpeak));
1688 fitalg->setProperty("InputWorkspace", dataws);
1689 fitalg->setProperty("WorkspaceIndex", 1);
1690 fitalg->setProperty("Minimizer", "Levenberg-MarquardtMD");
1691 fitalg->setProperty("CostFunction", "Least squares");
1692 fitalg->setProperty("MaxIterations", 1000);
1693 fitalg->setProperty("Output", "FitGaussianPeak");
1694
1695 // iv) Result
1696 bool successfulfit = fitalg->execute();
1697 if (!fitalg->isExecuted() || !successfulfit) {
1698 // Early return due to bad fit
1699 g_log.warning() << "Fitting Gaussian peak for peak around " << gaussianpeak->centre() << '\n';
1700 return false;
1701 }
1702
1703 double chi2;
1704 bool fitsuccess;
1705 std::string fitresult = parseFitResult(fitalg, chi2, fitsuccess);
1706 g_log.information() << "[Fit Gaussian Peak] Successful = " << fitsuccess << ", Result:\n" << fitresult << '\n';
1707
1708 // 4. Get result
1709 center = gaussianpeak->centre();
1710 height = gaussianpeak->height();
1711 double fwhm = gaussianpeak->fwhm();
1712 if (fwhm <= 0.0) {
1713 return false;
1714 }
1715 sigma = fwhm / 2.35;
1716
1717 // 5. Debug output
1718 API::ITableWorkspace_sptr paramws = fitalg->getProperty("OutputParameters");
1719 std::string infofit = parseFitParameterWorkspace(paramws);
1720 g_log.information() << "[DBx133] Fitted Gaussian Parameters: \n" << infofit << '\n';
1721
1722 return true;
1723}
1724
1725//----------------------------------------------------------------------------------------------
1731bool FitPowderDiffPeaks::fitOverlappedPeaks(vector<BackToBackExponential_sptr> peaks,
1732 const BackgroundFunction_sptr &backgroundfunction, double gfwhm) {
1733 // 1. Sort peak if necessary
1734 vector<pair<double, BackToBackExponential_sptr>> tofpeakpairs(peaks.size());
1735 for (size_t i = 0; i < peaks.size(); ++i)
1736 tofpeakpairs[i] = make_pair(peaks[i]->centre(), peaks[i]);
1737 sort(tofpeakpairs.begin(), tofpeakpairs.end());
1738
1739 // 2. Determine range of the data
1740 BackToBackExponential_sptr leftpeak = tofpeakpairs[0].second;
1741 BackToBackExponential_sptr rightpeak = tofpeakpairs.back().second;
1742 double peaksleftboundary, peaksrightboundary;
1743 if (gfwhm <= 0) {
1744 // Use peaks' preset value
1745 peaksleftboundary = leftpeak->centre() - 4 * leftpeak->fwhm();
1746 peaksrightboundary = rightpeak->centre() + 4 * rightpeak->fwhm();
1747 } else {
1748 // Use user input's guess fwhm
1749 peaksleftboundary = leftpeak->centre() - 4 * gfwhm;
1750 peaksrightboundary = rightpeak->centre() + 4 * gfwhm;
1751 }
1752
1753 // 3. Build partial data workspace
1754 Workspace2D_sptr peaksws = buildPartialWorkspace(m_dataWS, m_wsIndex, peaksleftboundary, peaksrightboundary);
1755
1756 // 4. Remove background
1757 estimateBackgroundCoarse(peaksws, backgroundfunction, 0, 2, 1);
1758
1759 // [DB] Debug output
1760 stringstream peakInfo;
1761 peakInfo << peaks.size() << "-Peaks Group Information: \n";
1762 for (size_t ipk = 0; ipk < tofpeakpairs.size(); ++ipk) {
1763 BackToBackExponential_sptr tmppeak = tofpeakpairs[ipk].second;
1764 peakInfo << "Peak " << ipk << " @ TOF = " << tmppeak->centre() << ", A = " << tmppeak->getParameter("A")
1765 << ", B = " << tmppeak->getParameter("B") << ", S = " << tmppeak->getParameter("S")
1766 << ", FWHM = " << tmppeak->fwhm() << '\n';
1767 }
1768 g_log.information() << "[DB1034] " << peakInfo.str();
1769
1770 stringstream peakWSData;
1771 peakWSData << "Partial workspace for peaks: \n";
1772 for (size_t i = 0; i < peaksws->x(0).size(); ++i)
1773 peakWSData << peaksws->x(1)[i] << "\t\t" << peaksws->y(1)[i] << "\t\t" << peaksws->e(1)[i] << "\t\t"
1774 << peaksws->y(0)[i] << '\n';
1775 g_log.information() << "[DB1042] " << peakWSData.str();
1776
1777 // 5. Estimate peak height according to pre-set peak value
1778 estimatePeakHeightsLeBail(peaksws, 1, peaks);
1779
1780 // 6. Set boundaries
1782
1783 // 7. Set up the composite function
1784 CompositeFunction_sptr peaksfunction(new CompositeFunction());
1785 for (auto &peak : peaks)
1786 peaksfunction->addFunction(peak);
1787
1788 // 8. Fit multiple peaks
1789 vector<double> chi2s;
1790 vector<bool> fitgoods;
1791 bool fitsuccess = doFitMultiplePeaks(peaksws, 1, peaksfunction, peaks, fitgoods, chi2s);
1792
1793 // 9. Plot peaks
1794 if (fitsuccess) {
1795 FunctionDomain1DVector domain(peaksws->x(1).rawData());
1796 plotFunction(peaksfunction, backgroundfunction, domain);
1797 }
1798
1799 return fitsuccess;
1800}
1801
1802//----------------------------------------------------------------------------------------------
1806 const CompositeFunction_sptr &peaksfunc,
1807 vector<BackToBackExponential_sptr> peakfuncs, vector<bool> &vecfitgood,
1808 vector<double> &vecchi2s) {
1809 // 0. Pre-debug output
1810 stringstream dbss;
1811 dbss << "[DBx529] Composite Function: " << peaksfunc->asString();
1812 g_log.notice(dbss.str());
1813
1814 // 1. Fit peaks intensities first
1815 const size_t numpeaks = peakfuncs.size();
1816 map<string, double> peaksfuncparams;
1817
1818 // a) Set up fit/fix
1819 vector<string> peakparnames = peakfuncs[0]->getParameterNames();
1820 for (size_t ipn = 0; ipn < peakparnames.size(); ++ipn) {
1821 bool isI = peakparnames[ipn] == "I";
1822
1823 for (size_t ipk = 0; ipk < numpeaks; ++ipk) {
1824 BackToBackExponential_sptr thispeak = peakfuncs[ipk];
1825 if (isI)
1826 thispeak->unfix(ipn);
1827 else
1828 thispeak->fix(ipn);
1829 }
1830 }
1831
1832 stringstream dbss1;
1833 dbss1 << "[DBx529A] Composite Function: " << peaksfunc->asString();
1834 g_log.notice(dbss1.str());
1835
1836 // b) Fit
1837 double chi2;
1838 bool fitgood = doFitNPeaksSimple(dataws, wsindex, peaksfunc, peakfuncs, "Levenberg-MarquardtMD", 1000, chi2);
1839 bool evergood = fitgood;
1840
1841 // c) Process result; possibly early return
1842 if (!fitgood) {
1843 vecfitgood.resize(numpeaks, false);
1844 vecchi2s.resize(numpeaks, -1.0);
1845 return false;
1846 } else {
1847 vecfitgood.resize(numpeaks, true);
1848 vecchi2s.resize(numpeaks, chi2);
1849 }
1850
1851 // 2. Fit A/B/S peak by peak
1852 for (size_t ipkfit = 0; ipkfit < numpeaks; ++ipkfit) {
1853 // a) Fix / unfix parameters
1854 for (size_t ipk = 0; ipk < numpeaks; ++ipk) {
1855 BackToBackExponential_sptr thispeak = peakfuncs[ipk];
1856 for (size_t iparam = 0; iparam < peakparnames.size(); ++iparam) {
1857 if (ipk == ipkfit) {
1858 // Peak to have parameters fit
1859 thispeak->unfix(iparam);
1860 } else {
1861 // Not the peak to fit, fix all
1862 thispeak->fix(iparam);
1863 }
1864 }
1865 }
1866
1867 // b) Fit
1868 storeFunctionParameters(peaksfunc, peaksfuncparams);
1869
1870 fitgood = doFitNPeaksSimple(dataws, wsindex, peaksfunc, peakfuncs, "Levenberg-MarquardtMD", 1000, chi2);
1871
1872 // not required. before loop starts, evergood=fitgood WITH fitgood==true
1873 // evergood = evergood || fitgood;
1874
1875 // c) Process the result
1876 if (!fitgood)
1877 restoreFunctionParameters(peaksfunc, peaksfuncparams);
1878 } // END OF FIT ALL SINGLE PEAKS
1879
1880 // 3. Fit all parameters (dangerous)
1881 for (size_t ipk = 0; ipk < numpeaks; ++ipk) {
1882 BackToBackExponential_sptr thispeak = peakfuncs[ipk];
1883 for (size_t iparam = 0; iparam < peakparnames.size(); ++iparam)
1884 thispeak->unfix(iparam);
1885 }
1886
1887 storeFunctionParameters(peaksfunc, peaksfuncparams);
1888 fitgood = doFitNPeaksSimple(dataws, wsindex, peaksfunc, peakfuncs, "Levenberg-MarquardtMD", 1000, chi2);
1889 // not required. before, evergood=fitgood WITH fitgood==true
1890 // evergood = evergood || fitgood;
1891
1892 if (!fitgood)
1893 restoreFunctionParameters(peaksfunc, peaksfuncparams);
1894
1895 // -1. Final debug output
1896 FunctionDomain1DVector domain(dataws->x(wsindex).rawData());
1897 FunctionValues values(domain);
1898 peaksfunc->function(domain, values);
1899 stringstream rss;
1900 for (size_t i = 0; i < domain.size(); ++i)
1901 rss << domain[i] << "\t\t" << values[i] << '\n';
1902 g_log.information() << "[T] Multiple peak fitting pattern:\n" << rss.str();
1903
1904 return evergood;
1905}
1906
1907//----------------------------------------------------------------------------------------------
1914 vector<BackToBackExponential_sptr> peaks) {
1915 // 1. Build data structures
1916 FunctionDomain1DVector domain(dataws->x(wsindex).rawData());
1917 FunctionValues values(domain);
1918 vector<vector<double>> peakvalues;
1919 for (size_t i = 0; i < (peaks.size() + 1); ++i) {
1920 vector<double> peakvalue(domain.size(), 0.0);
1921 peakvalues.emplace_back(peakvalue);
1922 }
1923
1924 // 2. Calcualte peak values
1925 size_t isum = peaks.size();
1926 for (size_t ipk = 0; ipk < peaks.size(); ++ipk) {
1927 BackToBackExponential_sptr thispeak = peaks[ipk];
1928 thispeak->setHeight(1.0);
1929 thispeak->function(domain, values);
1930 for (size_t j = 0; j < domain.size(); ++j) {
1931 peakvalues[ipk][j] = values[j];
1932 peakvalues[isum][j] += values[j];
1933 }
1934 }
1935
1936 // 3. Calculate peak height
1937 const auto &vecY = dataws->y(wsindex);
1938 for (size_t ipk = 0; ipk < peaks.size(); ++ipk) {
1939 double height = 0.0;
1940 for (size_t j = 0; j < domain.size() - 1; ++j) {
1941 if (vecY[j] > 0.0 && peakvalues[isum][j] > 1.0E-5) {
1942 double dtof = domain[j + 1] - domain[j];
1943 double temp = vecY[j] * peakvalues[ipk][j] / peakvalues[isum][j] * dtof;
1944 height += temp;
1945 }
1946 }
1947
1948 BackToBackExponential_sptr thispeak = peaks[ipk];
1949 thispeak->setHeight(height);
1950
1951 g_log.information() << "[DBx256] Peak @ " << thispeak->centre()
1952 << " Set Guessed Height (LeBail) = " << thispeak->height() << '\n';
1953 }
1954}
1955
1956//----------------------------------------------------------------------------------------------
1959void FitPowderDiffPeaks::setOverlappedPeaksConstraints(const vector<BackToBackExponential_sptr> &peaks) {
1960 for (const auto &thispeak : peaks) {
1961 // 1. Set constraint on X.
1962 double fwhm = thispeak->fwhm();
1963 double centre = thispeak->centre();
1964 double leftcentrebound = centre - 0.5 * fwhm;
1965 double rightcentrebound = centre + 0.5 * fwhm;
1966
1967 auto bc = std::make_unique<BoundaryConstraint>(thispeak.get(), "X0", leftcentrebound, rightcentrebound, false);
1968 thispeak->addConstraint(std::move(bc));
1969 }
1970}
1971
1972//----------------------------------------------------------------------------------------------
1976 const CompositeFunction_sptr &peaksfunc,
1977 const vector<BackToBackExponential_sptr> &peakfuncs,
1978 const string &minimizername, size_t maxiteration, double &chi2) {
1979 // 1. Debug output
1980 stringstream dbss0;
1981 dbss0 << "Starting Value: ";
1982 vector<string> paramNames = peaksfunc->getParameterNames();
1983 for (auto &paramName : paramNames)
1984 dbss0 << paramName << "= " << peaksfunc->getParameter(paramName) << ", \t";
1985 g_log.information() << "DBx430 " << dbss0.str() << '\n';
1986
1987 // 2. Create fit
1988 Algorithm_sptr fitalg = createChildAlgorithm("Fit", -1, -1, true);
1989 fitalg->initialize();
1990
1991 // 3. Set
1992 fitalg->setProperty("Function", std::dynamic_pointer_cast<API::IFunction>(peaksfunc));
1993 fitalg->setProperty("InputWorkspace", dataws);
1994 fitalg->setProperty("WorkspaceIndex", static_cast<int>(wsindex));
1995 fitalg->setProperty("Minimizer", minimizername);
1996 fitalg->setProperty("CostFunction", "Least squares");
1997 fitalg->setProperty("MaxIterations", static_cast<int>(maxiteration));
1998 fitalg->setProperty("Output", "FitPeak");
1999
2000 // 3. Execute and parse the result
2001 bool isexecute = fitalg->execute();
2002 bool fitsuccess = false;
2003 chi2 = DBL_MAX;
2004
2005 // 4. Output
2006 stringstream dbss;
2007 dbss << "Fit N-Peaks @ ";
2008 for (auto &peakfunc : peakfuncs)
2009 dbss << peakfunc->centre() << ", ";
2010
2011 if (isexecute) {
2012 // Figure out result
2013 std::string fitresult = parseFitResult(fitalg, chi2, fitsuccess);
2014
2015 dbss << " Result:" << fitsuccess << "\n"
2016 << "Detailed info = " << fitresult;
2017
2018 g_log.information() << "[DBx149A] " << dbss.str() << '\n';
2019
2020 // Debug information output
2021 API::ITableWorkspace_sptr paramws = fitalg->getProperty("OutputParameters");
2022 std::string infofit = parseFitParameterWorkspace(paramws);
2023 g_log.information() << "[DBx149B] Fitted Parameters: \n" << infofit << '\n';
2024 } else {
2025 dbss << ": Failed ";
2026 g_log.error() << "[DBx149C] " << dbss.str() << '\n';
2027 }
2028
2029 return fitsuccess;
2030}
2031
2032//=================================== Process Fit Result
2033//=====================================
2034//----------------------------------------------------------------------------------------------
2037std::string FitPowderDiffPeaks::parseFitResult(const API::IAlgorithm_sptr &fitalg, double &chi2, bool &fitsuccess) {
2038 stringstream rss;
2039
2040 chi2 = fitalg->getProperty("OutputChi2overDoF");
2041 string fitstatus = fitalg->getProperty("OutputStatus");
2042
2043 fitsuccess = (fitstatus == "success");
2044
2045 rss << " [Algorithm Fit]: Chi^2 = " << chi2 << "; Fit Status = " << fitstatus;
2046
2047 return rss.str();
2048}
2049
2050//----------------------------------------------------------------------------------------------
2054 // 1. Check
2055 if (!paramws) {
2056 g_log.warning() << "Input table workspace is NULL. \n";
2057 return "";
2058 }
2059
2060 // 2. Parse
2061 std::stringstream msgss;
2062 size_t numrows = paramws->rowCount();
2063 for (size_t i = 0; i < numrows; ++i) {
2064 API::TableRow row = paramws->getRow(i);
2065 std::string parname;
2066 double parvalue, parerror;
2067 row >> parname >> parvalue >> parerror;
2068
2069 msgss << parname << " = " << setw(10) << setprecision(5) << parvalue << " +/- " << setw(10) << setprecision(5)
2070 << parerror << '\n';
2071 }
2072
2073 return msgss.str();
2074}
2075
2076//================================ Process Input/Output
2077//======================================
2078//----------------------------------------------------------------------------------------------
2083 // 1. Check column orders
2084 std::vector<std::string> colnames = parameterWS->getColumnNames();
2085 if (colnames.size() < 2) {
2086 stringstream errss;
2087 errss << "Input parameter table workspace does not have enough number of "
2088 "columns. "
2089 << " Number of columns = " << colnames.size() << " >= 2 as required. ";
2090 g_log.error(errss.str());
2091 throw std::runtime_error(errss.str());
2092 }
2093
2094 if (colnames[0] != "Name" || colnames[1] != "Value") {
2095 stringstream errss;
2096 errss << "Input parameter table workspace does not have the columns in "
2097 "order as "
2098 << "Name, Value and etc. ";
2099 g_log.error(errss.str());
2100 throw std::runtime_error(errss.str());
2101 }
2102
2103 size_t numrows = parameterWS->rowCount();
2104
2105 g_log.notice() << "[DBx409] Import TableWorkspace " << parameterWS->getName() << " containing " << numrows
2106 << " instrument profile parameters\n";
2107
2108 // 2. Import data to maps
2109 std::string parname;
2110 double value;
2111 m_instrumentParmaeters.clear();
2112
2113 for (size_t ir = 0; ir < numrows; ++ir) {
2114 API::TableRow trow = parameterWS->getRow(ir);
2115 trow >> parname >> value;
2116 m_instrumentParmaeters.emplace(parname, value);
2117 g_log.notice() << "[DBx211] Import parameter " << parname << ": " << value << '\n';
2118 }
2119}
2120
2123void FitPowderDiffPeaks::parseBraggPeakTable(const TableWorkspace_sptr &peakws, vector<map<string, double>> &parammaps,
2124 vector<map<string, int>> &hklmaps) {
2125 // 1. Get columns' types and names
2126 vector<string> paramnames = peakws->getColumnNames();
2127 size_t numcols = paramnames.size();
2128 vector<string> coltypes(numcols);
2129 for (size_t i = 0; i < numcols; ++i) {
2130 Column_sptr col = peakws->getColumn(i);
2131 string coltype = col->type();
2132 coltypes[i] = coltype;
2133 }
2134
2135 // 2. Parse table rows
2136 size_t numrows = peakws->rowCount();
2137 for (size_t irow = 0; irow < numrows; ++irow) {
2138 // 1. Create map
2139 map<string, int> intmap;
2140 map<string, double> doublemap;
2141
2142 // 2. Parse
2143 for (size_t icol = 0; icol < numcols; ++icol) {
2144 string coltype = coltypes[icol];
2145 string colname = paramnames[icol];
2146
2147 if (coltype == "int") {
2148 // Integer
2149 int temp = peakws->cell<int>(irow, icol);
2150 intmap.emplace(colname, temp);
2151 } else if (coltype == "double") {
2152 // Double
2153 double temp = peakws->cell<double>(irow, icol);
2154 doublemap.emplace(colname, temp);
2155 }
2156
2157 } // ENDFOR Column
2158
2159 parammaps.emplace_back(doublemap);
2160 hklmaps.emplace_back(intmap);
2161
2162 } // ENDFOR Row
2163
2164 g_log.information() << "Import " << hklmaps.size() << " entries from Bragg peak TableWorkspace " << peakws->getName()
2165 << '\n';
2166}
2167
2168//----------------------------------------------------------------------------
2173 int workspaceindex) {
2174 // 1. Init
2175 const auto &X = m_dataWS->x(workspaceindex);
2176
2177 if (pattern.size() != X.size()) {
2178 stringstream errmsg;
2179 errmsg << "Input pattern (" << pattern.size() << ") and algorithm's input workspace (" << X.size()
2180 << ") have different size. ";
2181 g_log.error() << errmsg.str() << '\n';
2182 throw std::logic_error(errmsg.str());
2183 }
2184
2185 // 2. Create data workspace
2186 Workspace2D_sptr dataws = std::dynamic_pointer_cast<Workspace2D>(
2187 WorkspaceFactory::Instance().create("Workspace2D", 5, pattern.size(), pattern.size()));
2188
2189 // 3. Set up
2190 for (size_t iw = 0; iw < 5; ++iw) {
2191 dataws->setSharedX(iw, m_dataWS->sharedX(workspaceindex));
2192 }
2193
2194 dataws->setSharedY(0, m_dataWS->sharedY(workspaceindex));
2195 dataws->mutableY(1) = pattern;
2196 dataws->mutableY(2) = m_dataWS->y(workspaceindex) - pattern;
2197
2198 return dataws;
2199}
2200
2201//----------------------------------------------------------------------------
2205 // 1. Check and prepare
2206 if (m_vecPeakFunctions.size() != m_peakFitChi2.size()) {
2207 throw runtime_error("Wrong definition of m_peakFitChi2");
2208 }
2209
2210 size_t numpeaks = m_vecPeakFunctions.size();
2211
2212 // 2. Collect parameters of peak fitted good
2213 vector<double> vecdh, vectofh, vecalpha, vecbeta, vecsigma, vecchi2;
2214 for (size_t i = 0; i < numpeaks; ++i) {
2215 double &chi2 = m_peakFitChi2[i];
2216 if (chi2 > 0) {
2217 // a) Get values
2218 double dh = m_vecPeakFunctions[i].first;
2219 BackToBackExponential_sptr peak = m_vecPeakFunctions[i].second.second;
2220
2221 double p_a = peak->getParameter("A");
2222 double p_b = peak->getParameter("B");
2223 double p_x = peak->getParameter("X0");
2224 double p_s = peak->getParameter("S");
2225
2226 // b) To vectors
2227 vecchi2.emplace_back(chi2);
2228 vecdh.emplace_back(dh);
2229 vectofh.emplace_back(p_x);
2230 vecalpha.emplace_back(p_a);
2231 vecbeta.emplace_back(p_b);
2232 vecsigma.emplace_back(p_s);
2233 }
2234 } // ENDFOR i(peak)
2235
2236 // 3. Create workspace2D
2237 size_t numgoodpeaks = vecdh.size();
2238 Workspace2D_sptr paramws = std::dynamic_pointer_cast<Workspace2D>(
2239 WorkspaceFactory::Instance().create("Workspace2D", 4, numgoodpeaks, numgoodpeaks));
2240 for (size_t j = 0; j < 4; ++j) {
2241 paramws->mutableX(j) = vecdh;
2242 paramws->mutableE(j) = vecchi2;
2243 }
2244 paramws->mutableY(0) = vectofh;
2245 paramws->mutableY(1) = vecalpha;
2246 paramws->mutableY(2) = vecbeta;
2247 paramws->mutableY(3) = vecsigma;
2248
2249 // 4. Set Axis label
2250 paramws->getAxis(0)->setUnit("dSpacing");
2251
2252 auto tAxis = std::make_unique<TextAxis>(4);
2253 tAxis->setLabel(0, "X0");
2254 tAxis->setLabel(1, "A");
2255 tAxis->setLabel(2, "B");
2256 tAxis->setLabel(3, "S");
2257
2258 paramws->replaceAxis(1, std::move(tAxis));
2259
2260 return paramws;
2261}
2262
2263//----------------------------------------------------------------------------
2268pair<TableWorkspace_sptr, TableWorkspace_sptr> FitPowderDiffPeaks::genPeakParametersWorkspace() {
2269 // 1. Debug/Test Output
2270 for (size_t i = 0; i < m_vecPeakFunctions.size(); ++i) {
2271 g_log.debug() << "Peak @ d = " << m_vecPeakFunctions[i].first << ": Chi^2 = " << m_peakFitChi2[i] << '\n';
2272 }
2273
2274 if (m_vecPeakFunctions.size() != m_peakFitChi2.size()) {
2275 throw runtime_error("Wrong definition of m_peakFitChi2");
2276 }
2277
2278 size_t numpeaks = m_vecPeakFunctions.size();
2279 vector<double> vectofh(numpeaks), vecalpha(numpeaks), vecbeta(numpeaks), vecsigma(numpeaks);
2280
2281 // 2. Generate the TableWorkspace for peak parameters
2282 auto tablewsptr = new TableWorkspace();
2283 TableWorkspace_sptr tablews = TableWorkspace_sptr(tablewsptr);
2284
2285 tablews->addColumn("int", "H");
2286 tablews->addColumn("int", "K");
2287 tablews->addColumn("int", "L");
2288
2289 tablews->addColumn("double", "d_h");
2290 tablews->addColumn("double", "TOF_h");
2291 tablews->addColumn("double", "Height");
2292 tablews->addColumn("double", "Alpha");
2293 tablews->addColumn("double", "Beta");
2294 tablews->addColumn("double", "Sigma");
2295 tablews->addColumn("double", "Chi2");
2296
2297 /*
2298 stringstream outbuf;
2299 outbuf << setw(10) << "H" << setw(10) << "K" << setw(10) << "L"
2300 << setw(10) << "d_h" << setw(10) << "X0" << setw(10) << "I"
2301 << setw(10) << "A" << setw(10) << "B" << setw(10) << "S" << setw(10) <<
2302 "Chi2\n";
2303 */
2304
2305 for (size_t i = 0; i < numpeaks; ++i) {
2306 const double &chi2 = m_peakFitChi2[i];
2307 if (chi2 > 0) {
2308 // Bad fit peak has chi^2 < 0;
2309 double dh = m_vecPeakFunctions[i].first;
2310 const vector<int> &hkl = m_vecPeakFunctions[i].second.first;
2311 BackToBackExponential_sptr peak = m_vecPeakFunctions[i].second.second;
2312
2313 TableRow newrow = tablews->appendRow();
2314
2315 // i. H, K, L, d_h
2316 newrow << hkl[0] << hkl[1] << hkl[2] << dh;
2317
2318 // ii. A, B, I, S, X0
2319 double p_a = peak->getParameter("A");
2320 double p_b = peak->getParameter("B");
2321 double p_i = peak->getParameter("I");
2322 double p_x = peak->getParameter("X0");
2323 double p_s = peak->getParameter("S");
2324 newrow << p_x << p_i << p_a << p_b << p_s;
2325
2326 // iii. Chi^2
2327 newrow << chi2;
2328
2329 // iv. Prepare for Z-score
2330 // vectofh.emplace_back(p_x);
2331 // vecalpha.emplace_back(p_a);
2332 // vecbeta.emplace_back(p_b);
2333 // vecsigma.emplace_back(p_s);
2334 vectofh[i] = p_x;
2335 vecalpha[i] = p_a;
2336 vecbeta[i] = p_b;
2337 vecsigma[i] = p_s;
2338
2339 /*
2340 outbuf << setw(10) << setprecision(5) << chi2 << '\n';
2341 outbuf << setw(10) << hkl[0] << setw(10) << hkl[1] << setw(10) << hkl[2]
2342 << setw(10) << setprecision(5) << dh;
2343 outbuf << setw(10) << setprecision(5) << p_x
2344 << setw(10) << setprecision(5) << p_i
2345 << setw(10) << setprecision(5) << p_a
2346 << setw(10) << setprecision(5) << p_b
2347 << setw(10) << setprecision(5) << p_s;
2348 */
2349 }
2350 } // ENDFOR i(peak)
2351
2352 // 3. Z-score table
2353 // i. Calculate Z-scores
2354 vector<double> zcentres = Kernel::getZscore(vectofh);
2355 vector<double> zalphas = Kernel::getZscore(vecalpha);
2356 vector<double> zbetas = Kernel::getZscore(vecbeta);
2357 vector<double> zsigma = Kernel::getZscore(vecsigma);
2358
2359 // ii. Build table workspace for Z scores
2360 auto ztablewsptr = new TableWorkspace();
2361 TableWorkspace_sptr ztablews = TableWorkspace_sptr(ztablewsptr);
2362
2363 ztablews->addColumn("int", "H");
2364 ztablews->addColumn("int", "K");
2365 ztablews->addColumn("int", "L");
2366
2367 ztablews->addColumn("double", "d_h");
2368 ztablews->addColumn("double", "Z_TOF_h");
2369 ztablews->addColumn("double", "Z_Alpha");
2370 ztablews->addColumn("double", "Z_Beta");
2371 ztablews->addColumn("double", "Z_Sigma");
2372
2373 // iii. Set values
2374 for (size_t i = 0; i < m_vecPeakFunctions.size(); ++i) {
2375 double chi2 = m_peakFitChi2[i];
2376 if (chi2 > 0) {
2377 // A good fit has chi^2 larger than 0
2378 double dh = m_vecPeakFunctions[i].first;
2379 const vector<int> &hkl = m_vecPeakFunctions[i].second.first;
2380
2381 TableRow newrow = ztablews->appendRow();
2382 newrow << hkl[0] << hkl[1] << hkl[2] << dh;
2383
2384 // ii. Z scores
2385 double p_x = zcentres[i];
2386 double p_a = zalphas[i];
2387 double p_b = zbetas[i];
2388 double p_s = zsigma[i];
2389
2390 newrow << p_x << p_a << p_b << p_s;
2391 }
2392 } // ENDFOR i(peak)
2393
2394 return make_pair(tablews, ztablews);
2395}
2396
2397//----------------------------------------------------------------------------------------------
2403 // Check and clear input and output
2404 if (!peakparamws) {
2405 stringstream errss;
2406 errss << "Input tableworkspace for peak parameters is invalid!";
2407 g_log.error(errss.str());
2408 throw std::invalid_argument(errss.str());
2409 }
2410
2411 m_vecPeakFunctions.clear();
2412
2413 // Give name to peak parameters
2414 BackToBackExponential tempeak;
2415 tempeak.initialize();
2417 mPeakParameterNames.emplace_back("S2");
2418
2419 // Parse TableWorkspace
2420 vector<map<std::string, double>> peakparametermaps;
2421 vector<map<std::string, int>> peakhkls;
2422 this->parseBraggPeakTable(peakparamws, peakparametermaps, peakhkls);
2423
2424 // Create a map to convert the Bragg peak Table paramter name to Back to back
2425 // exponential+pseudo-voigt
2426 map<string, string> bk2bk2braggmap{{"A", "Alpha"}, {"B", "Beta"}, {"X0", "TOF_h"},
2427 {"I", "Height"}, {"S", "Sigma"}, {"S2", "Sigma2"}};
2428
2429 // Generate Peaks
2430 size_t numrows = peakparamws->rowCount();
2431 for (size_t ir = 0; ir < numrows; ++ir) {
2432 double d_h;
2433 vector<int> hkl;
2434 bool good;
2435 BackToBackExponential_sptr newpeak = genPeak(peakhkls[ir], peakparametermaps[ir], bk2bk2braggmap, good, hkl, d_h);
2436
2437 if (good) {
2438 m_vecPeakFunctions.emplace_back(d_h, make_pair(hkl, newpeak));
2439 }
2440 } // ENDFOR Each potential peak
2441
2442 // Sort by peaks' centre in d-space
2443 sort(m_vecPeakFunctions.begin(), m_vecPeakFunctions.end());
2444
2445 // Remove all peaks outside of tof_min and tof_max
2446 double tofmin = m_dataWS->x(m_wsIndex)[0];
2447 double tofmax = m_dataWS->x(m_wsIndex).back();
2448
2449 stringstream dbss;
2450 dbss << "Specified range for peaks in TOF: " << tofmin << ", " << tofmax << "\n";
2451
2452 vector<pair<double, pair<vector<int>, BackToBackExponential_sptr>>>::iterator deliter;
2453 for (deliter = m_vecPeakFunctions.begin(); deliter != m_vecPeakFunctions.end(); ++deliter) {
2454 double d_h = deliter->first;
2455 const vector<int> &hkl = deliter->second.first;
2456 BackToBackExponential_sptr peak = deliter->second.second;
2457 double tofh = peak->getParameter("X0");
2458
2459 dbss << "\tPeak (" << hkl[0] << ", " << hkl[1] << ", " << hkl[2] << ") @ d = " << d_h
2460 << " (calculated TOF) = " << tofh << " ";
2461
2462 if (tofh < tofmin || tofh > tofmax) {
2463 // Remove peak from vector
2464 deliter = m_vecPeakFunctions.erase(deliter);
2465 dbss << ": Delete due to out of range. \n";
2466 if (deliter == m_vecPeakFunctions.end()) {
2467 break;
2468 }
2469 } else {
2470 dbss << ": Kept. \n";
2471 }
2472 }
2473
2474 g_log.notice() << "[DBx453] " << dbss.str() << "\n";
2475
2476 // Remove peaks lower than minimum
2477 if (m_minimumHKL.size() == 3) {
2478 // Only keep peaks from and above minimum HKL
2479 for (deliter = m_vecPeakFunctions.begin(); deliter != m_vecPeakFunctions.end(); ++deliter) {
2480 vector<int> hkl = deliter->second.first;
2481 if (hkl == m_minimumHKL) {
2482 break;
2483 }
2484 }
2485 if (deliter != m_vecPeakFunctions.end()) {
2486 // Find the real minum
2487 int indminhkl = static_cast<int>(deliter - m_vecPeakFunctions.begin());
2488 int ind1stpeak = indminhkl - m_numPeaksLowerToMin;
2489 if (ind1stpeak > 0) {
2490 deliter = m_vecPeakFunctions.begin() + ind1stpeak;
2491 m_vecPeakFunctions.erase(m_vecPeakFunctions.begin(), deliter);
2492 }
2493 } else {
2494 // Miminum HKL peak does not exist
2495 vector<int> hkl = m_minimumHKL;
2496 g_log.warning() << "Minimum peak " << hkl[0] << ", " << hkl[1] << ", " << hkl[2] << " does not exit. \n";
2497 }
2498 }
2499
2500 // Keep some input information
2501 stringstream dbout;
2502 for (deliter = m_vecPeakFunctions.begin(); deliter != m_vecPeakFunctions.end(); ++deliter) {
2503 vector<int> hkl = deliter->second.first;
2504 double d_h = deliter->first;
2505 double tof_h = deliter->second.second->centre();
2506 dbout << "Peak (" << hkl[0] << ", " << hkl[1] << ", " << hkl[2] << ") @ d = " << d_h << ", TOF = " << tof_h << '\n';
2507 }
2508 g_log.information() << "[DBx531] Peaks To Fit: Number of peaks = " << m_vecPeakFunctions.size() << "\n"
2509 << dbout.str();
2510}
2511
2512//-----------------------------------------------------------------------------------------
2523BackToBackExponential_sptr FitPowderDiffPeaks::genPeak(map<string, int> hklmap, map<string, double> parammap,
2524 map<string, string> bk2bk2braggmap, bool &good, vector<int> &hkl,
2525 double &d_h) {
2526 // Generate a peak function
2527 BackToBackExponential_sptr newpeakptr = std::make_shared<BackToBackExponential>();
2528 newpeakptr->initialize();
2529
2530 // Check miller index (HKL) is a valid value in a miller indexes pool (hklmap)
2531 good = getHKLFromMap(std::move(hklmap), hkl);
2532 if (!good) {
2533 // Ignore and return
2534 return newpeakptr;
2535 }
2536
2537 // Set the peak parameters either by calculating or from table workspace
2538 string peakcalmode;
2540 // a) Use Bragg peak table's (HKL) and calculate the peak parameters
2541 double alph0 = getParameter("Alph0");
2542 double alph1 = getParameter("Alph1");
2543 double alph0t = getParameter("Alph0t");
2544 double alph1t = getParameter("Alph1t");
2545 double beta0 = getParameter("Beta0");
2546 double beta1 = getParameter("Beta1");
2547 double beta0t = getParameter("Beta0t");
2548 double beta1t = getParameter("Beta1t");
2549 double sig0 = getParameter("Sig0");
2550 double sig1 = getParameter("Sig1");
2551 double sig2 = getParameter("Sig2");
2552 double tcross = getParameter("Tcross");
2553 double width = getParameter("Width");
2554 double dtt1 = getParameter("Dtt1");
2555 double dtt1t = getParameter("Dtt1t");
2556 double dtt2t = getParameter("Dtt2t");
2557 double zero = getParameter("Zero");
2558 double zerot = getParameter("Zerot");
2559
2560 // b) Check validity and make choice
2561 if (tcross == EMPTY_DBL() || width == EMPTY_DBL() || dtt1 == EMPTY_DBL() || dtt1t == EMPTY_DBL() ||
2562 dtt2t == EMPTY_DBL() || zero == EMPTY_DBL() || zerot == EMPTY_DBL()) {
2563 stringstream errss;
2564 errss << "In input InstrumentParameterTable, one of the following is not "
2565 "given. Unable to process. \n";
2566 errss << "Tcross = " << tcross << "; Width = " << width << ", Dtt1 = " << dtt1 << ", Dtt1t = " << dtt1t << '\n';
2567 errss << "Dtt2t = " << dtt2t << ", Zero = " << zero << ", Zerot = " << zerot;
2568
2569 g_log.error(errss.str());
2570 throw runtime_error(errss.str());
2571 }
2572
2573 bool caltofonly = false;
2574 if (alph0 == EMPTY_DBL() || alph1 == EMPTY_DBL() || alph0t == EMPTY_DBL() || alph1t == EMPTY_DBL() ||
2575 beta0 == EMPTY_DBL() || beta1 == EMPTY_DBL() || beta0t == EMPTY_DBL() || beta1t == EMPTY_DBL() ||
2576 sig0 == EMPTY_DBL() || sig1 == EMPTY_DBL() || sig2 == EMPTY_DBL()) {
2577 caltofonly = true;
2578 g_log.warning() << "[DBx343] At least one of the input instrument-peak "
2579 "profile parameters is not given. "
2580 << "Use (HKL) only!\n";
2581 g_log.warning() << "Alph0 = " << alph0 << ", Alph1 = " << alph1 << ", Alph0t = " << alph0t
2582 << ", Alph1t = " << alph1t << '\n';
2583 g_log.warning() << "Beta0 = " << beta0 << ", Beta1 = " << beta1 << ", Beta0t = " << beta0t
2584 << ", Beta1t = " << beta1t << '\n';
2585 g_log.warning() << "Sig0 = " << sig0 << ", Sig1 = " << sig1 << ", Sig2 = " << sig2 << '\n';
2586 }
2587
2588 if (caltofonly) {
2589 // c) Calcualte d->TOF only
2590 // Calcualte d-spacing
2591 d_h = m_unitCell.d(hkl[0], hkl[1], hkl[2]);
2592 // d_h = calCubicDSpace(latticesize, hkl[0], hkl[1],
2593 // hkl[2]);
2594 if ((d_h != d_h) || (d_h < -DBL_MAX) || (d_h > DBL_MAX)) {
2595 stringstream warnss;
2596 warnss << "Peak with Miller Index = " << hkl[0] << ", " << hkl[1] << ", " << hkl[2]
2597 << " has unphysical d-spacing value = " << d_h;
2598 g_log.warning(warnss.str());
2599 good = false;
2600 return newpeakptr;
2601 }
2602
2603 // Calcualte TOF_h
2604 double tof_h = calThermalNeutronTOF(d_h, dtt1, dtt1t, dtt2t, zero, zerot, width, tcross);
2605 newpeakptr->setCentre(tof_h);
2606 } else {
2607 // d) Calculate a lot of peak parameters
2608 // Initialize the function
2610 tnb2bfunc.initialize();
2611 tnb2bfunc.setMillerIndex(hkl[0], hkl[1], hkl[2]);
2612
2613 vector<string> tnb2bfuncparnames = tnb2bfunc.getParameterNames();
2614
2615 // Set peak parameters
2616 for (const auto &parname : tnb2bfuncparnames) {
2617 if (parname != "Height") {
2618 auto miter = m_instrumentParmaeters.find(parname);
2619 if (miter == m_instrumentParmaeters.end()) {
2620 stringstream errss;
2621 errss << "Cannot find peak parameter " << parname << " in input instrument parameter "
2622 << "TableWorkspace. This mode is unable to execute. Quit!";
2623 g_log.error(errss.str());
2624 throw runtime_error(errss.str());
2625 }
2626 double parvalue = miter->second;
2627 tnb2bfunc.setParameter(parname, parvalue);
2628 }
2629 }
2630
2631 // Calculate peak parameters A, B, S, and X0
2632 tnb2bfunc.calculateParameters(false);
2633 d_h = tnb2bfunc.getPeakParameter("d_h");
2634 double alpha = tnb2bfunc.getPeakParameter("Alpha");
2635 double beta = tnb2bfunc.getPeakParameter("Beta");
2636 double sigma2 = tnb2bfunc.getPeakParameter("Sigma2");
2637 double tof_h = tnb2bfunc.centre();
2638
2639 newpeakptr->setParameter("A", alpha);
2640 newpeakptr->setParameter("B", beta);
2641 newpeakptr->setParameter("S", sqrt(sigma2));
2642 newpeakptr->setParameter("X0", tof_h);
2643 }
2644
2645 peakcalmode = "Calculate all parameters by thermal neutron peak function.";
2646
2647 } // Calculate peak variable from Profile function
2649 // e) Import from input table workspace
2650 vector<string>::iterator nameiter;
2651 for (nameiter = mPeakParameterNames.begin(); nameiter != mPeakParameterNames.end(); ++nameiter) {
2652 // BackToBackExponential parameter name
2653 string b2bexpname = *nameiter;
2654 // Map to instrument parameter
2655 map<string, string>::iterator refiter;
2656 refiter = bk2bk2braggmap.find(b2bexpname);
2657 if (refiter == bk2bk2braggmap.end())
2658 throw runtime_error("Programming error!");
2659 string instparname = refiter->second;
2660
2661 // Search in Bragg peak table
2662 map<string, double>::iterator miter;
2663 miter = parammap.find(instparname);
2664 if (miter != parammap.end()) {
2665 // Parameter exist in input
2666 double parvalue = miter->second;
2667 if (b2bexpname == "S2") {
2668 newpeakptr->setParameter("S", sqrt(parvalue));
2669 } else {
2670 newpeakptr->setParameter(b2bexpname, parvalue);
2671 }
2672 }
2673 } // FOR PEAK TABLE CELL
2674
2675 peakcalmode = "Import from Bragg peaks table";
2676 }
2677
2678 // Debug output
2679 stringstream infoss;
2680 string peakinfo = getFunctionInfo(std::dynamic_pointer_cast<IFunction>(newpeakptr));
2681
2682 infoss << "Generate Peak (" << hkl[0] << ", " << hkl[1] << ", " << hkl[2] << ") Of Mode " << peakcalmode << ".\n";
2683 g_log.notice(infoss.str());
2684 g_log.information() << "[DBx426B] " << peakinfo << ".\n";
2685
2686 good = true;
2687
2688 return newpeakptr;
2689}
2690
2691//----------------------------------------------------------------------------------------------
2700 const FunctionDomain1DVector &domain) {
2701 // 1. Determine range
2702 const auto &vecX = m_dataWS->x(m_wsIndex);
2703 double x0 = domain[0];
2704 auto viter = lower_bound(vecX.cbegin(), vecX.cend(), x0);
2705 int ix0 = static_cast<int>(std::distance(vecX.cbegin(), viter));
2706
2707 // Check boundary
2708 if ((static_cast<int>(domain.size()) + ix0) > static_cast<int>(m_peakData.size()))
2709 throw runtime_error("Plot single peak out of boundary error!");
2710
2711 // 2. Calculation of peaks
2712 FunctionValues values1(domain);
2713 peakfunction->function(domain, values1);
2714
2715 for (int i = 0; i < static_cast<int>(domain.size()); ++i)
2716 m_peakData[i + ix0] = values1[i];
2717
2718 // 3.Calculation of background
2719 FunctionValues values2(domain);
2720 background->function(domain, values2);
2721
2722 for (int i = 0; i < static_cast<int>(domain.size()); ++i)
2723 m_peakData[i + ix0] += values2[i];
2724}
2725
2726//===================================== Auxiliary Functions
2727//===================================
2728//----------------------------------------------------------------------------------------------
2734bool FitPowderDiffPeaks::getHKLFromMap(map<string, int> intmap, vector<int> &hkl) {
2735 vector<string> strhkl(3);
2736 strhkl[0] = "H";
2737 strhkl[1] = "K";
2738 strhkl[2] = "L";
2739
2740 hkl.clear();
2741
2742 map<string, int>::iterator miter;
2743 for (size_t i = 0; i < 3; ++i) {
2744 string parname = strhkl[i];
2745 miter = intmap.find(parname);
2746 if (miter == intmap.end())
2747 return false;
2748
2749 hkl.emplace_back(miter->second);
2750 }
2751
2752 return true;
2753}
2754
2755//----------------------------------------------------------------------------------------------
2760void FitPowderDiffPeaks::cropWorkspace(double tofmin, double tofmax) {
2761 auto cropalg = createChildAlgorithm("CropWorkspace", -1, -1, true);
2762 cropalg->initialize();
2763
2764 cropalg->setProperty("InputWorkspace", m_dataWS);
2765 cropalg->setPropertyValue("OutputWorkspace", "MyData");
2766 cropalg->setProperty("XMin", tofmin);
2767 cropalg->setProperty("XMax", tofmax);
2768
2769 bool cropstatus = cropalg->execute();
2770
2771 std::stringstream errmsg;
2772 if (!cropstatus) {
2773 errmsg << "DBx309 Cropping workspace unsuccessful. Fatal Error. Quit!";
2774 g_log.error() << errmsg.str() << '\n';
2775 throw std::runtime_error(errmsg.str());
2776 }
2777
2778 m_dataWS = cropalg->getProperty("OutputWorkspace");
2779 if (!m_dataWS) {
2780 errmsg << "Unable to retrieve a Workspace2D object from ChildAlgorithm Crop.";
2781 g_log.error(errmsg.str());
2782 throw std::runtime_error(errmsg.str());
2783 } else {
2784 g_log.information() << "[DBx211] Cropped Workspace Range: " << m_dataWS->x(m_wsIndex)[0] << ", "
2785 << m_dataWS->x(m_wsIndex).back() << '\n';
2786 }
2787}
2788
2789//----------------------------------------------------------------------------------------------
2794double FitPowderDiffPeaks::getParameter(const string &parname) {
2795 map<string, double>::iterator mapiter;
2796 mapiter = m_instrumentParmaeters.find(parname);
2797
2798 if (mapiter == m_instrumentParmaeters.end()) {
2799 stringstream errss;
2800 errss << "Instrument parameter map (having " << m_instrumentParmaeters.size() << " entries) "
2801 << "does not have parameter " << parname << ". ";
2802 g_log.debug(errss.str());
2803 return (EMPTY_DBL());
2804 }
2805
2806 return mapiter->second;
2807}
2808
2809//----------------------------------------------------------------------------------------------
2817 size_t workspaceindex, double leftbound, double rightbound) {
2818 // 1. Check
2819 const auto &X = sourcews->x(workspaceindex);
2820 const auto &Y = sourcews->y(workspaceindex);
2821 const auto &E = sourcews->e(workspaceindex);
2822
2823 if (leftbound >= rightbound) {
2824 stringstream errmsg;
2825 errmsg << "[BuildPartialWorkspace] Input left boundary = " << leftbound << " is larger than input right boundary "
2826 << rightbound << ". It is not allowed. ";
2827 throw std::invalid_argument(errmsg.str());
2828 }
2829 if (leftbound >= X.back() || rightbound <= X[0]) {
2830 throw std::invalid_argument("Boundary is out side of the input data set. ");
2831 }
2832
2833 // 2. Determine the size of the "partial" workspace
2834 int ileft = static_cast<int>(std::lower_bound(X.begin(), X.end(), leftbound) - X.begin());
2835 if (ileft > 0)
2836 --ileft;
2837 int iright = static_cast<int>(std::lower_bound(X.begin(), X.end(), rightbound) - X.begin());
2838 if (iright >= static_cast<int>(X.size()))
2839 iright = static_cast<int>(X.size() - 1);
2840
2841 auto wssize = static_cast<size_t>(iright - ileft + 1);
2842
2843 // 3. Build the partial workspace
2844 size_t nspec = 6;
2845 Workspace2D_sptr partws = std::dynamic_pointer_cast<Workspace2D>(
2846 API::WorkspaceFactory::Instance().create("Workspace2D", nspec, wssize, wssize));
2847
2848 // 4. Put data there
2849 // TODO: Try to use vector copier
2850 for (size_t iw = 0; iw < partws->getNumberHistograms(); ++iw) {
2851 MantidVec &nX = partws->dataX(iw);
2852 for (size_t i = 0; i < wssize; ++i) {
2853 nX[i] = X[i + ileft];
2854 }
2855 }
2856 auto &nY = partws->mutableY(0);
2857 auto &nE = partws->mutableE(0);
2858 for (size_t i = 0; i < wssize; ++i) {
2859 nY[i] = Y[i + ileft];
2860 nE[i] = E[i + ileft];
2861 }
2862
2863 return partws;
2864}
2865
2866//----------------------------------------------------------------------------------------------
2870string getFunctionInfo(const IFunction_sptr &function) {
2871 stringstream outss;
2872 vector<string> parnames = function->getParameterNames();
2873 size_t numpars = parnames.size();
2874 outss << "Number of Parameters = " << numpars << '\n';
2875 for (size_t i = 0; i < numpars; ++i)
2876 outss << parnames[i] << " = " << function->getParameter(i) << ", \t\tFitted = " << function->isActive(i) << '\n';
2877
2878 return outss.str();
2879}
2880
2881//----------------------------------------------------------------------------------------------
2894 size_t wsindexraw, size_t wsindexbkgd, size_t wsindexpeak) {
2895 // 1. Get prepared
2896 if (dataws->getNumberHistograms() < 3) {
2897 stringstream errss;
2898 errss << "Function estimateBackgroundCoase() requires input Workspace2D "
2899 "has at least 3 spectra."
2900 << "Present input has " << dataws->getNumberHistograms() << " spectra.";
2901 throw runtime_error(errss.str());
2902 }
2903 const auto &X = dataws->x(wsindexraw);
2904 const auto &Y = dataws->y(wsindexraw);
2905
2906 // TODO: This is a magic number!
2907 size_t numsamplepts = 2;
2908 if (X.size() <= 10) {
2909 // Make it at minimum to estimate background
2910 numsamplepts = 1;
2911 }
2912
2913 // 2. Average the first and last three data points
2914 double y0 = 0;
2915 double x0 = 0;
2916
2917 for (size_t i = 0; i < numsamplepts; ++i) {
2918 x0 += X[i];
2919 y0 += Y[i];
2920 }
2921 x0 = x0 / static_cast<double>(numsamplepts);
2922 y0 = y0 / static_cast<double>(numsamplepts);
2923
2924 double xf = 0;
2925 double yf = 0;
2926 for (size_t i = X.size() - numsamplepts; i < X.size(); ++i) {
2927 xf += X[i];
2928 yf += Y[i];
2929 }
2930 xf = xf / static_cast<double>(numsamplepts);
2931 yf = yf / static_cast<double>(numsamplepts);
2932
2933 // 3. Calculate B(x) = B0 + B1*x
2934 double b1 = (yf - y0) / (xf - x0);
2935 double b0 = yf - b1 * xf;
2936
2937 background->setParameter("A0", b0);
2938 background->setParameter("A1", b1);
2939
2940 // 4. Calcualte background
2941 FunctionDomain1DVector domain(X.rawData());
2942 FunctionValues values(domain);
2943 background->function(domain, values);
2944
2945 dataws->mutableY(wsindexbkgd) = values.toVector();
2946 dataws->mutableY(wsindexpeak) = dataws->y(wsindexraw) - dataws->y(wsindexbkgd);
2947 dataws->mutableE(wsindexpeak) = dataws->e(wsindexraw);
2948}
2949
2950//-----------------------------------------------------------------------------------------------------------
2964bool observePeakParameters(const Workspace2D_sptr &dataws, size_t wsindex, double &centre, double &height, double &fwhm,
2965 string &errmsg) {
2966 // 1. Get the value of the Max Height
2967 const auto &X = dataws->x(wsindex);
2968 const auto &Y = dataws->y(wsindex);
2969
2970 // 2. The highest peak should be the centre
2971 size_t icentre = findMaxValue(Y.rawData());
2972 centre = X[icentre];
2973 height = Y[icentre];
2974
2975 if (icentre <= 1 || icentre > X.size() - 2) {
2976 stringstream errss;
2977 errss << "Peak center = " << centre << " is at the edge of the input workspace [" << X[0] << ", " << X.back()
2978 << ". It is unable to proceed the estimate of FWHM. Quit with error!.";
2979 errmsg = errss.str();
2980 return false;
2981 }
2982 if (height <= 0) {
2983 stringstream errss;
2984 errss << "Max height = " << height << " in input workspace [" << X[0] << ", " << X.back()
2985 << " is negative. Fatal error is design of the algorithm.";
2986 errmsg = errss.str();
2987 return false;
2988 }
2989
2990 // 3. Calculate FWHM
2991 double halfMax = height * 0.5;
2992
2993 // a) Deal with left side
2994 bool continueloop = true;
2995 size_t index = icentre - 1;
2996 while (continueloop) {
2997 if (Y[index] <= halfMax) {
2998 // Located the data points
2999 continueloop = false;
3000 } else if (index == 0) {
3001 // Reach the end of the boundary, but haven't found. return with error.
3002 stringstream errss;
3003 errss << "The peak is not complete (left side) in the given data range.";
3004 errmsg = errss.str();
3005 return false;
3006 } else {
3007 // Contineu to locate
3008 --index;
3009 }
3010 }
3011 double x0 = X[index];
3012 double xf = X[index + 1];
3013 double y0 = Y[index];
3014 double yf = Y[index + 1];
3015
3016 // Formular for linear iterpolation: X = [(xf-x0)*Y - (xf*y0-x0*yf)]/(yf-y0)
3017 double xl = linearInterpolateX(x0, xf, y0, yf, halfMax);
3018
3019 double lefthalffwhm = centre - xl;
3020
3021 // 3. Deal with right side
3022 continueloop = true;
3023 index = icentre + 1;
3024 while (continueloop) {
3025 if (Y[index] <= halfMax) {
3026 // Located the data points
3027 continueloop = false;
3028 } else if (index == Y.size() - 1) {
3029 // Reach the end of the boundary, but haven't found. return with error.
3030 stringstream errss;
3031 errss << "The peak is not complete (right side) in the given data range.";
3032 errmsg = errss.str();
3033 return false;
3034 } else {
3035 ++index;
3036 }
3037 }
3038 x0 = X[index - 1];
3039 xf = X[index];
3040 y0 = Y[index - 1];
3041 yf = Y[index];
3042
3043 // Formula for linear iterpolation: X = [(xf-x0)*Y - (xf*y0-x0*yf)]/(yf-y0)
3044 double xr = linearInterpolateX(x0, xf, y0, yf, halfMax);
3045
3046 double righthalffwhm = xr - centre;
3047
3048 // Final
3049 fwhm = lefthalffwhm + righthalffwhm;
3050
3051 return true;
3052}
3053
3054//-----------------------------------------------------------------------------------------------------------
3059size_t findMaxValue(const std::vector<double> &Y) {
3060
3061 auto maxIt = std::max_element(Y.begin(), Y.end());
3062 return std::distance(Y.begin(), maxIt);
3063}
3064
3065//----------------------------------------------------------------------------------------------
3073size_t findMaxValue(const MatrixWorkspace_sptr &dataws, size_t wsindex, double leftbound, double rightbound) {
3074 const auto &X = dataws->x(wsindex);
3075 const auto &Y = dataws->y(wsindex);
3076
3077 // 1. Determine xmin, xmax range
3078 std::vector<double>::const_iterator viter;
3079
3080 viter = std::lower_bound(X.begin(), X.end(), leftbound);
3081 size_t ixmin = size_t(viter - X.begin());
3082 if (ixmin != 0)
3083 --ixmin;
3084 viter = std::lower_bound(X.begin(), X.end(), rightbound);
3085 size_t ixmax = size_t(viter - X.begin());
3086
3087 // 2. Search imax
3088 size_t imax = ixmin;
3089 double maxY = Y[ixmin];
3090 for (size_t i = ixmin + 1; i <= ixmax; ++i) {
3091 if (Y[i] > maxY) {
3092 maxY = Y[i];
3093 imax = i;
3094 }
3095 }
3096
3097 return imax;
3098}
3099
3100} // namespace Mantid::CurveFitting::Algorithms
#define DECLARE_ALGORITHM(classname)
Definition Algorithm.h:538
double value
The value of the point.
Definition FitMW.cpp:51
#define PEAKFITRANGEFACTOR
Factor on FWHM for fitting a peak.
#define PEAKBOUNDARYFACTOR
Factor on FWHM for defining a peak's range.
double sigma
Definition GetAllEi.cpp:156
double height
Definition GetAllEi.cpp:155
std::map< DeltaEMode::Type, std::string > index
#define UNUSED_ARG(x)
Function arguments are sometimes unused in certain implmentations but are required for documentation ...
Definition System.h:48
void declareProperty(std::unique_ptr< Kernel::Property > p, const std::string &doc="") override
Add a property to the list of managed properties.
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.
Kernel::Logger & g_log
Definition Algorithm.h:422
A composite function is a function containing other functions.
Implements FunctionDomain1D with its own storage in form of a std::vector.
size_t size() const override
Return the number of arguments in the domain.
A class to store values calculated by a function.
const std::vector< double > & toVector() const
Return the calculated values as a vector.
virtual void initialize()
Iinialize the function.
Definition IFunction.h:425
std::vector< std::string > getParameterNames() const
Return a vector with all parameter names.
virtual void setMillerIndex(int h, int k, int l)
Set Miller Indicies.
virtual double centre() const
Overwrite IFunction base class methods.
TableRow represents a row in a TableWorkspace.
Definition TableRow.h:39
size_t size() const
Returns the number of rows in the TableWorkspace.
Definition TableRow.h:45
A property class for workspaces.
FitPowderDiffPeaks : Fit peaks in powder diffraction pattern.
double m_rightmostPeakLeftBound
Right most peak's left boundary.
void init() override
Implement abstract Algorithm methods.
void restoreFunctionParameters(const API::IFunction_sptr &function, std::map< std::string, double > parammap)
Restore the function's parameter values from a map.
bool doFit1PeakSequential(const DataObjects::Workspace2D_sptr &dataws, size_t workspaceindex, const Functions::BackToBackExponential_sptr &peakfunction, std::vector< std::string > minimzernames, std::vector< size_t > maxiterations, const std::vector< double > &dampfactors, double &chi2)
Fit 1 peak by using a sequential of minimizer.
void observePeakRange(const Functions::BackToBackExponential_sptr &thispeak, const Functions::BackToBackExponential_sptr &rightpeak, double refpeakshift, double &peakleftbound, double &peakrightbound)
Observe peak range with hint from right peak's properties.
void storeFunctionParameters(const API::IFunction_sptr &function, std::map< std::string, double > &parammaps)
Store the function's parameter values to a map.
double getParameter(const std::string &parname)
Get parameter value from m_instrumentParameters.
bool doFitMultiplePeaks(const DataObjects::Workspace2D_sptr &dataws, size_t wsindex, const API::CompositeFunction_sptr &peaksfunc, std::vector< Functions::BackToBackExponential_sptr > peakfuncs, std::vector< bool > &vecfitgood, std::vector< double > &vecchi2s)
Fit multiple (overlapped) peaks.
bool doFitGaussianPeak(const DataObjects::Workspace2D_sptr &dataws, size_t workspaceindex, double in_center, double leftfwhm, double rightfwhm, double &center, double &sigma, double &height)
Fit background-removed peak by Gaussian.
std::vector< int > m_rightmostPeakHKL
Right most peak HKL.
std::string parseFitResult(const API::IAlgorithm_sptr &fitalg, double &chi2, bool &fitsuccess)
Parse the fitting result.
enum Mantid::CurveFitting::Algorithms::FitPowderDiffPeaks::@0 m_fitMode
Fit mode.
std::pair< DataObjects::TableWorkspace_sptr, DataObjects::TableWorkspace_sptr > genPeakParametersWorkspace()
Generate output peak parameters workspace.
DataObjects::Workspace2D_sptr genOutputFittedPatternWorkspace(const std::vector< double > &pattern, int workspaceindex)
Create a Workspace2D for fitted peaks (pattern)
bool getHKLFromMap(std::map< std::string, int > intmap, std::vector< int > &hkl)
Get (HKL) from a map; Return false if the information is incomplete.
double m_rightmostPeakRightBound
Right most peak's right boundary.
std::vector< bool > m_goodFit
Peak fitting status.
std::vector< double > m_peakFitChi2
Peak fitting information.
double m_minPeakHeight
Minimum peak height for peak to be refined.
bool doFit1PeakSimple(const DataObjects::Workspace2D_sptr &dataws, size_t workspaceindex, const Functions::BackToBackExponential_sptr &peakfunction, const std::string &minimzername, size_t maxiteration, double &chi2)
Fit 1 peak by 1 minimizer of 1 call of minimzer (simple version)
void exec() override
Implement abstract Algorithm methods.
bool m_fitPeakBackgroundComposite
Fit peak + background as the last step.
void importInstrumentParameterFromTable(const DataObjects::TableWorkspace_sptr &parameterWS)
Import instrument parameters from (input) table workspace.
void setOverlappedPeaksConstraints(const std::vector< Functions::BackToBackExponential_sptr > &peaks)
Set constraints on a group of overlapped peaks for fitting.
DataObjects::TableWorkspace_sptr m_profileTable
Instrument profile parameter table workspace.
void genPeaksFromTable(const DataObjects::TableWorkspace_sptr &peakparamws)
Generate peaks from input table workspace.
void calculatePeakFitBoundary(size_t ileftpeak, size_t irightpeak, double &peakleftboundary, double &peakrightboundary)
Calculate the range to fit peak/peaks group.
Geometry::UnitCell m_unitCell
Unit cell of powder crystal.
std::pair< bool, double > doFitPeak(const DataObjects::Workspace2D_sptr &dataws, const Functions::BackToBackExponential_sptr &peakfunction, double guessedfwhm)
Fit peak with flexiblity in multiple steps Prerequisit:
std::vector< double > m_peakData
Data for each individual peaks. (HKL)^2, vector index, function values.
void parseBraggPeakTable(const DataObjects::TableWorkspace_sptr &peakws, std::vector< std::map< std::string, double > > &parammaps, std::vector< std::map< std::string, int > > &hklmaps)
Import Bragg peak table workspace.
void cropWorkspace(double tofmin, double tofmax)
Crop data workspace.
enum Mantid::CurveFitting::Algorithms::FitPowderDiffPeaks::@1 m_genPeakStartingValue
Choice to generate peak profile paramter starting value.
DataObjects::Workspace2D_sptr genPeakParameterDataWorkspace()
Create data workspace for X0, A, B and S of peak with good fit.
DataObjects::TableWorkspace_sptr m_peakParamTable
Bragg peak parameter.
Functions::BackToBackExponential_sptr genPeak(std::map< std::string, int > hklmap, std::map< std::string, double > parammap, std::map< std::string, std::string > bk2bk2braggmap, bool &good, std::vector< int > &hkl, double &d_h)
Generate a peak.
std::vector< std::string > mPeakParameterNames
Peak parmeter names.
bool fitOverlappedPeaks(std::vector< Functions::BackToBackExponential_sptr > peaks, const Functions::BackgroundFunction_sptr &backgroundfunction, double gfwhm)
Fit peaks with confidence in fwhm and etc.
std::map< std::string, double > m_instrumentParmaeters
Map for function (instrument parameter)
std::vector< std::pair< double, std::pair< std::vector< int >, Functions::BackToBackExponential_sptr > > > m_vecPeakFunctions
Sorted vector for peaks. double = d_h, vector = (HKL), peak.
int m_numPeaksLowerToMin
Number of peaks to fit lower to minimum HKL.
bool doFit1PeakBackground(const DataObjects::Workspace2D_sptr &dataws, size_t wsindex, const Functions::BackToBackExponential_sptr &peak, const Functions::BackgroundFunction_sptr &backgroundfunction, double &chi2)
Fit 1 peak and background.
void plotFunction(const API::IFunction_sptr &peakfunction, const Functions::BackgroundFunction_sptr &background, const API::FunctionDomain1DVector &domain)
Plot a single peak to output vector.
std::string parseFitParameterWorkspace(const API::ITableWorkspace_sptr &paramws)
Parse Fit() output parameter workspace.
DataObjects::Workspace2D_sptr buildPartialWorkspace(const API::MatrixWorkspace_sptr &sourcews, size_t workspaceindex, double leftbound, double rightbound)
Build partial workspace for fitting.
bool m_useGivenTOFh
Flag to use given Bragg peaks' centre in TOF.
bool fitSinglePeakRobust(const Functions::BackToBackExponential_sptr &peak, const Functions::BackgroundFunction_sptr &backgroundfunction, double peakleftbound, double peakrightbound, const std::map< std::string, double > &rightpeakparammap, double &finalchi2)
Fit single peak in robust mode (no hint)
void estimatePeakHeightsLeBail(const DataObjects::Workspace2D_sptr &dataws, size_t wsindex, std::vector< Functions::BackToBackExponential_sptr > peaks)
Use Le Bail method to estimate and set the peak heights.
bool fitSinglePeakConfident(const Functions::BackToBackExponential_sptr &peak, const Functions::BackgroundFunction_sptr &backgroundfunction, double leftbound, double rightbound, double &chi2, bool &annhilatedpeak)
Fit peak with trustful peak parameters.
int m_wsIndex
TOF vector of data workspace to process with.
bool doFitNPeaksSimple(const DataObjects::Workspace2D_sptr &dataws, size_t wsindex, const API::CompositeFunction_sptr &peaksfunc, const std::vector< Functions::BackToBackExponential_sptr > &peakfuncs, const std::string &minimizername, size_t maxiteration, double &chi2)
Fit N overlapped peaks in a simple manner.
bool fitSinglePeakSimulatedAnnealing(const Functions::BackToBackExponential_sptr &peak, const std::vector< std::string > &paramtodomc)
Fit signle peak by Monte Carlo/simulated annealing.
Provide BackToBackExponential peak shape function interface to IPeakFunction.
ThermalNeutronBk2BkExpConvPVoigt : Back-to-back exponential convoluted with pseudo Voigt for thermal ...
double getPeakParameter(const std::string &) override
Overwrite IPeakFunction base class methods.
void setParameter(size_t i, const double &value, bool explicitlySet=true) override
Core function to calcualte peak values for whole region.
void calculateParameters(bool explicitoutput) const override
Calculate peak parameters (alpha, beta, sigma2..)
TableWorkspace is an implementation of Workspace in which the data are organised in columns of same s...
void set(double _a, double _b, double _c, double _alpha, double _beta, double _gamma, const int angleunit=angDegrees)
Set lattice parameters.
Definition UnitCell.cpp:303
double d(double h, double k, double l) const
Return d-spacing ( ) for a given h,k,l coordinate.
Definition UnitCell.cpp:700
IPropertyManager * setProperty(const std::string &name, const T &value)
Templated method to set the value of a PropertyWithValue.
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 information(const std::string &msg)
Logs at information level.
Definition Logger.cpp:136
static T & Instance()
Return a reference to the Singleton instance, creating it if it does not already exist Creation is do...
std::shared_ptr< IAlgorithm > IAlgorithm_sptr
shared pointer to Mantid::API::IAlgorithm
std::shared_ptr< Column > Column_sptr
Definition Column.h:232
std::shared_ptr< ITableWorkspace > ITableWorkspace_sptr
shared pointer to Mantid::API::ITableWorkspace
std::shared_ptr< IFunction > IFunction_sptr
shared pointer to the function base class
Definition IFunction.h:743
std::shared_ptr< Algorithm > Algorithm_sptr
Typedef for a shared pointer to an Algorithm.
Definition Algorithm.h:52
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
double linearInterpolateX(double x0, double xf, double y0, double yf, double y)
Formular for linear iterpolation: X = [(xf-x0)*Y - (xf*y0-x0*yf)]/(yf-y0)
size_t findMaxValue(const std::vector< double > &Y)
Find maximum value.
bool observePeakParameters(const DataObjects::Workspace2D_sptr &dataws, size_t wsindex, double &centre, double &height, double &fwhm, std::string &errmsg)
Estimate peak parameters;.
void estimateBackgroundCoarse(const DataObjects::Workspace2D_sptr &dataws, const Functions::BackgroundFunction_sptr &background, size_t wsindexraw, size_t wsindexbkgd, size_t wsindexpeak)
Estimate background for a pattern in a coarse mode.
std::string getFunctionInfo(const API::IFunction_sptr &function)
Get function parameter name, value and etc information in string.
double calThermalNeutronTOF(double dh, double dtt1, double dtt1t, double dtt2t, double zero, double zerot, double width, double tcross)
Calcualte TOF from d-spacing value for thermal neutron.
std::shared_ptr< Polynomial > Polynomial_sptr
Definition Polynomial.h:53
std::shared_ptr< BackToBackExponential > BackToBackExponential_sptr
std::shared_ptr< BackgroundFunction > BackgroundFunction_sptr
std::shared_ptr< Workspace2D > Workspace2D_sptr
shared pointer to Mantid::DataObjects::Workspace2D
std::shared_ptr< TableWorkspace > TableWorkspace_sptr
shared pointer to Mantid::DataObjects::TableWorkspace
std::unique_ptr< T > create(const P &parent, const IndexArg &indexArg, const HistArg &histArg)
This is the create() method that all the other create() methods call.
std::vector< double > getZscore(const std::vector< TYPE > &data)
Return the Z score values for a dataset.
Helper class which provides the Collimation Length for SANS instruments.
std::vector< double > MantidVec
typedef for the data storage used in Mantid matrix workspaces
Definition cow_ptr.h:172
constexpr double EMPTY_DBL() noexcept
Returns what we consider an "empty" double within a property.
Definition EmptyValues.h:42
STL namespace.
@ InOut
Both an input & output workspace.
Definition Property.h:55
@ Input
An input workspace.
Definition Property.h:53
@ Output
An output workspace.
Definition Property.h:54