Mantid
Loading...
Searching...
No Matches
PDCalibration.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 +
14#include "MantidAPI/Run.h"
16#include "MantidAPI/TableRow.h"
27#include "MantidHistogramData/Histogram.h"
36#include "MantidKernel/Unit.h"
37
38#include <algorithm>
39#include <gsl/gsl_multifit_nlin.h>
40#include <gsl/gsl_multimin.h>
41#include <limits>
42#include <numeric>
43
44namespace Mantid::Algorithms {
45
46using namespace Mantid::DataObjects;
47using namespace Mantid::HistogramData;
69using std::vector;
70
71// Register the algorithm into the AlgorithmFactory
72DECLARE_ALGORITHM(PDCalibration)
73
74namespace { // anonymous
75const auto isNonZero = [](const double value) { return value != 0.; };
76} // namespace
77
80public:
90 FittedPeaks(const API::MatrixWorkspace_const_sptr &wksp, const std::size_t wkspIndex) {
91 this->wkspIndex = wkspIndex;
92
93 // convert workspace index into detector id
94 const auto &spectrum = wksp->getSpectrum(wkspIndex);
95 this->detid = spectrum.getDetectorIDs();
96
97 const auto &X = spectrum.x();
98 const auto &Y = spectrum.y();
99 tofMin = X.front();
100 tofMax = X.back();
101
102 // determine tof min supported by the workspace (non-zero counts)
103 size_t minIndex = 0; // want to store value
104 for (; minIndex < Y.size(); ++minIndex) {
105 if (isNonZero(Y[minIndex])) {
106 tofMin = X[minIndex];
107 break; // we found the first bin with non-zero counts
108 }
109 }
110
111 // determine tof max supported by the workspace (non-zero counts)
112 size_t maxIndex = Y.size() - 1;
113 for (; maxIndex > minIndex; --maxIndex) {
114 if (isNonZero(Y[maxIndex])) {
115 tofMax = X[maxIndex];
116 break; // we found the last bin with non-zero counts
117 }
118 }
119 }
120
131 void setPositions(const std::vector<double> &peaksInD, const std::vector<double> &peaksInDWindows, const double difa,
132 const double difc, const double tzero) {
133 // clear out old values
134 inDPos.clear();
135 inTofPos.clear();
136 inTofWindows.clear();
137
138 // assign things
139 inDPos.assign(peaksInD.begin(), peaksInD.end());
140 inTofPos.assign(peaksInD.begin(), peaksInD.end());
141 inTofWindows.assign(peaksInDWindows.begin(), peaksInDWindows.end());
142
143 // convert the bits that matter to TOF
144 Kernel::Units::dSpacing dSpacingUnit;
145 std::vector<double> yunused;
146 dSpacingUnit.toTOF(
147 inTofPos, yunused, -1, 0,
149 dSpacingUnit.toTOF(
150 inTofWindows, yunused, -1, 0,
152 }
153
154 std::size_t wkspIndex;
155 std::set<detid_t> detid;
156 double tofMin; // TOF of bin with minimum TOF and non-zero counts
157 double tofMax; // TOF of bin with maximum TOF and non-zero counts
158 std::vector<double> inTofPos; // peak centers, in TOF
159 // left and right fit ranges for each peak center, in TOF
160 std::vector<double> inTofWindows;
161 std::vector<double> inDPos; // peak centers, in d-spacing
162};
163
164//----------------------------------------------------------------------------------------------
168
169//----------------------------------------------------------------------------------------------
173
174//----------------------------------------------------------------------------------------------
175
177const std::string PDCalibration::name() const { return "PDCalibration"; }
178
180int PDCalibration::version() const { return 1; }
181
183const std::string PDCalibration::category() const { return "Diffraction\\Calibration"; }
184
186const std::string PDCalibration::summary() const {
187 return "This algorithm determines the diffractometer constants that convert diffraction spectra "
188 "from time-of-flight (TOF) to d-spacing. The results are output to a calibration table.";
189}
190
191//----------------------------------------------------------------------------------------------
195 declareProperty(std::make_unique<WorkspaceProperty<MatrixWorkspace>>("InputWorkspace", "", Direction::Input),
196 "Input workspace containing spectra as a function of TOF measured on a standard sample.");
197
198 auto mustBePositive = std::make_shared<BoundedValidator<int>>();
199 mustBePositive->setLower(0);
200 declareProperty("StartWorkspaceIndex", 0, mustBePositive, "Starting workspace index for fit");
202 "StopWorkspaceIndex", EMPTY_INT(),
203 "Last workspace index for fit is the smaller of this value and the workspace index of last spectrum.");
204
206 std::make_unique<ArrayProperty<double>>("TofBinning", std::make_shared<RebinParamsValidator>()),
207 "Min, Step, and Max of TOF bins. "
208 "Logarithmic binning is used if Step is negative. Chosen binning should ensure sufficient datapoints across "
209 "the peaks to be fitted, considering the number of parameters required by PeakFunction and BackgroundType.");
210
211 const std::vector<std::string> exts2{".h5", ".cal"};
212 declareProperty(std::make_unique<FileProperty>("PreviousCalibrationFile", "", FileProperty::OptionalLoad, exts2),
213 "An existent file to be loaded into a calibration table, which is used for converting PeakPositions "
214 "from d-spacing to TOF.");
216 "PreviousCalibrationTable", "", Direction::Input, API::PropertyMode::Optional),
217 "An existent calibration table used for converting PeakPositions from d-spacing to TOF. "
218 "This property has precedence over PreviousCalibrationFile.");
219
220 // properties about peak positions to fit
221 std::vector<std::string> peaktypes{"BackToBackExponential", "Gaussian", "Lorentzian", "PseudoVoigt",
222 "IkedaCarpenterPV"};
223 declareProperty("PeakFunction", "Gaussian", std::make_shared<StringListValidator>(peaktypes),
224 "Function to fit input peaks.");
225 vector<std::string> bkgdtypes{"Flat", "Linear", "Quadratic"};
226 declareProperty("BackgroundType", "Linear", std::make_shared<StringListValidator>(bkgdtypes),
227 "Function to fit input peaks background.");
228
229 auto peaksValidator = std::make_shared<CompositeValidator>();
230 auto mustBePosArr = std::make_shared<Kernel::ArrayBoundedValidator<double>>();
231 mustBePosArr->setLower(0.0);
232 peaksValidator->add(mustBePosArr);
233 peaksValidator->add(std::make_shared<MandatoryValidator<std::vector<double>>>());
234 declareProperty(std::make_unique<ArrayProperty<double>>("PeakPositions", peaksValidator),
235 "Comma-delimited positions (d-spacing) of reference peaks. Care should be taken to avoid using peaks "
236 "whose predicted positions are too close considering their peak windows.");
237
238 auto windowValidator = std::make_shared<CompositeValidator>();
239 windowValidator->add(mustBePosArr);
240 auto lengthValidator = std::make_shared<Kernel::ArrayLengthValidator<double>>();
241 lengthValidator->setLengthMin(1);
242 windowValidator->add(lengthValidator);
243 windowValidator->add(std::make_shared<MandatoryValidator<std::vector<double>>>());
244 declareProperty(std::make_unique<ArrayProperty<double>>("PeakWindow", "0.1", windowValidator),
245 "Width of the window (d-spacing) over which to fit a peak. If a single value is supplied, it will be "
246 "used as half "
247 "the window width for all peaks. Otherwise, the expected input is a comma-delimited list of 2*N "
248 "window boundaries, "
249 "where N is the number of values in PeakPositions. The order of the window boundaries should match "
250 "the order in PeakPositions.");
251
252 auto min = std::make_shared<BoundedValidator<double>>();
253 min->setLower(1e-3);
254 declareProperty("PeakWidthPercent", EMPTY_DBL(), min,
255 "Used for estimating peak width (an initial parameter for peak fitting) by multiplying peak's value "
256 "in PeakPositions by this factor.");
257
258 declareProperty("MinimumPeakHeight", 2.,
259 "Used for validating peaks before and after fitting. If a peak's observed/estimated or "
260 "fitted height is under this value, the peak will be marked as an error.");
261
262 declareProperty("MaxChiSq", 100.,
263 "Used for validating peaks after fitting. If the chi-squared value is higher than this value, "
264 "the peak will be excluded from calibration. The recommended value is between 2 and 10.");
265
266 declareProperty("ConstrainPeakPositions", false,
267 "If true, a peak center being fit will be constrained by the estimated peak position "
268 "+/- 0.5 * \"estimated peak width\", where the estimated peak position is the highest Y-value "
269 "position in the window, "
270 " and the estimated peak width is FWHM calculated over the window data.");
271
272 declareProperty("HighBackground", false,
273 "Flag whether the input data has high background compared to peaks' heights. "
274 "This option is recommended for data with peak-to-background ratios under ~5.");
275
276 declareProperty("MinimumSignalToNoiseRatio", 0.,
277 "Used for validating peaks before fitting. If the signal-to-noise ratio is under this value, "
278 "the peak will be excluded from fitting and calibration. This check does not apply to peaks for "
279 "which the noise cannot be estimated. "
280 "The minimum recommended value is 12.");
281
282 std::vector<std::string> modes{"DIFC", "DIFC+TZERO", "DIFC+TZERO+DIFA"};
283 declareProperty("CalibrationParameters", "DIFC", std::make_shared<StringListValidator>(modes),
284 "Select which diffractometer constants (GSAS convention) to determine.");
285
286 declareProperty(std::make_unique<ArrayProperty<double>>("TZEROrange"),
287 "Range for allowable calibrated TZERO value. Default: no restriction.");
288
289 declareProperty(std::make_unique<ArrayProperty<double>>("DIFArange"),
290 "Range for allowable calibrated DIFA value. Default: no restriction.");
291
292 declareProperty("UseChiSq", false,
293 "Defines the weighting scheme used in the least-squares fit of the extracted peak centers "
294 "that determines the diffractometer constants. If true, the peak weight will be "
295 "the inverse square of the error on the fitted peak center. If false, the peak weight will be "
296 "the square of the fitted peak height.");
297
299 std::make_unique<WorkspaceProperty<API::ITableWorkspace>>("OutputCalibrationTable", "", Direction::Output),
300 "Output table workspace containing the calibration.");
301
302 // Mantid's python API _requires_ a non empty-string name for any Output workspace, even when 'PropertyMode::Optional'
303 // is specified.
304 declareProperty(std::make_unique<WorkspaceProperty<MaskWorkspace>>("MaskWorkspace", "_empty_", Direction::Output,
306 "Mask workspace (optional input / output workspace):"
307 " when specified, if the workspace already exists, any incoming masked detectors will be combined"
308 " with any additional outgoing masked detectors detected by the algorithm");
309
311 std::make_unique<WorkspaceProperty<API::WorkspaceGroup>>("DiagnosticWorkspaces", "", Direction::Output),
312 "Auxiliary workspaces containing extended information on the calibration results.");
313
314 declareProperty("MinimumPeakTotalCount", EMPTY_DBL(),
315 "Used for validating peaks before fitting. If the total peak window Y-value count "
316 "is under this value, the peak will be excluded from fitting and calibration.");
317
318 declareProperty("MinimumSignalToSigmaRatio", 0.,
319 "Used for validating peaks after fitting. If the signal-to-sigma ratio is under this value, "
320 "the peak will be excluded from fitting and calibration.");
321
322 // make group for Input properties
323 std::string inputGroup("Input Options");
324 setPropertyGroup("InputWorkspace", inputGroup);
325 setPropertyGroup("StartWorkspaceIndex", inputGroup);
326 setPropertyGroup("StopWorkspaceIndex", inputGroup);
327 setPropertyGroup("TofBinning", inputGroup);
328 setPropertyGroup("PreviousCalibrationFile", inputGroup);
329 setPropertyGroup("PreviousCalibrationTable", inputGroup);
330
331 std::string funcgroup("Function Types");
332 setPropertyGroup("PeakFunction", funcgroup);
333 setPropertyGroup("BackgroundType", funcgroup);
334
335 // make group for FitPeaks properties
336 std::string fitPeaksGroup("Peak Fitting");
337 setPropertyGroup("PeakPositions", fitPeaksGroup);
338 setPropertyGroup("PeakWindow", fitPeaksGroup);
339 setPropertyGroup("PeakWidthPercent", fitPeaksGroup);
340 setPropertyGroup("MinimumPeakHeight", fitPeaksGroup);
341 setPropertyGroup("MinimumSignalToNoiseRatio", fitPeaksGroup);
342 setPropertyGroup("MinimumPeakTotalCount", fitPeaksGroup);
343 setPropertyGroup("MinimumSignalToSigmaRatio", fitPeaksGroup);
344 setPropertyGroup("HighBackground", fitPeaksGroup);
345 setPropertyGroup("MaxChiSq", fitPeaksGroup);
346 setPropertyGroup("ConstrainPeakPositions", fitPeaksGroup);
347
348 // make group for type of calibration
349 std::string calGroup("Calibration Type");
350 setPropertyGroup("CalibrationParameters", calGroup);
351 setPropertyGroup("TZEROrange", calGroup);
352 setPropertyGroup("DIFArange", calGroup);
353 setPropertyGroup("UseChiSq", calGroup);
354}
355
356std::map<std::string, std::string> PDCalibration::validateInputs() {
357 std::map<std::string, std::string> messages;
358
359 if (MaskWorkspace_const_sptr maskWS = getProperty("MaskWorkspace")) {
360 MatrixWorkspace_const_sptr inputWS = getProperty("InputWorkspace");
361
362 // detectors which are monitors are not included in the mask
363 if (maskWS->getInstrument()->getNumberDetectors(true) != inputWS->getInstrument()->getNumberDetectors(true)) {
364 messages["MaskWorkspace"] = "incoming mask workspace must have the same instrument as the input workspace";
365 } else if (maskWS->getNumberHistograms() != inputWS->getInstrument()->getNumberDetectors(true)) {
366 messages["MaskWorkspace"] = "incoming mask workspace must have one spectrum per detector";
367 }
368 }
369
370 vector<double> tzeroRange = getProperty("TZEROrange");
371 if (!tzeroRange.empty()) {
372 if (tzeroRange.size() != 2) {
373 messages["TZEROrange"] = "Require two values [min,max]";
374 } else if (tzeroRange[0] >= tzeroRange[1]) {
375 messages["TZEROrange"] = "min must be less than max";
376 }
377 }
378
379 vector<double> difaRange = getProperty("DIFArange");
380 if (!difaRange.empty()) {
381 if (difaRange.size() != 2) {
382 messages["DIFArange"] = "Require two values [min,max]";
383 } else if (difaRange[0] >= difaRange[1]) {
384 messages["DIFArange"] = "min must be less than max";
385 }
386 }
387
388 vector<double> peakWindow = getProperty("PeakWindow");
389 vector<double> peakCentres = getProperty("PeakPositions");
390 if (peakWindow.size() > 1 && peakWindow.size() != 2 * peakCentres.size()) {
391 messages["PeakWindow"] = "PeakWindow must be a vector with either 1 element (interpreted as half the width of "
392 "the window) or twice the number of peak centres provided.";
393 }
394
395 return messages;
396}
397
398namespace {
399
401bool hasDasIDs(const API::ITableWorkspace_const_sptr &table) {
402 const auto columnNames = table->getColumnNames();
403 return (std::find(columnNames.begin(), columnNames.end(), std::string("dasid")) != columnNames.end());
404}
405
413double getWidthToFWHM(const std::string &peakshape) {
414 // could we use actual peak function here?
415 if (peakshape == "Gaussian") {
416 return 2 * std::sqrt(2. * std::log(2.));
417 } else if (peakshape == "Lorentzian") {
418 return 2.;
419 } else if (peakshape == "BackToBackExponential") {
420 return 1.; // TODO the conversion isn't document in the function
421 } else {
422 return 1.;
423 }
424}
425
426} // end of anonymous namespace
427
428//----------------------------------------------------------------------------------------------
432 vector<double> tofBinningParams = getProperty("TofBinning");
433 m_tofMin = tofBinningParams.front();
434 m_tofMax = tofBinningParams.back();
435
436 vector<double> tzeroRange = getProperty("TZEROrange");
437 if (tzeroRange.size() == 2) {
438 m_tzeroMin = tzeroRange[0];
439 m_tzeroMax = tzeroRange[1];
440
441 std::stringstream msg;
442 msg << "Using tzero range of " << m_tzeroMin << " <= "
443 << "TZERO <= " << m_tzeroMax;
444 g_log.information(msg.str());
445 } else {
446 g_log.information("Using all TZERO values");
447
448 m_tzeroMin = std::numeric_limits<double>::lowest();
449 m_tzeroMax = std::numeric_limits<double>::max();
450 }
451
452 vector<double> difaRange = getProperty("DIFArange");
453 if (difaRange.size() == 2) {
454 m_difaMin = difaRange[0];
455 m_difaMax = difaRange[1];
456
457 std::stringstream msg;
458 msg << "Using difa range of " << m_difaMin << " <= "
459 << "DIFA <= " << m_difaMax;
460 g_log.information(msg.str());
461 } else {
462 g_log.information("Using all DIFA values");
463
464 m_difaMin = std::numeric_limits<double>::lowest();
465 m_difaMax = std::numeric_limits<double>::max();
466 }
467
468 m_peaksInDspacing = getProperty("PeakPositions");
469 std::vector<double> peakWindow = getProperty("PeakWindow");
470 const std::size_t NUMPEAKS = m_peaksInDspacing.size();
471 // if peak windows given for each peak center, sort all by peak center
472 if (peakWindow.size() > 1) {
473 // load windows into ordered map keyed by peak center
474 // this preserves association of centers and window edges
475 std::map<double, std::pair<double, double>> peakEdgeAndCenter;
476 for (std::size_t i = 0; i < m_peaksInDspacing.size(); i++) {
477 peakEdgeAndCenter[m_peaksInDspacing[i]] = {peakWindow[2 * i], peakWindow[2 * i + 1]};
478 }
479 m_peaksInDspacing.clear();
480 m_peaksInDspacing.reserve(NUMPEAKS);
481 peakWindow.clear();
482 peakWindow.reserve(2 * NUMPEAKS);
483 // retrieve ordered center, windows
484 for (auto it = peakEdgeAndCenter.begin(); it != peakEdgeAndCenter.end(); it++) {
485 m_peaksInDspacing.push_back(it->first);
486 peakWindow.push_back(it->second.first);
487 peakWindow.push_back(it->second.second);
488 }
489 } else {
490 // Sort peak positions, required for correct peak window calculations
491 std::sort(m_peaksInDspacing.begin(), m_peaksInDspacing.end());
492 }
493
494 const double minPeakHeight = getProperty("MinimumPeakHeight");
495 const double minPeakTotalCount = getProperty("MinimumPeakTotalCount");
496 const double minSignalToNoiseRatio = getProperty("MinimumSignalToNoiseRatio");
497 const double minSignalToSigmaRatio = getProperty("MinimumSignalToSigmaRatio");
498 const double maxChiSquared = getProperty("MaxChiSq");
499
500 const std::string calParams = getPropertyValue("CalibrationParameters");
501 if (calParams == std::string("DIFC"))
503 else if (calParams == std::string("DIFC+TZERO"))
505 else if (calParams == std::string("DIFC+TZERO+DIFA"))
507 else
508 throw std::runtime_error("Encountered impossible CalibrationParameters value");
509
511 setProperty("InputWorkspace", m_uncalibratedWS);
512
513 m_startWorkspaceIndex = getProperty("StartWorkspaceIndex");
514 m_stopWorkspaceIndex = isDefault("StopWorkspaceIndex") ? static_cast<int>(m_uncalibratedWS->getNumberHistograms() - 1)
515 : getProperty("StopWorkspaceIndex");
516
517 auto uncalibratedEWS = std::dynamic_pointer_cast<EventWorkspace>(m_uncalibratedWS);
518 auto isEvent = bool(uncalibratedEWS);
519
520 // Load Previous Calibration or create calibration table from signal file
521 if ((!static_cast<std::string>(getProperty("PreviousCalibrationFile")).empty()) ||
522 (!getPropertyValue("PreviousCalibrationTable").empty())) { //"PreviousCalibrationTable"
524 } else {
525 createCalTableNew(); // calculates "difc" values from instrument geometry
526 }
528
529 // Use the incoming mask workspace, or start a new one if the workspace does not exist.
530 MaskWorkspace_sptr maskWS;
531 if (!isDefault("MaskWorkspace")) {
532 maskWS = getProperty("MaskWorkspace");
533 }
534 if (!maskWS) {
535 g_log.debug() << "[PDCalibration]: CREATING new MaskWorkspace.\n";
536 // A new mask is completely cleared at creation.
537 maskWS = std::make_shared<MaskWorkspace>(m_uncalibratedWS->getInstrument());
538 } else {
539 g_log.debug() << "[PDCalibration]: Using EXISTING MaskWorkspace.\n";
540 }
541 // Include any incoming masked detector flags in the mask-workspace values.
542 maskWS->combineFromDetectorMasks(m_uncalibratedWS->detectorInfo());
543
544 const std::string peakFunction = getProperty("PeakFunction");
545 const double WIDTH_TO_FWHM = getWidthToFWHM(peakFunction);
546 if (WIDTH_TO_FWHM == 1.) {
547 g_log.notice() << "Unknown conversion for \"" << peakFunction
548 << "\", found peak widths and resolution should not be "
549 "directly compared to delta-d/d";
550 }
551 auto NUMHIST = m_stopWorkspaceIndex - m_startWorkspaceIndex + 1;
552
553 // Create a pair of workspaces, one containing the nominal peak centers in TOF units,
554 // the other containing the left and right fitting ranges around each nominal
555 // peak center, also in TOF units. This for each pixel of the instrument
556 auto matrix_pair = createTOFPeakCenterFitWindowWorkspaces(m_uncalibratedWS, peakWindow);
557 API::MatrixWorkspace_sptr tof_peak_center_ws = matrix_pair.first;
558 API::MatrixWorkspace_sptr tof_peak_window_ws = matrix_pair.second;
559
560 double peak_width_percent = getProperty("PeakWidthPercent");
561
562 const std::string diagnostic_prefix = getPropertyValue("DiagnosticWorkspaces");
563
564 // Refine the position of the peak centers starting from the nominal peak
565 // centers and fitting them against a peak fit function (e.g. a Gaussian)
566 auto algFitPeaks = createChildAlgorithm("FitPeaks", .2, .7);
567 algFitPeaks->setLoggingOffset(3);
568
569 algFitPeaks->setProperty("InputWorkspace", m_uncalibratedWS);
570
571 // limit the spectra to fit
572 algFitPeaks->setProperty("StartWorkspaceIndex", static_cast<int>(m_startWorkspaceIndex));
573 algFitPeaks->setProperty("StopWorkspaceIndex", static_cast<int>(m_stopWorkspaceIndex));
574
575 // theoretical peak center
576 algFitPeaks->setProperty("PeakCentersWorkspace", tof_peak_center_ws);
577
578 // peak and background functions
579 algFitPeaks->setProperty<std::string>("PeakFunction", peakFunction);
580 algFitPeaks->setProperty<std::string>("BackgroundType", getProperty("BackgroundType"));
581 // peak range setup
582 algFitPeaks->setProperty("FitPeakWindowWorkspace", tof_peak_window_ws);
583 algFitPeaks->setProperty("PeakWidthPercent", peak_width_percent);
584 algFitPeaks->setProperty("MinimumPeakHeight", minPeakHeight);
585 algFitPeaks->setProperty("MinimumPeakTotalCount", minPeakTotalCount);
586 algFitPeaks->setProperty("MinimumSignalToNoiseRatio", minSignalToNoiseRatio);
587 algFitPeaks->setProperty("MinimumSignalToSigmaRatio", minSignalToSigmaRatio);
588 // some fitting strategy
589 algFitPeaks->setProperty("FitFromRight", true);
590 const bool highBackground = getProperty("HighBackground");
591 algFitPeaks->setProperty("HighBackground", highBackground);
592 bool constrainPeakPosition = getProperty("ConstrainPeakPositions");
593 algFitPeaks->setProperty("ConstrainPeakPositions",
594 constrainPeakPosition); // TODO Pete: need to test this option
595 // optimization setup // TODO : need to test LM or LM-MD
596 algFitPeaks->setProperty("Minimizer", "Levenberg-Marquardt");
597 algFitPeaks->setProperty("CostFunction", "Least squares");
598
599 // FitPeaks will abstract the peak parameters if you ask (if using chisq then
600 // need FitPeaks to output fitted params rather than height, width)
601 const bool useChiSq = getProperty("UseChiSq");
602 algFitPeaks->setProperty("RawPeakParameters", useChiSq);
603
604 // Analysis output
605 // If using a Gaussian peak shape plus a constant background, then
606 // OutputPeakParametersWorkspace is a table with columns:
607 // wsindex_0 peakindex_0 centre width height intensity A0 chi2
608 // ...
609 // wsindex_0 peakindex_N centre width height intensity A0 chi2
610 // ...
611 // ...
612 // wsindex_M peakindex_N centre width height intensity
613 algFitPeaks->setPropertyValue("OutputPeakParametersWorkspace", diagnostic_prefix + "_fitparam");
614 // contains the same intensities as input m_uncalibratedWS except within
615 // the fitting range of each successfully fitted peak. Within this range,
616 // the actual intensities are replaced with the values resulting from
617 // evaluating the peak function (e.g. a Gaussian peak function)
618 algFitPeaks->setPropertyValue("FittedPeaksWorkspace", diagnostic_prefix + "_fitted");
619 if (useChiSq) {
620 algFitPeaks->setPropertyValue("OutputParameterFitErrorsWorkspace", diagnostic_prefix + "_fiterrors");
621 }
622
623 // run and get the result
624 algFitPeaks->executeAsChildAlg();
625 g_log.information("finished FitPeaks");
626
627 // get the fit result
628 API::ITableWorkspace_sptr fittedTable = algFitPeaks->getProperty("OutputPeakParametersWorkspace");
629 API::MatrixWorkspace_sptr calculatedWS = algFitPeaks->getProperty("FittedPeaksWorkspace");
630 API::ITableWorkspace_sptr errorTable; // or nullptr as in FitPeaks L1997
631 if (useChiSq) {
632 errorTable = algFitPeaks->getProperty("OutputParameterFitErrorsWorkspace");
633 }
634
635 // check : for Pete
636 if (!fittedTable)
637 throw std::runtime_error("FitPeaks does not have output OutputPeakParametersWorkspace.");
638 if (fittedTable->rowCount() != NUMHIST * m_peaksInDspacing.size())
639 throw std::runtime_error("The number of rows in OutputPeakParametersWorkspace is not correct!");
640
641 // END-OF (FitPeaks)
642 const std::string backgroundType = getPropertyValue("BackgroundType");
643
644 API::Progress prog(this, 0.7, 1.0, NUMHIST);
645
646 // calculate fitting ranges to the left and right of each nominal peak
647 // center, in d-spacing units
648 const auto windowsInDSpacing = dSpacingWindows(m_peaksInDspacing, peakWindow);
649
650 // get spectrum info to check workspace index correpsonds to a valid spectrum
651 const auto &spectrumInfo = m_uncalibratedWS->spectrumInfo();
652
653 // Scan the table containing the fit parameters for every peak, retrieve the
654 // parameters for peaks that were successfully fitting, then use this info
655 // to obtain difc, difa, and tzero for each pixel
656
657 PRAGMA_OMP(parallel for schedule(dynamic, 1))
658 for (int wkspIndex = m_startWorkspaceIndex; wkspIndex <= m_stopWorkspaceIndex; ++wkspIndex) {
660 if ((isEvent && uncalibratedEWS->getSpectrum(wkspIndex).empty()) || !spectrumInfo.hasDetectors(wkspIndex) ||
661 spectrumInfo.isMonitor(wkspIndex) ||
662 maskWS->isMasked(m_uncalibratedWS->getSpectrum(wkspIndex).getDetectorIDs())) {
663 prog.report();
664
665 if (spectrumInfo.hasDetectors(wkspIndex) && !spectrumInfo.isMonitor(wkspIndex)) {
666 if (isEvent && uncalibratedEWS->getSpectrum(wkspIndex).empty()) {
667 maskWS->setMasked(m_uncalibratedWS->getSpectrum(wkspIndex).getDetectorIDs(), true);
668 g_log.debug() << "FULLY masked spectrum, index: " << wkspIndex << "\n";
669 }
670 }
671
672 continue;
673 }
674
675 // object to hold information about the peak positions, detid, and workspace index
677 const auto [dif_c, dif_a, tzero] =
678 getDSpacingToTof(peaks.detid); // doesn't matter which one - all have same difc etc.
679 peaks.setPositions(m_peaksInDspacing, windowsInDSpacing, dif_a, dif_c, tzero);
680
681 // includes peaks that aren't used in the fit
682 // The following data structures will hold information for the peaks
683 // found in the current spectrum
684 const size_t numPeaks = m_peaksInDspacing.size();
685 // TOF of fitted peak centers, default `nan` for failed fitted peaks
686 std::vector<double> tof_vec_full(numPeaks, std::nan(""));
687 std::vector<double> d_vec; // nominal peak centers of fitted peaks
688 std::vector<double> tof_vec; // TOF of fitted peak centers only
689 // width of fitted peak centers, default `nan` for failed fitted peaks
690 std::vector<double> width_vec_full(numPeaks, std::nan(""));
691 // height of fitted peak centers, default `nan` for failed fitted peaks
692 std::vector<double> height_vec_full(numPeaks, std::nan(""));
693 std::vector<double> weights; // weights for diff const fits
694 // row where first peak occurs
695 const size_t rowNumInFitTableOffset = (wkspIndex - m_startWorkspaceIndex) * numPeaks;
696 // We assumed that the current spectrum contains peaks near the nominal
697 // peak centers. Now we check how many peaks we actually found
698 for (size_t peakIndex = 0; peakIndex < numPeaks; ++peakIndex) {
699 size_t rowIndexInFitTable = rowNumInFitTableOffset + peakIndex;
700
701 // check indices in PeaksTable
702 if (fittedTable->getRef<int>("wsindex", rowIndexInFitTable) != wkspIndex)
703 throw std::runtime_error("workspace index mismatch!");
704 if (fittedTable->getRef<int>("peakindex", rowIndexInFitTable) != static_cast<int>(peakIndex))
705 throw std::runtime_error("peak index mismatch but workspace index matched");
706
707 const double chi2 = fittedTable->getRef<double>("chi2", rowIndexInFitTable);
708 double centre = 0.0;
709 double centre_error = 0.0; // only used if useChiSq true
710 double width = 0.0;
711 double height = 0.0;
712 if (!useChiSq) {
713 // get the effective peak parameters from FitPeaks output
714 centre = fittedTable->getRef<double>("centre", rowIndexInFitTable);
715 width = fittedTable->getRef<double>("width", rowIndexInFitTable);
716 height = fittedTable->getRef<double>("height", rowIndexInFitTable);
717 } else {
718 // FitPeaks outputs actual fitting parameters
719 // extract these from the peak function (not efficient)
720 auto peakfunc = std::dynamic_pointer_cast<API::IPeakFunction>(
721 API::FunctionFactory::Instance().createFunction(peakFunction));
722 // set peak functio nparameters from fit
723 for (size_t ipar = 0; ipar < peakfunc->nParams(); ipar++) {
724 peakfunc->setParameter(ipar, fittedTable->getRef<double>(peakfunc->parameterName(ipar), rowIndexInFitTable));
725 }
726 centre = peakfunc->centre();
727 width = peakfunc->fwhm();
728 height = peakfunc->height();
729 centre_error = errorTable->getRef<double>(peakfunc->getCentreParameterName(), rowIndexInFitTable);
730 }
731
732 // check chi-square
733 if (chi2 > maxChiSquared || chi2 < 0.) {
734 g_log.debug("failure to fit: chi2 > maximum");
735 continue; // peak fit deemed as failure
736 }
737
738 // rule out of peak with wrong position. `centre` should be within its
739 // left and right window ranges
740 if (peaks.inTofWindows[2 * peakIndex] >= centre || peaks.inTofWindows[2 * peakIndex + 1] <= centre) {
741 g_log.debug("failure to fit: peak center is out-of-range");
742 continue; // peak fit deemed as failure
743 }
744
745 // check height: make sure 0 is smaller than 0
746 if (height < minPeakHeight + 1.E-15) {
747 g_log.debug("failure to fit: peak height is less than minimum");
748 continue; // peak fit deemed as failure
749 }
750
751 // the peak fit was a success. Collect info
752 g_log.getLogStream(Logger::Priority::PRIO_TRACE) << "successful fit: peak centered at " << centre << "\n";
753
754 d_vec.emplace_back(m_peaksInDspacing[peakIndex]);
755 tof_vec.emplace_back(centre);
756 if (!useChiSq) {
757 weights.emplace_back(height * height);
758 } else {
759 weights.emplace_back(1 / (centre_error * centre_error));
760 }
761 tof_vec_full[peakIndex] = centre;
762 width_vec_full[peakIndex] = width;
763 height_vec_full[peakIndex] = height;
764 }
765
766 if (d_vec.size() < 2) {
767 // If less than two peaks were fitted successfully, indicate failure by
768 // masking all of the detectors contributing to the spectrum.
769 maskWS->setMasked(peaks.detid, true);
770
771 g_log.debug() << "MASKING:\n";
772 for (const auto &det : peaks.detid) {
773 g_log.debug() << " " << det << "\n";
774 }
775 g_log.debug() << "\n";
776
777 continue;
778 } else {
779 // obtain difc, difa, and t0 by fitting the nominal peak center
780 // positions, in d-spacing against the fitted peak center positions, in
781 // TOF units.
782 double difc = 0., t0 = 0., difa = 0.;
783 fitDIFCtZeroDIFA_LM(d_vec, tof_vec, weights, difc, t0, difa);
784 for (auto iter = peaks.detid.begin(); iter != peaks.detid.end(); ++iter) {
785 auto det = *iter;
786 const auto rowIndexOutputPeaks = m_detidToRow[det];
787 // chisq represent the deviations between the nominal peak positions
788 // and the peak positions using the GSAS formula with optimized difc,
789 // difa, and tzero
790 double chisq = 0.;
792 dSpacingUnit.initialize(-1., 0,
796 for (std::size_t i = 0; i < numPeaks; ++i) {
797 if (std::isnan(tof_vec_full[i]))
798 continue;
799 // Find d-spacing using the GSAS formula with optimized difc, difa,
800 // t0 for the TOF of the current peak's center.
801 const double dspacing = dSpacingUnit.singleFromTOF(tof_vec_full[i]);
802 // `temp` is residual between the nominal position in d-spacing for
803 // the current peak, and the fitted position in d-spacing
804 const double temp = m_peaksInDspacing[i] - dspacing;
805 chisq += (temp * temp);
806 m_peakPositionTable->cell<double>(rowIndexOutputPeaks, i + 1) = dspacing;
807 m_peakWidthTable->cell<double>(rowIndexOutputPeaks, i + 1) =
808 WIDTH_TO_FWHM * (width_vec_full[i] / (2 * difa * dspacing + difc));
809 m_peakHeightTable->cell<double>(rowIndexOutputPeaks, i + 1) = height_vec_full[i];
810 }
811 m_peakPositionTable->cell<double>(rowIndexOutputPeaks, m_peaksInDspacing.size() + 1) = chisq;
812 m_peakPositionTable->cell<double>(rowIndexOutputPeaks, m_peaksInDspacing.size() + 2) =
813 chisq / static_cast<double>(numPeaks - 1);
814
815 setCalibrationValues(det, difc, difa, t0);
816 }
817 }
818 prog.report();
819
821 }
823
824 // sort the calibration tables by increasing detector ID
826 setProperty("OutputCalibrationTable", m_calibrationTable);
827
828 // Return the mask workspace only if it was specified as a parameter.
829 if (!isDefault("MaskWorkspace")) {
830 // Align the detector mask flags of the mask workspace with the workspace values:
831 maskWS->combineToDetectorMasks();
832 setProperty("MaskWorkspace", maskWS);
833 }
834
835 // fix-up the diagnostic workspaces
839
840 // a derived table from the position and width
841 auto resolutionWksp = calculateResolutionTable();
842
843 // set the diagnostic workspaces out
844 auto diagnosticGroup = std::make_shared<API::WorkspaceGroup>();
845 // add workspaces calculated by FitPeaks
846 API::AnalysisDataService::Instance().addOrReplace(diagnostic_prefix + "_fitparam", fittedTable);
847 diagnosticGroup->addWorkspace(fittedTable);
848 API::AnalysisDataService::Instance().addOrReplace(diagnostic_prefix + "_fitted", calculatedWS);
849 diagnosticGroup->addWorkspace(calculatedWS);
850 if (useChiSq) {
851 API::AnalysisDataService::Instance().addOrReplace(diagnostic_prefix + "_fiterror", errorTable);
852 diagnosticGroup->addWorkspace(errorTable);
853 }
854
855 // add workspaces calculated by PDCalibration
856 API::AnalysisDataService::Instance().addOrReplace(diagnostic_prefix + "_dspacing", m_peakPositionTable);
857 diagnosticGroup->addWorkspace(m_peakPositionTable);
858 API::AnalysisDataService::Instance().addOrReplace(diagnostic_prefix + "_width", m_peakWidthTable);
859 diagnosticGroup->addWorkspace(m_peakWidthTable);
860 API::AnalysisDataService::Instance().addOrReplace(diagnostic_prefix + "_height", m_peakHeightTable);
861 diagnosticGroup->addWorkspace(m_peakHeightTable);
862 API::AnalysisDataService::Instance().addOrReplace(diagnostic_prefix + "_resolution", resolutionWksp);
863 diagnosticGroup->addWorkspace(resolutionWksp);
864 setProperty("DiagnosticWorkspaces", diagnosticGroup);
865}
866
867namespace { // anonymous namespace
875double gsl_costFunction(const gsl_vector *v, void *peaks) {
876 // this array is [numPeaks, numParams, vector<tof>, vector<dspace>,
877 // vector<weights>]
878 // index as [0, 1, 2, , 2+n , 2+2n]
879 const std::vector<double> *peakVec = reinterpret_cast<std::vector<double> *>(peaks);
880 // number of peaks being fit
881 const auto numPeaks = static_cast<size_t>(peakVec->at(0));
882 // number of parameters
883 const auto numParams = static_cast<size_t>(peakVec->at(1));
884
885 // isn't strictly necessary, but makes reading the code much easier
886 const std::vector<double> tofObs(peakVec->begin() + 2, peakVec->begin() + 2 + numPeaks);
887 const std::vector<double> dspace(peakVec->begin() + (2 + numPeaks), peakVec->begin() + (2 + 2 * numPeaks));
888 const std::vector<double> weights(peakVec->begin() + (2 + 2 * numPeaks), peakVec->begin() + (2 + 3 * numPeaks));
889
890 // create the function to convert tof to dspacing
891 double difc = gsl_vector_get(v, 0);
892 double tzero = 0.;
893 double difa = 0.;
894 if (numParams > 1) {
895 tzero = gsl_vector_get(v, 1);
896 if (numParams > 2)
897 difa = gsl_vector_get(v, 2);
898 }
900 dSpacingUnit.initialize(-1., 0,
904
905 // calculate the sum of the residuals from observed peaks
906 double errsum = 0.0;
907 for (size_t i = 0; i < numPeaks; ++i) {
908 const double tofCalib = dSpacingUnit.singleToTOF(dspace[i]);
909 const double errsum_i = std::pow(tofObs[i] - tofCalib, 2) * weights[i];
910 errsum += errsum_i;
911 }
912
913 return errsum;
914}
915
934double fitDIFCtZeroDIFA(std::vector<double> &peaks, double &difc, double &t0, double &difa) {
935 const auto numParams = static_cast<size_t>(peaks[1]);
936
937 // initial starting point as [DIFC, 0, 0]
938 gsl_vector *fitParams = gsl_vector_alloc(numParams);
939 gsl_vector_set_all(fitParams, 0.0); // set all parameters to zero
940 gsl_vector_set(fitParams, 0, difc);
941 if (numParams > 1) {
942 gsl_vector_set(fitParams, 1, t0);
943 if (numParams > 2) {
944 gsl_vector_set(fitParams, 2, difa);
945 }
946 }
947
948 // Set initial step sizes
949 gsl_vector *stepSizes = gsl_vector_alloc(numParams);
950 gsl_vector_set_all(stepSizes, 0.1);
951 gsl_vector_set(stepSizes, 0, 0.01);
952
953 // Initialize method and iterate
954 gsl_multimin_function minex_func;
955 minex_func.n = numParams;
956 minex_func.f = &gsl_costFunction;
957 minex_func.params = &peaks;
958
959 // Set up GSL minimzer - simplex is overkill
960 const gsl_multimin_fminimizer_type *minimizerType = gsl_multimin_fminimizer_nmsimplex2;
961 gsl_multimin_fminimizer *minimizer = gsl_multimin_fminimizer_alloc(minimizerType, numParams);
962 gsl_multimin_fminimizer_set(minimizer, &minex_func, fitParams, stepSizes);
963
964 // Finally do the fitting
965 size_t iter = 0; // number of iterations
966 const size_t MAX_ITER = 75 * numParams;
967 int status = 0;
968 do {
969 iter++;
970 status = gsl_multimin_fminimizer_iterate(minimizer);
971 if (status)
972 break;
973
974 double size = gsl_multimin_fminimizer_size(minimizer);
975 status = gsl_multimin_test_size(size, 1e-4);
976
977 } while (status == GSL_CONTINUE && iter < MAX_ITER);
978
979 // only update calibration values on successful fit
980 double errsum = 0.; // return 0. if fit didn't work
981 std::string status_msg = gsl_strerror(status);
982 if (status_msg == "success") {
983 difc = gsl_vector_get(minimizer->x, 0);
984 if (numParams > 1) {
985 t0 = gsl_vector_get(minimizer->x, 1);
986 if (numParams > 2) {
987 difa = gsl_vector_get(minimizer->x, 2);
988 }
989 }
990 // return from gsl_costFunction can be accessed as fval
991 errsum = minimizer->fval;
992 }
993
994 // free memory
995 gsl_vector_free(fitParams);
996 gsl_vector_free(stepSizes);
997 gsl_multimin_fminimizer_free(minimizer);
998 return errsum;
999}
1000
1001} // end of anonymous namespace
1002
1017void PDCalibration::fitDIFCtZeroDIFA_LM(const std::vector<double> &d, const std::vector<double> &tof,
1018 const std::vector<double> &weights, double &difc, double &t0, double &difa) {
1019 const size_t numPeaks = d.size();
1020 if (numPeaks <= 1) {
1021 return; // don't do anything
1022 }
1023 // number of fit parameters 1=[DIFC], 2=[DIFC,TZERO], 3=[DIFC,TZERO,DIFA]
1024 // set the maximum number of parameters that will be used
1025 // statistics doesn't support having too few peaks
1026 size_t maxParams = std::min<size_t>(numPeaks - 1, m_numberMaxParams);
1027
1028 // this must have the same layout as the unpacking in gsl_costFunction above
1029 // `peaks` has the following structure for a fit session with three peaks
1030 // and two fit parameters:
1031 // 3, 2, tof1, tof2, tof3, d1, d2, d3, h21, h22, h23
1032 std::vector<double> peaks(numPeaks * 3 + 2, 0.);
1033 peaks[0] = static_cast<double>(d.size());
1034 peaks[1] = 1.; // number of parameters to fit. Initialize to just one
1035 for (size_t i = 0; i < numPeaks; ++i) {
1036 peaks[i + 2] = tof[i];
1037 peaks[i + 2 + numPeaks] = d[i];
1038 peaks[i + 2 + 2 * numPeaks] = weights[i];
1039 }
1040
1041 // calculate a starting DIFC
1042 double difc_start = difc;
1043 if (difc_start == 0.) {
1044 const double d_sum = std::accumulate(d.begin(), d.end(), 0.);
1045 const double tof_sum = std::accumulate(tof.begin(), tof.end(), 0.);
1046 difc_start = tof_sum / d_sum; // number of peaks falls out of division
1047 }
1048
1049 // save the best values so far
1050 double best_errsum = std::numeric_limits<double>::max();
1051 double best_difc = difc_start;
1052 double best_t0 = 0.;
1053 double best_difa = 0.;
1054
1055 // loop over possible number of parameters, doing up to three sequential fits.
1056 // We first start with equation TOF = DIFC * d and obtain
1057 // optimized DIFC which serves as initial guess for next fit with equation
1058 // TOF = DIFC * d + TZERO, obtaining optimized DIFC and TZERO which serve as
1059 // initial guess for final fit TOF = DIFC * d + DIFA * d^2 + TZERO.
1060 for (size_t numParams = 1; numParams <= maxParams; ++numParams) {
1061 peaks[1] = static_cast<double>(numParams);
1062 double difc_local = best_difc;
1063 double t0_local = best_t0;
1064 double difa_local = best_difa;
1065 double errsum = fitDIFCtZeroDIFA(peaks, difc_local, t0_local, difa_local);
1066 if (errsum > 0.) {
1067 // normalize by degrees of freedom
1068 errsum = errsum / static_cast<double>(numPeaks - numParams);
1069 // save the best and forget the rest
1070 // the selected (DIFC, DIFA, TZERO) correspond to those of the fit that
1071 // minimizes errsum. It doesn't have to be the last fit.
1072 if (errsum < best_errsum) {
1073 if (difa_local > m_difaMax || difa_local < m_difaMin)
1074 continue; // unphysical fit
1075 if (t0_local > m_tzeroMax || t0_local < m_tzeroMin)
1076 continue; // unphysical fit
1077 best_errsum = errsum;
1078 best_difc = difc_local;
1079 best_t0 = t0_local;
1080 best_difa = difa_local;
1081 }
1082 }
1083 }
1084
1085 difc = best_difc;
1086 // check that something actually fit and set to the best result
1087 if (best_difc > 0. && best_errsum < std::numeric_limits<double>::max()) {
1088 t0 = best_t0;
1089 difa = best_difa;
1090 }
1091}
1092
1102vector<double> PDCalibration::dSpacingWindows(const std::vector<double> &centres,
1103 const std::vector<double> &windows_in) {
1104
1105 if (!(windows_in.size() == 1 || windows_in.size() / 2 == centres.size()))
1106 throw std::logic_error("the peak-window vector must contain either a single peak-width value, or a pair of values "
1107 "for each peak center specified");
1108
1109 const std::size_t numPeaks = centres.size();
1110
1111 // assumes distance between peaks can be used for window sizes
1112 if (!(numPeaks >= 2))
1113 throw std::logic_error("at least two peak centres must be specified: the distance between these centres will be "
1114 "used to estimate the peak widths");
1115
1116 vector<double> windows_out(2 * numPeaks);
1117 double left;
1118 double right;
1119 for (std::size_t i = 0; i < numPeaks; ++i) {
1120 if (windows_in.size() == 1) {
1121 left = centres[i] - windows_in[0];
1122 right = centres[i] + windows_in[0];
1123 } else {
1124 left = windows_in[2 * i];
1125 right = windows_in[2 * i + 1];
1126 }
1127 // check boundaries don't exceed half dist to adjacent peaks
1128 if (i > 0) {
1129 left = std::max(left, centres[i] - 0.5 * (centres[i] - centres[i - 1]));
1130 }
1131 if (i < numPeaks - 1) {
1132 right = std::min(right, centres[i] + 0.5 * (centres[i + 1] - centres[i]));
1133 }
1134 // set the windows
1135 windows_out[2 * i] = left;
1136 windows_out[2 * i + 1] = right;
1137 }
1138 return windows_out;
1139}
1140
1148std::tuple<double, double, double> PDCalibration::getDSpacingToTof(const std::set<detid_t> &detIds) {
1149 // to start this is the old calibration values
1150 double difc = 0.;
1151 double difa = 0.;
1152 double tzero = 0.;
1153 for (auto detId : detIds) {
1154 auto rowNum = m_detidToRow[detId];
1155 difc += m_calibrationTable->getRef<double>("difc", rowNum);
1156 difa += m_calibrationTable->getRef<double>("difa", rowNum);
1157 tzero += m_calibrationTable->getRef<double>("tzero", rowNum);
1158 }
1159 if (detIds.size() > 1) {
1160 double norm = 1. / static_cast<double>(detIds.size());
1161 difc = norm * difc;
1162 difa = norm * difa;
1163 tzero = norm * tzero;
1164 }
1165
1166 return {difc, difa, tzero};
1167}
1168
1169void PDCalibration::setCalibrationValues(const detid_t detid, const double difc, const double difa,
1170 const double tzero) {
1171 // don't set values that aren't in the table
1172 const auto rowIter = m_detidToRow.find(detid);
1173 if (rowIter == m_detidToRow.end())
1174 return;
1175
1176 // get the row number
1177 auto rowNum = rowIter->second;
1178
1179 // detid is already there
1180 m_calibrationTable->cell<double>(rowNum, 1) = difc;
1181 m_calibrationTable->cell<double>(rowNum, 2) = difa;
1182 m_calibrationTable->cell<double>(rowNum, 3) = tzero;
1183
1184 size_t hasDasIdsOffset = 0; // because it adds a column
1185 if (m_hasDasIds)
1186 hasDasIdsOffset++;
1187
1188 const auto tofMinMax = getTOFminmax(difc, difa, tzero);
1189 m_calibrationTable->cell<double>(rowNum, 4 + hasDasIdsOffset) = tofMinMax[0];
1190 m_calibrationTable->cell<double>(rowNum, 5 + hasDasIdsOffset) = tofMinMax[1];
1191}
1192
1201vector<double> PDCalibration::getTOFminmax(const double difc, const double difa, const double tzero) {
1202 vector<double> tofminmax(2);
1203
1204 Kernel::Units::dSpacing dSpacingUnit;
1205 tofminmax[0] = dSpacingUnit.calcTofMin(difc, difa, tzero, m_tofMin);
1206 tofminmax[1] = dSpacingUnit.calcTofMax(difc, difa, tzero, m_tofMax);
1207
1208 return tofminmax;
1209}
1210MatrixWorkspace_sptr PDCalibration::load(const std::string &filename) {
1211 // TODO this assumes that all files are event-based
1212 const double maxChunkSize = getProperty("MaxChunkSize");
1213 const double filterBadPulses = getProperty("FilterBadPulses");
1214
1215 auto alg = createChildAlgorithm("LoadEventAndCompress");
1216 alg->setLoggingOffset(1);
1217 alg->setProperty("Filename", filename);
1218 alg->setProperty("MaxChunkSize", maxChunkSize);
1219 alg->setProperty("FilterByTofMin", m_tofMin);
1220 alg->setProperty("FilterByTofMax", m_tofMax);
1221 alg->setProperty("FilterBadPulses", filterBadPulses);
1222 alg->setProperty("LoadMonitors", false);
1223 alg->executeAsChildAlg();
1224 API::Workspace_sptr workspace = alg->getProperty("OutputWorkspace");
1225
1226 return std::dynamic_pointer_cast<MatrixWorkspace>(workspace);
1227}
1228
1234
1236 g_log.information("Binning data in time-of-flight");
1237 auto rebin = createChildAlgorithm("Rebin");
1238 rebin->setLoggingOffset(1);
1239 rebin->setProperty("InputWorkspace", wksp);
1240 rebin->setProperty("OutputWorkspace", wksp);
1241 rebin->setProperty("Params", getPropertyValue("TofBinning"));
1242 rebin->setProperty("PreserveEvents", true);
1243 rebin->executeAsChildAlg();
1244 wksp = rebin->getProperty("OutputWorkspace");
1245
1246 return wksp;
1247}
1248
1250 std::set<detid_t> detids;
1251
1252 // return early since everything is being used
1253 if (isDefault("StartWorkspaceIndex") && isDefault("StopWorkspaceIndex"))
1254 return detids;
1255
1256 // get the indices to loop over
1257 std::size_t startIndex = static_cast<std::size_t>(m_startWorkspaceIndex);
1258 std::size_t stopIndex = static_cast<std::size_t>(m_stopWorkspaceIndex);
1259
1260 for (std::size_t i = startIndex; i <= stopIndex; ++i) {
1261 const auto detidsForSpectrum = m_uncalibratedWS->getSpectrum(i).getDetectorIDs();
1262 for (const auto &detid : detidsForSpectrum) {
1263 detids.emplace(detid);
1264 }
1265 }
1266 return detids;
1267}
1268
1270 // create a new workspace
1271 m_calibrationTable = std::make_shared<DataObjects::TableWorkspace>();
1272 // TODO m_calibrationTable->setTitle("");
1273 m_calibrationTable->addColumn("int", "detid");
1274 m_calibrationTable->addColumn("double", "difc");
1275 m_calibrationTable->addColumn("double", "difa");
1276 m_calibrationTable->addColumn("double", "tzero");
1277 if (m_hasDasIds)
1278 m_calibrationTable->addColumn("int", "dasid");
1279 m_calibrationTable->addColumn("double", "tofmin");
1280 m_calibrationTable->addColumn("double", "tofmax");
1281}
1282
1297 API::ITableWorkspace_sptr calibrationTableOld = getProperty("PreviousCalibrationTable");
1298 if (calibrationTableOld == nullptr) {
1299 // load from file
1300 std::string filename = getProperty("PreviousCalibrationFile");
1301 auto alg = createChildAlgorithm("LoadDiffCal");
1302 alg->setLoggingOffset(1);
1303 alg->setProperty("Filename", filename);
1304 alg->setProperty("WorkspaceName", "NOMold"); // TODO
1305 alg->setProperty("MakeGroupingWorkspace", false);
1306 alg->setProperty("MakeMaskWorkspace", false);
1307 alg->setProperty("TofMin", m_tofMin);
1308 alg->setProperty("TofMax", m_tofMax);
1309 alg->executeAsChildAlg();
1310 calibrationTableOld = alg->getProperty("OutputCalWorkspace");
1311 }
1312
1313 m_hasDasIds = hasDasIDs(calibrationTableOld);
1314
1315 // create a new workspace
1316 this->createCalTableHeader();
1317
1318 const auto includedDetids = detIdsForTable();
1319 const bool useAllDetids = includedDetids.empty();
1320
1321 // generate the map of detid -> row
1322 API::ColumnVector<int> detIDs = calibrationTableOld->getVector("detid");
1323 std::size_t rowNum = 0;
1324 for (std::size_t i = 0; i < detIDs.size(); ++i) {
1325 // only add rows for detids that exist in input workspace
1326 if ((useAllDetids) || (includedDetids.count(detIDs[i]) > 0)) {
1327 m_detidToRow[static_cast<detid_t>(detIDs[i])] = rowNum++;
1328 API::TableRow newRow = m_calibrationTable->appendRow();
1329 int detid = calibrationTableOld->getRef<int>("detid", i);
1330 double difc = calibrationTableOld->getRef<double>("difc", i);
1331 double difa = calibrationTableOld->getRef<double>("difa", i);
1332 double tzero = calibrationTableOld->getRef<double>("tzero", i);
1333
1334 newRow << detid << difc << difa << tzero;
1335 if (m_hasDasIds)
1336 newRow << calibrationTableOld->getRef<int>("dasid", i);
1337
1338 // adjust tofmin and tofmax for this pixel
1339 const auto tofMinMax = getTOFminmax(difc, difa, tzero);
1340 newRow << tofMinMax[0] << tofMinMax[1]; // tofmin/tofmax
1341 }
1342 }
1343}
1344
1352 // create new calibraion table for when an old one isn't loaded
1353 // using the signal workspace and CalculateDIFC
1354 auto alg = createChildAlgorithm("CalculateDIFC");
1355 alg->setLoggingOffset(1);
1356 alg->setProperty("InputWorkspace", m_uncalibratedWS);
1357 alg->executeAsChildAlg();
1358 API::MatrixWorkspace_const_sptr difcWS = alg->getProperty("OutputWorkspace");
1359
1360 // create a new workspace
1361 this->createCalTableHeader();
1362
1363 const detid2index_map allDetectors = difcWS->getDetectorIDToWorkspaceIndexMap(false);
1364
1365 const auto includedDetids = detIdsForTable();
1366 const bool useAllDetids = includedDetids.empty();
1367
1368 // copy over the values
1369 size_t i = 0;
1370 for (auto it = allDetectors.begin(); it != allDetectors.end(); ++it) {
1371 const detid_t detID = it->first;
1372 // only add rows for detids that exist in input workspace
1373 if (useAllDetids || (includedDetids.count(detID) > 0)) {
1374 m_detidToRow[detID] = i++;
1375 const size_t wi = it->second;
1376 API::TableRow newRow = m_calibrationTable->appendRow();
1377 newRow << detID;
1378 newRow << difcWS->y(wi)[0];
1379 newRow << 0.; // difa
1380 newRow << 0.; // tzero
1381 newRow << 0.; // tofmin
1382 newRow << DBL_MAX; // tofmax
1383 }
1384 }
1385}
1386
1393 // table for the fitted location of the various peaks, in d-spacing units
1394 m_peakPositionTable = std::make_shared<DataObjects::TableWorkspace>();
1395 m_peakWidthTable = std::make_shared<DataObjects::TableWorkspace>();
1396 m_peakHeightTable = std::make_shared<DataObjects::TableWorkspace>();
1397
1398 m_peakPositionTable->addColumn("int", "detid");
1399 m_peakWidthTable->addColumn("int", "detid");
1400 m_peakHeightTable->addColumn("int", "detid");
1401
1402 for (double dSpacing : m_peaksInDspacing) {
1403 std::stringstream namess;
1404 namess << "@" << std::setprecision(5) << dSpacing;
1405 m_peakPositionTable->addColumn("double", namess.str());
1406 m_peakWidthTable->addColumn("double", namess.str());
1407 m_peakHeightTable->addColumn("double", namess.str());
1408 }
1409 m_peakPositionTable->addColumn("double", "chisq");
1410 m_peakPositionTable->addColumn("double", "normchisq");
1411 // residuals aren't needed for FWHM or height
1412
1413 // convert the map of m_detidToRow to be a vector of detector ids
1414 std::vector<detid_t> detIds(m_detidToRow.size());
1415 for (const auto &it : m_detidToRow) {
1416 detIds[it.second] = it.first;
1417 }
1418
1419 // copy the detector ids from the main table and add lots of NaNs
1420 for (const auto &detId : detIds) {
1421 API::TableRow newPosRow = m_peakPositionTable->appendRow();
1422 API::TableRow newWidthRow = m_peakWidthTable->appendRow();
1423 API::TableRow newHeightRow = m_peakHeightTable->appendRow();
1424
1425 newPosRow << detId;
1426 newWidthRow << detId;
1427 newHeightRow << detId;
1428
1429 for (double dSpacing : m_peaksInDspacing) {
1430 UNUSED_ARG(dSpacing);
1431 newPosRow << std::nan("");
1432 newWidthRow << std::nan("");
1433 newHeightRow << std::nan("");
1434 }
1435 }
1436}
1437
1440 std::make_shared<DataObjects::SpecialWorkspace2D>(m_uncalibratedWS->getInstrument());
1441 resolutionWksp->setTitle("average width/height");
1442
1443 // assume both tables have the same number of rows b/c the algorithm created
1444 // both
1445 // they are also in the same order
1446 // accessing cells is done by (row, col)
1447 const size_t numRows = m_peakPositionTable->rowCount();
1448 const size_t numPeaks = m_peaksInDspacing.size();
1449 std::vector<double> resolution; // vector of non-nan resolutions
1450 for (size_t rowIndex = 0; rowIndex < numRows; ++rowIndex) {
1451 resolution.clear();
1452 // first column is detid
1453 const auto detId = static_cast<detid_t>(m_peakPositionTable->Int(rowIndex, 0));
1454 for (size_t peakIndex = 1; peakIndex < numPeaks + 1; ++peakIndex) {
1455 const double pos = m_peakPositionTable->Double(rowIndex, peakIndex);
1456 if (std::isnormal(pos)) {
1457 resolution.emplace_back(m_peakWidthTable->Double(rowIndex, peakIndex) / pos);
1458 }
1459 }
1460 if (resolution.empty()) {
1461 resolutionWksp->setValue(detId, 0., 0.); // instview doesn't like nan
1462 } else {
1463 // calculate the mean
1464 const double mean =
1465 std::accumulate(resolution.begin(), resolution.end(), 0.) / static_cast<double>(resolution.size());
1466 double stddev = 0.;
1467 std::for_each(resolution.cbegin(), resolution.cend(),
1468 [&stddev, mean](const auto value) { stddev += (value - mean) * (value * mean); });
1469 stddev = std::sqrt(stddev / static_cast<double>(resolution.size() - 1));
1470 resolutionWksp->setValue(detId, mean, stddev);
1471 }
1472 }
1473
1474 return resolutionWksp;
1475}
1476
1479 auto alg = createChildAlgorithm("SortTableWorkspace");
1480 alg->setLoggingOffset(1);
1481 alg->setProperty("InputWorkspace", table);
1482 alg->setProperty("OutputWorkspace", table);
1483 alg->setProperty("Columns", "detid");
1484 alg->executeAsChildAlg();
1485 table = alg->getProperty("OutputWorkspace");
1486
1487 return table;
1488}
1489
1504std::pair<API::MatrixWorkspace_sptr, API::MatrixWorkspace_sptr>
1506 const std::vector<double> &peakWindow) {
1507
1508 // calculate fitting ranges to the left and right of each nominal peak
1509 // center, in d-spacing units
1510 const auto windowsInDSpacing = dSpacingWindows(m_peaksInDspacing, peakWindow);
1511
1512 g_log.information() << "DSPACING WINDOWS\n";
1513 for (std::size_t i = 0; i < m_peaksInDspacing.size(); ++i) {
1514 g_log.information() << "[" << i << "] " << windowsInDSpacing[2 * i] << " < " << m_peaksInDspacing[i] << " < "
1515 << windowsInDSpacing[2 * i + 1] << "\n";
1516 }
1517
1518 // create workspaces for nominal peak centers and fit ranges
1519 size_t numspec = dataws->getNumberHistograms();
1520 size_t numpeaks = m_peaksInDspacing.size();
1521 MatrixWorkspace_sptr peak_pos_ws = create<Workspace2D>(numspec, Points(numpeaks));
1522 MatrixWorkspace_sptr peak_window_ws = create<Workspace2D>(numspec, Points(numpeaks * 2));
1523
1524 const auto NUM_HIST = m_stopWorkspaceIndex - m_startWorkspaceIndex + 1;
1525 API::Progress prog(this, 0., .2, NUM_HIST);
1526
1527 g_log.information() << "TOF WINDOWS\n";
1528 PRAGMA_OMP(parallel for schedule(dynamic, 1) )
1529 for (int64_t iiws = m_startWorkspaceIndex; iiws <= static_cast<int64_t>(m_stopWorkspaceIndex); iiws++) {
1531 std::size_t iws = static_cast<std::size_t>(iiws);
1532 // calculatePositionWindowInTOF
1533 PDCalibration::FittedPeaks peaks(dataws, iws);
1534 // toTof is a function that converts from d-spacing to TOF for a particular
1535 // pixel
1536 auto [difc, difa, tzero] = getDSpacingToTof(peaks.detid);
1537 // setpositions initializes peaks.inTofPos and peaks.inTofWindows
1538 peaks.setPositions(m_peaksInDspacing, windowsInDSpacing, difa, difc, tzero);
1539 peak_pos_ws->setPoints(iws, peaks.inTofPos);
1540 peak_window_ws->setPoints(iws, peaks.inTofWindows);
1541 prog.report();
1542
1543 for (std::size_t i = 0; i < peaks.inTofPos.size(); i++) {
1544 g_log.information() << "[" << iws << "," << i << "] " << peaks.inTofWindows[2 * i] << " < " << peaks.inTofPos[i]
1545 << " < " << peaks.inTofWindows[2 * i + 1] << "\n";
1546 }
1547
1549 }
1551
1552 return std::make_pair(peak_pos_ws, peak_window_ws);
1553}
1554
1555} // namespace Mantid::Algorithms
#define DECLARE_ALGORITHM(classname)
Definition Algorithm.h:538
double value
The value of the point.
Definition FitMW.cpp:51
double height
Definition GetAllEi.cpp:155
IPeaksWorkspace_sptr workspace
double left
double right
#define PARALLEL_START_INTERRUPT_REGION
Begins a block to skip processing is the algorithm has been interupted Note the end of the block if n...
#define PARALLEL_END_INTERRUPT_REGION
Ends a block to skip processing is the algorithm has been interupted Note the start of the block if n...
#define PRAGMA_OMP(expression)
#define PARALLEL_CHECK_INTERRUPT_REGION
Adds a check after a Parallel region to see if it was interupted.
#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.
std::string getPropertyValue(const std::string &name) const override
Get the value of a property as a string.
TypedValue getProperty(const std::string &name) const override
Get the value of a property.
virtual std::shared_ptr< Algorithm > createChildAlgorithm(const std::string &name, const double startProgress=-1., const double endProgress=-1., const bool enableLogging=true, const int &version=-1)
Create a Child Algorithm.
Kernel::Logger & g_log
Definition Algorithm.h:422
bool isDefault(const std::string &name) const
ColumnVector gives access to the column elements without alowing its resizing.
size_t size()
Size of the vector.
A specialized class for dealing with file properties.
@ OptionalLoad
to specify a file to read but the file doesn't have to exist
Base MatrixWorkspace Abstract Class.
Helper class for reporting progress from algorithms.
Definition Progress.h:25
TableRow represents a row in a TableWorkspace.
Definition TableRow.h:39
A property class for workspaces.
void setPositions(const std::vector< double > &peaksInD, const std::vector< double > &peaksInDWindows, const double difa, const double difc, const double tzero)
Store the positions of the peak centers, as well as the left and right fit ranges for each peak,...
FittedPeaks(const API::MatrixWorkspace_const_sptr &wksp, const std::size_t wkspIndex)
Find the bins with non-zero counts and with minimum and maximum TOF.
API::ITableWorkspace_sptr m_peakWidthTable
std::map< detid_t, size_t > m_detidToRow
API::ITableWorkspace_sptr m_peakHeightTable
void setCalibrationValues(const detid_t detid, const double difc, const double difa, const double tzero)
API::MatrixWorkspace_sptr m_uncalibratedWS
API::MatrixWorkspace_sptr rebin(API::MatrixWorkspace_sptr wksp)
int m_stopWorkspaceIndex
stop index (workspace index of the last spectrum included)
API::ITableWorkspace_sptr sortTableWorkspace(API::ITableWorkspace_sptr &table)
sort the calibration table according increasing values in column "detid"
std::vector< double > dSpacingWindows(const std::vector< double > &centres, const std::vector< double > &widthMax)
Fitting ranges to the left and right of peak center (the window cannot exceed half the distance to th...
const std::string category() const override
Algorithm's category for identification.
std::pair< API::MatrixWorkspace_sptr, API::MatrixWorkspace_sptr > createTOFPeakCenterFitWindowWorkspaces(const API::MatrixWorkspace_sptr &dataws, const std::vector< double > &peakWindowMaxInDSpacing)
NEW: convert peak positions in dSpacing to peak centers workspace.
API::ITableWorkspace_sptr m_calibrationTable
void init() override
Initialize the algorithm's properties.
API::MatrixWorkspace_sptr loadAndBin()
load input workspace and rebin with parameters "TofBinning" provided by User
std::vector< double > getTOFminmax(const double difc, const double difa, const double tzero)
Adjustment of TofMin and TofMax values, to ensure positive values of d-spacing when converting from T...
void createInformationWorkspaces()
Table workspaces where the first column is the detector ID and subsequent columns are termed "@x....
~PDCalibration() override
Destructor.
std::map< std::string, std::string > validateInputs() override
Method checking errors on ALL the inputs, before execution.
const std::string name() const override
Algorithms name for identification.
void fitDIFCtZeroDIFA_LM(const std::vector< double > &d, const std::vector< double > &tof, const std::vector< double > &height2, double &difc, double &t0, double &difa)
Fit the nominal peak center positions, in d-spacing against the fitted peak center positions,...
double m_tofMin
first bin boundary when rebinning in TOF (user input)
std::tuple< double, double, double > getDSpacingToTof(const std::set< detid_t > &detIds)
Return a function that converts from d-spacing to TOF units for a particular pixel,...
API::MatrixWorkspace_sptr load(const std::string &filename)
void createCalTableNew()
Initialize the calibration table workspace.
API::MatrixWorkspace_sptr calculateResolutionTable()
int version() const override
Algorithm's version for identification.
const std::string summary() const override
Algorithm's summary for use in the GUI and help.
API::ITableWorkspace_sptr m_peakPositionTable
std::vector< double > m_peaksInDspacing
void exec() override
Execute the algorithm.
void createCalTableFromExisting()
Read a calibration table workspace provided by user, or load from a file provided by User.
double m_tofMax
last bin boundary when rebinning in TOF (user input)
This class is intended to fulfill the design specified in <https://github.com/mantidproject/documents...
Kernel/ArrayBoundedValidator.h.
Support for a property that holds an array of values.
BoundedValidator is a validator that requires the values to be between upper or lower bounds,...
IPropertyManager * setProperty(const std::string &name, const T &value)
Templated method to set the value of a PropertyWithValue.
void setPropertyGroup(const std::string &name, const std::string &group)
Set the group for a given property.
ListValidator is a validator that requires the value of a property to be one of a defined list of pos...
The Logger class is in charge of the publishing messages from the framework through various channels.
Definition Logger.h:51
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
std::ostream & getLogStream(const Priority &priority)
gets the correct log stream for a priority
Definition Logger.cpp:404
void information(const std::string &msg)
Logs at information level.
Definition Logger.cpp:136
Validator to check that a property is not left empty.
void report()
Increments the loop counter by 1, then sends the progress notification on behalf of its algorithm.
Validator to check the format of a vector providing the rebin parameters to an algorithm.
void initialize(const double &_l1, const int &_emode, const UnitParametersMap &params)
Initialize the unit to perform conversion using singleToTof() and singleFromTof()
Definition Unit.cpp:133
void toTOF(std::vector< double > &xdata, std::vector< double > const &ydata, const double &_l1, const int &_emode, std::initializer_list< std::pair< const UnitParams, double > > params)
Convert from the concrete unit to time-of-flight.
Definition Unit.cpp:148
d-Spacing in Angstrom
Definition Unit.h:351
double calcTofMax(const double difc, const double difa, const double tzero, const double tofmax=0.)
Definition Unit.cpp:782
double calcTofMin(const double difc, const double difa, const double tzero, const double tofmin=0.)
Definition Unit.cpp:775
double singleFromTOF(const double tof) const override
DIFA * d^2 + DIFC * d + T0 - TOF = 0.
Definition Unit.cpp:696
double singleToTOF(const double x) const override
Convert a single X value to TOF.
Definition Unit.cpp:670
std::shared_ptr< ITableWorkspace > ITableWorkspace_sptr
shared pointer to Mantid::API::ITableWorkspace
std::shared_ptr< Workspace > Workspace_sptr
shared pointer to Mantid::API::Workspace
Mantid::Kernel::SingletonHolder< AnalysisDataServiceImpl > AnalysisDataService
std::shared_ptr< const ITableWorkspace > ITableWorkspace_const_sptr
shared pointer to Mantid::API::ITableWorkspace (const version)
std::shared_ptr< const MatrixWorkspace > MatrixWorkspace_const_sptr
shared pointer to the matrix workspace base class (const version)
std::shared_ptr< MatrixWorkspace > MatrixWorkspace_sptr
shared pointer to the matrix workspace base class
static double gsl_costFunction(const gsl_vector *v, void *params)
The gsl_costFunction is optimized by GSL simplex.
std::shared_ptr< const MaskWorkspace > MaskWorkspace_const_sptr
shared pointer to a const MaskWorkspace
std::shared_ptr< MaskWorkspace > MaskWorkspace_sptr
shared pointer to the MaskWorkspace class
std::shared_ptr< SpecialWorkspace2D > SpecialWorkspace2D_sptr
shared pointer to the SpecialWorkspace2D class
std::shared_ptr< const Instrument > Instrument_const_sptr
Shared pointer to an const instrument object.
std::unordered_map< UnitParams, double > UnitParametersMap
Definition Unit.h:30
constexpr int EMPTY_INT() noexcept
Returns what we consider an "empty" integer within a property.
Definition EmptyValues.h:24
int32_t detid_t
Typedef for a detector ID.
std::unordered_map< detid_t, size_t > detid2index_map
Map with key = detector ID, value = workspace index.
constexpr double EMPTY_DBL() noexcept
Returns what we consider an "empty" double within a property.
Definition EmptyValues.h:42
Enumeration for a mandatory/optional property.
Describes the direction (within an algorithm) of a Property.
Definition Property.h:50
@ Input
An input workspace.
Definition Property.h:53
@ Output
An output workspace.
Definition Property.h:54