Mantid
Loading...
Searching...
No Matches
StripPeaks.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
11#include "MantidAPI/TableRow.h"
17
18namespace Mantid::Algorithms {
19
20// Register the algorithm into the AlgorithmFactory
21DECLARE_ALGORITHM(StripPeaks)
22
23using namespace Kernel;
24using namespace API;
25
27 declareProperty(std::make_unique<WorkspaceProperty<>>("InputWorkspace", "", Direction::Input),
28 "The name of the input workspace.");
29 declareProperty(std::make_unique<WorkspaceProperty<>>("OutputWorkspace", "", Direction::Output),
30 "The name to use for the output workspace.");
31
32 auto min = std::make_shared<BoundedValidator<int>>();
33 min->setLower(1);
34 // The estimated width of a peak in terms of number of channels
35 declareProperty("FWHM", 7, min,
36 "The number of points covered, on average, "
37 "by the fwhm of a peak (default 7).\nPassed "
38 "through to [[FindPeaks]].");
39 // The tolerance allowed in meeting the conditions
40 declareProperty("Tolerance", 4, min,
41 "A measure of the strictness desired in "
42 "meeting the condition on peak "
43 "candidates,\n"
44 "Mariscotti recommends 2 (default 4)");
45
46 declareProperty(std::make_unique<ArrayProperty<double>>("PeakPositions"),
47 "Optional: enter a comma-separated list of the expected "
48 "X-position of the centre of the peaks. Only peaks near "
49 "these positions will be fitted.");
50
51 declareProperty("PeakPositionTolerance", 0.01,
52 "Tolerance on the found peaks' positions against the input "
53 "peak positions. Non-positive value indicates that this "
54 "option is turned off.");
55
56 std::vector<std::string> bkgdtypes{"Linear", "Quadratic"};
57 declareProperty("BackgroundType", "Linear", std::make_shared<StringListValidator>(bkgdtypes),
58 "Type of Background. Present choices include 'Linear' and 'Quadratic'");
59
60 declareProperty("HighBackground", true, "Flag whether the input data has high background compared to peak heights.");
61
62 auto mustBePositive = std::make_shared<BoundedValidator<int>>();
63 mustBePositive->setLower(0);
64 declareProperty("WorkspaceIndex", EMPTY_INT(), mustBePositive,
65 "If set, will remove peaks only in the given spectrum of the "
66 "workspace. Otherwise, all spectra will be searched.");
67
68 auto mustBePositiveDbl = std::make_shared<BoundedValidator<double>>();
69 mustBePositiveDbl->setLower(0.);
70 declareProperty("MaximumChisq", 100., mustBePositiveDbl,
71 "The maximum chisq value for fits to remove the peak. Default 100.");
72}
73
75 // Retrieve the input workspace
76 MatrixWorkspace_sptr inputWS = getProperty("InputWorkspace");
77 m_maxChiSq = getProperty("MaximumChisq");
78
79 // Call FindPeaks as a Child Algorithm
80 ITableWorkspace_sptr peakslist = this->findPeaks(inputWS);
81
82 MatrixWorkspace_sptr outputWS;
83 if (peakslist->rowCount() > 0) {
84 // Remove the peaks that were successfully fitted
85 outputWS = this->removePeaks(inputWS, peakslist);
86 } else {
87 // No peaks found!!! - just assign the input workspace pointer to the output
88 // one
89 g_log.warning("No peaks to remove!");
90 outputWS = inputWS;
91 }
92
93 // Set the output workspace
94 setProperty("OutputWorkspace", outputWS);
95}
96
102 g_log.debug("Calling FindPeaks as a Child Algorithm");
103
104 // Read from properties
105 int fwhm = getProperty("FWHM");
106 int tolerance = getProperty("Tolerance");
107 int wsindex = getProperty("WorkspaceIndex");
108 std::string backgroundtype = getProperty("BackgroundType");
109 bool highbackground = getProperty("HighBackground");
110 std::vector<double> peakpositions = getProperty("PeakPositions");
111 double peakpostol = getProperty("PeakPositionTolerance");
112 if (peakpostol < 0.)
113 peakpostol = EMPTY_DBL();
114
115 // Set up and execute algorithm
116 bool showlog = true;
117 auto findpeaks = createChildAlgorithm("FindPeaks", 0.0, 0.2, showlog);
118 findpeaks->setProperty("InputWorkspace", WS);
119 findpeaks->setProperty<int>("FWHM", fwhm);
120 findpeaks->setProperty<int>("Tolerance", tolerance);
121 findpeaks->setProperty<int>("WorkspaceIndex", wsindex);
122 // Get the specified peak positions, which is optional
123 findpeaks->setProperty<std::vector<double>>("PeakPositions", peakpositions);
124 findpeaks->setProperty<std::string>("BackgroundType", backgroundtype);
125 findpeaks->setProperty<bool>("HighBackground", highbackground);
126 findpeaks->setProperty<double>("PeakPositionTolerance", peakpostol);
127 findpeaks->setProperty<bool>("RawPeakParameters", true);
128
129 findpeaks->executeAsChildAlg();
130
131 ITableWorkspace_sptr peaklistws = findpeaks->getProperty("PeaksList");
132 if (!peaklistws)
133 throw std::runtime_error("Algorithm FindPeaks() returned a NULL pointer "
134 "for table workspace 'PeaksList'.");
135
136 std::stringstream infoss;
137 infoss << "Call FindPeaks() to find " << peakpositions.size() << " peaks with parameters: \n";
138 infoss << "\t FWHM = " << fwhm << ";\n";
139 infoss << "\t Tolerance = " << tolerance << ";\n";
140 infoss << "\t HighBackground = " << highbackground << ";\n";
141 infoss << "\t BackgroundType = " << backgroundtype << ";\n";
142 infoss << "\t Peak positions = ";
143 for (size_t i = 0; i < peakpositions.size(); ++i) {
144 infoss << peakpositions[i];
145 if (i < peakpositions.size() - 1)
146 infoss << ", ";
147 }
148 infoss << ")\n"
149 << "Returned number of fitted peak = " << peaklistws->rowCount() << ".";
150 g_log.information(infoss.str());
151
152 return peaklistws;
153}
154
163 const API::ITableWorkspace_sptr &peakslist) {
164 g_log.debug("Subtracting peaks");
165 // Create an output workspace - same size as input one
166 MatrixWorkspace_sptr outputWS = WorkspaceFactory::Instance().create(input);
167 // Copy the data over from the input to the output workspace
168 const size_t hists = input->getNumberHistograms();
169 // progress 0.2 to 0.3 in this loop
170 double prg = 0.2;
171 for (size_t k = 0; k < hists; ++k) {
172 outputWS->setHistogram(k, input->histogram(k));
173 prg += (0.1 / static_cast<double>(hists));
174 progress(prg);
175 }
176
177 // progress from 0.3 to 1.0 here
178 prg = 0.3;
179 // Loop over the list of peaks
180 for (size_t i = 0; i < peakslist->rowCount(); ++i) {
181 // Get references to the data - X can be const
182 // but Y needs to be mutable
183 auto &X = outputWS->x(peakslist->getRef<int>("spectrum", i));
184 auto &Y = outputWS->mutableY(peakslist->getRef<int>("spectrum", i));
185 // Get back the gaussian parameters
186 const double height = peakslist->getRef<double>("Height", i);
187 const double centre = peakslist->getRef<double>("PeakCentre", i);
188 const double sigma = peakslist->getRef<double>("Sigma", i);
189 const double chisq = peakslist->getRef<double>("chi2", i);
190 // These are some heuristic rules to discard bad fits.
191 // Hope to be able to remove them when we have better fitting routine
192 if (height < 0) {
193 g_log.error() << "Find Peak with Negative Height\n";
194 continue; // Height must be positive
195 }
196 if (chisq > m_maxChiSq) {
197 g_log.information() << "Error for fit peak at " << centre << " is " << chisq << ", which exceeds allowed value "
198 << m_maxChiSq << "\n";
199 if (chisq != DBL_MAX)
200 g_log.error() << "StripPeaks(): Peak Index = " << i << " @ X = " << centre
201 << " Error: Peak fit with too high of chisq " << chisq << " > " << m_maxChiSq << "\n";
202 continue;
203 } else if (chisq < 0.) {
204 g_log.warning() << "StripPeaks(): Peak Index = " << i << " @ X = " << centre
205 << ". Error: Peak fit with too wide peak width" << sigma << " denoted by chi^2 = " << chisq
206 << " <= 0. \n";
207 }
208 {
209 // find the bin width at the center of the peak - average of the bins on
210 // either side
211 const auto center_iter = lower_bound(X.begin(), X.end(), centre);
212 const double bin_width = .5 * (*(center_iter + 1) - (*(center_iter - 1)));
213
214 // seek wikipedia if you don't believe the conversion factor
215 const double fwhm = sigma * 2. * std::sqrt(2. * std::log(2.));
216
217 if ((fwhm / bin_width) < 1.5) {
218 g_log.warning() << "StripPeaks(): Peak Index = " << i << " @ X = " << centre
219 << " Error: Peak fit with too narrow of peak"
220 << " fwhm = " << fwhm << " bin size = " << bin_width
221 << " num bins in peak = " << 2. * (fwhm / bin_width) << " <3\n";
222 continue;
223 }
224 }
225
226 g_log.information() << "Subtracting peak " << i << " from spectrum " << peakslist->getRef<int>("spectrum", i)
227 << " at x = " << centre << " h = " << height << " s = " << sigma << " chi2 = " << chisq << "\n";
228
229 { // log the background function
230 double a0 = 0.;
231 double a1 = 0.;
232 double a2 = 0.;
233 const std::vector<std::string> columnNames = peakslist->getColumnNames();
234 if (std::find(columnNames.begin(), columnNames.end(), "A0") != columnNames.end())
235 a0 = peakslist->getRef<double>("A0", i);
236 if (std::find(columnNames.begin(), columnNames.end(), "A1") != columnNames.end())
237 a1 = peakslist->getRef<double>("A1", i);
238 if (std::find(columnNames.begin(), columnNames.end(), "A2") != columnNames.end())
239 a2 = peakslist->getRef<double>("A2", i);
240 g_log.information() << " background = " << a0 << " + " << a1 << " x + " << a2 << " x^2\n";
241 }
242
243 // Get central bin values using points
244 auto binCenters = outputWS->points(peakslist->getRef<int>("spectrum", i));
245 const auto spectrumLength = static_cast<int>(Y.size());
246 for (int j = 0; j < spectrumLength; ++j) {
247 double x = binCenters[j];
248 // Skip if not anywhere near this peak
249 if (x < centre - 5.0 * sigma)
250 continue;
251 if (x > centre + 5.0 * sigma)
252 break;
253 // Calculate the value of the Gaussian function at this point
254 const double funcVal = height * exp(-0.5 * pow((x - centre) / sigma, 2));
255 // Subtract the calculated value from the data
256 Y[j] -= funcVal;
257 }
258 prg += (0.7 / static_cast<double>(peakslist->rowCount()));
259 progress(prg);
260 }
261
262 return outputWS;
263}
264
265} // namespace Mantid::Algorithms
#define DECLARE_ALGORITHM(classname)
Definition Algorithm.h:538
double sigma
Definition GetAllEi.cpp:156
double height
Definition GetAllEi.cpp:155
double tolerance
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
void progress(double p, const std::string &msg="", double estimatedTime=0.0, int progressPrecision=0)
Sends ProgressNotification.
A property class for workspaces.
void exec() override
Execution code.
API::MatrixWorkspace_sptr removePeaks(const API::MatrixWorkspace_const_sptr &input, const API::ITableWorkspace_sptr &peakslist)
If a peak was successfully fitted, it is subtracted from the data.
void init() override
Initialisation code.
API::ITableWorkspace_sptr findPeaks(const API::MatrixWorkspace_sptr &WS)
Calls FindPeaks as a Child Algorithm.
Support for a property that holds an array of values.
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 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< ITableWorkspace > ITableWorkspace_sptr
shared pointer to Mantid::API::ITableWorkspace
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
constexpr int EMPTY_INT() noexcept
Returns what we consider an "empty" integer within a property.
Definition EmptyValues.h:24
constexpr double EMPTY_DBL() noexcept
Returns what we consider an "empty" double within a property.
Definition EmptyValues.h:42
@ Input
An input workspace.
Definition Property.h:53
@ Output
An output workspace.
Definition Property.h:54