Mantid
Loading...
Searching...
No Matches
InterpolatingRebin.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#include "MantidAPI/Axis.h"
13#include "MantidHistogramData/Histogram.h"
16#include "MantidKernel/Spline.h"
18
19#include <gsl/gsl_errno.h>
20
21#include <boost/lexical_cast.hpp>
22#include <span>
23
24namespace Mantid::Algorithms {
25
26// Register the class into the algorithm factory
27DECLARE_ALGORITHM(InterpolatingRebin)
28
29using namespace Kernel;
30using namespace API;
31using namespace HistogramData;
32using namespace DataObjects;
33
38 declareProperty(std::make_unique<WorkspaceProperty<>>("InputWorkspace", "", Direction::Input),
39 "Workspace containing the input data");
40 declareProperty(std::make_unique<WorkspaceProperty<>>("OutputWorkspace", "", Direction::Output),
41 "The name to give the output workspace");
42
43 declareProperty(std::make_unique<ArrayProperty<double>>("Params", std::make_shared<RebinParamsValidator>()),
44 "A comma separated list of first bin boundary, width, last bin boundary. "
45 "Optionally "
46 "this can be followed by a comma and more widths and last boundary "
47 "pairs. "
48 "Optionally this can also be a single number, which is the bin width. "
49 "In this case, the boundary of binning will be determined by minimum and "
50 "maximum TOF "
51 "values among all events, or previous binning boundary, in case of event "
52 "Workspace, or "
53 "non-event Workspace, respectively. Negative width values indicate "
54 "logarithmic binning. ");
55}
56
62std::map<std::string, std::string> InterpolatingRebin::validateInputs() {
63 std::map<std::string, std::string> helpMessages;
64
65 // as part of validation, force the binwidth to be compatible with set bin mode
66 MatrixWorkspace_sptr inputWS = getProperty("InputWorkspace");
67 std::vector<double> rbParams = getProperty("Params");
68 try {
69 Rebin::rebinParamsFromInput(rbParams, *inputWS, g_log);
70 } catch (std::exception &err) {
71 helpMessages["Params"] = err.what();
72 }
73 return helpMessages;
74}
75
81 // Get the input workspace
82 MatrixWorkspace_sptr inputW = getProperty("InputWorkspace");
83
84 // retrieve the properties
85 std::vector<double> rb_params = Rebin::rebinParamsFromInput(getProperty("Params"), *inputW, g_log);
86 HistogramData::BinEdges XValues_new(0);
87 // create new output X axis
88 const int ntcnew = VectorHelper::createAxisFromRebinParams(rb_params, XValues_new.mutableRawData());
89
90 const auto nHists = static_cast<int>(inputW->getNumberHistograms());
91 // make output Workspace the same type as the input but with the new axes
92 MatrixWorkspace_sptr outputW = create<MatrixWorkspace>(*inputW, BinEdges(ntcnew));
93 // Copy over the 'vertical' axis
94 if (inputW->axes() > 1)
95 outputW->replaceAxis(1, std::unique_ptr<Axis>(inputW->getAxis(1)->clone(outputW.get())));
96 outputW->setDistribution(true);
97
98 // this calculation requires a distribution workspace but deal with the
99 // situation when we don't get this
100 bool distCon(false);
101 if (!inputW->isDistribution()) {
102 g_log.debug() << "Converting the input workspace to a distribution\n";
104 distCon = true;
105 }
106
107 try {
108 // evaluate the rebinned data
109 outputYandEValues(inputW, XValues_new, outputW);
110 } catch (std::exception &) {
111
112 if (distCon) {
113 // we need to return the input workspace to the state we found it in
115 }
116 throw;
117 }
118
119 // check if there was a convert to distribution done previously
120 if (distCon) {
121 g_log.debug() << "Converting the input and output workspaces _from_ distributions\n";
123 // the calculation produces a distribution workspace but if they passed a
124 // non-distribution workspace they maybe not expect it, so convert back to
125 // the same form that was given
127 outputW->setDistribution(false);
128 }
129
130 // Now propagate any masking correctly to the output workspace
131 // More efficient to have this in a separate loop because
132 // MatrixWorkspace::maskBins blocks multi-threading
133 for (int i = 0; i < nHists; ++i) {
134 if (inputW->hasMaskedBins(i)) // Does the current spectrum have any masked bins?
135 {
136 this->propagateMasks(inputW, outputW, i);
137 }
138 }
139
140 for (int i = 0; i < outputW->axes(); ++i) {
141 outputW->getAxis(i)->unit() = inputW->getAxis(i)->unit();
142 }
143
144 // Assign it to the output workspace property
145 setProperty("OutputWorkspace", outputW);
146}
154 const HistogramData::BinEdges &XValues_new,
155 const API::MatrixWorkspace_sptr &outputW) {
156 g_log.debug() << "Preparing to calculate y-values using splines and estimate errors\n";
157
158 // prepare to use GSL functions but don't let them terminate Mantid
159 gsl_error_handler_t *old_handler = gsl_set_error_handler(nullptr);
160
161 const auto histnumber = static_cast<int>(inputW->getNumberHistograms());
162 Progress prog(this, 0.0, 1.0, histnumber);
163 for (int hist = 0; hist < histnumber; ++hist) {
164
165 try {
166 // output data arrays are implicitly filled by function
167 outputW->setHistogram(hist, cubicInterpolation(inputW->histogram(hist), XValues_new));
168 } catch (std::exception &ex) {
169 g_log.error() << "Error in rebin function: " << ex.what() << '\n';
170 throw;
171 }
172
173 // Populate the output workspace X values
174 outputW->setBinEdges(hist, XValues_new);
175
176 prog.report();
177 }
178
179 gsl_set_error_handler(old_handler);
180}
181
208Histogram InterpolatingRebin::cubicInterpolation(const Histogram &oldHistogram, const BinEdges &xNew) const {
209 const auto &yOld = oldHistogram.y();
210
211 const size_t size_old = yOld.size();
212 if (size_old == 0)
213 throw std::runtime_error("Empty spectrum found, aborting!");
214
215 const size_t size_new = xNew.size() - 1; // -1 because BinEdges
216
217 // get the bin centres of the input data
218 auto xCensOld = oldHistogram.points();
219 VectorHelper::convertToBinCentre(oldHistogram.x().rawData(), xCensOld.mutableRawData());
220 // the centres of the output data
221 Points xCensNew(size_new);
222 VectorHelper::convertToBinCentre(xNew.rawData(), xCensNew.mutableRawData());
223
224 // find the range of input values whose x-values just suround the output
225 // x-values
226 size_t oldIn1 = std::lower_bound(xCensOld.begin(), xCensOld.end(), xCensNew.front()) - xCensOld.begin();
227 if (oldIn1 == 0) { // the lowest interpolation value might be out of range but
228 // if it is almost on the boundary let it through
229 if (std::abs(xCensOld.front() - xCensNew.front()) < 1e-8 * (xCensOld.back() - xCensOld.front())) {
230 oldIn1 = 1;
231 // make what should be a very small correction
232 xCensNew.mutableRawData().front() = xCensOld.front();
233 }
234 }
235
236 size_t oldIn2 = std::lower_bound(xCensOld.begin(), xCensOld.end(), xCensNew.back()) - xCensOld.begin();
237 if (oldIn2 == size_old) { // the highest point is nearly out of range of the
238 // input data but if it's very near the boundary let
239 // it through
240 if (std::abs(xCensOld.back() - xCensNew.back()) < 1e-8 * (xCensOld.back() - xCensOld.front())) {
241 oldIn2 = size_old - 1;
242 // make what should be a very small correction
243 xCensNew.mutableRawData().back() = xCensOld.back();
244 }
245 }
246
247 // check that the intepolation points fit well enough within the data for
248 // reliable intepolation to be done
249 bool goodRangeLow(false), goodRangeHigh(false), canInterpol(false);
250 if (oldIn1 > 1) { // set the range of the fit, including more input data to
251 // improve accuracy
252 oldIn1 -= 2;
253 goodRangeLow = true;
254 canInterpol = true;
255 } else {
256 if (oldIn1 > 0) {
257 canInterpol = true;
258 oldIn1--;
259 }
260 }
261
262 if (oldIn2 < size_old - 1) {
263 oldIn2++;
264 goodRangeHigh = true;
265 } else {
266 if (oldIn2 >= size_old) {
267 canInterpol = false;
268 }
269 }
270
271 const auto &xOld = oldHistogram.x();
272 const auto &eOld = oldHistogram.e();
273
274 // No Interpolation branch
275 if (!canInterpol) {
276 if (VectorHelper::isConstantValue(yOld.rawData())) {
277 // this copies the single y-value into the output array, errors are still
278 // calculated from the nearest input data points
279 // this is as much as we need to do in this (trival) case
280 return noInterpolation(oldHistogram, xNew);
281 } else { // some points are two close to the edge of the data
282
283 throw std::invalid_argument(std::string("At least one x-value to interpolate to is outside the "
284 "range of the original data.\n") +
285 "original data range: " + boost::lexical_cast<std::string>(xOld.front()) + " to " +
286 boost::lexical_cast<std::string>(xOld.back()) + "\n" +
287 "range to try to interpolate to " + boost::lexical_cast<std::string>(xNew.front()) +
288 " to " + boost::lexical_cast<std::string>(xNew.back()));
289 }
290 }
291
292 // Can Interpolate
293 // Create Histogram
294 Histogram newHistogram{xNew, Frequencies(xNew.size() - 1)};
295
296 auto &yNew = newHistogram.mutableY();
297 auto &eNew = newHistogram.mutableE();
298
299 if ((!goodRangeLow) || (!goodRangeHigh)) {
300 g_log.information() << "One or more points in the interpolation are near "
301 "the boundary of the input data, these points will "
302 "have slightly less accuracy\n";
303 }
304
305 // get the GSL to allocate the memory
306 const size_t nPoints = oldIn2 - oldIn1 + 1;
307 std::span<double const> xSpan(xCensOld.cbegin() + oldIn1, nPoints);
308 std::span<double const> ySpan(yOld.cbegin() + oldIn1, nPoints);
309 std::span<double const> xNewSpan(xCensNew.cbegin(), xCensNew.cend());
311 std::transform(xCensNew.cbegin(), xCensNew.cend(), eNew.begin(),
312 [this, &xCensOld, &eOld](double x) { return this->estimateError(xCensOld, eOld, x); });
313 return newHistogram;
314}
324Histogram InterpolatingRebin::noInterpolation(const Histogram &oldHistogram, const BinEdges &xNew) const {
325 Histogram newHistogram{xNew, Frequencies(xNew.size() - 1)};
326 auto &yNew = newHistogram.mutableY();
327 auto &eNew = newHistogram.mutableE();
328
329 yNew.assign(yNew.size(), oldHistogram.y().front());
330
331 const auto &xPointData = oldHistogram.points();
332 const auto &eOld = oldHistogram.e();
333
334 // -1 because xNew.size is 1 bigger than eNew
335 std::transform(xNew.cbegin(), xNew.cend() - 1, eNew.begin(),
336 [&](double x) { return estimateError(xPointData, eOld, x); });
337
338 return newHistogram;
339}
351double InterpolatingRebin::estimateError(const Points &xsOld, const HistogramE &esOld, const double xNew) const {
352
353 // find the first point in the array that has a higher value of x, we'll base
354 // some of the error estimate on the error on this point
355
356 const size_t indAbove = std::lower_bound(xsOld.begin(), xsOld.end(), xNew) - xsOld.begin();
357
358 // if the point's x-value is out of the range covered by the x-values in the
359 // input data return the error value at the end of the range
360 if (indAbove == 0) {
361 return esOld.front();
362 }
363 // xsOld is 1 longer than xsOld
364 if (indAbove >= esOld.size()) { // cubicInterpolation() checks that that there
365 // are no empty histograms
366 return esOld.back();
367 }
368
369 const double error1 = esOld[indAbove];
370 // ratio of weightings will be inversely proportional to the distance between
371 // the points
372 double weight1 = xsOld[indAbove] - xNew;
373 // check if the points are close enough agnoring any spurious effects that can
374 // occur with exact comparisons of floating point numbers
375 if (weight1 < 1e-100) {
376 // the point is on an input point, all the weight is on this point ignore
377 // the other
378 return error1;
379 }
380 weight1 = 1 / weight1;
381
382 // if p were zero lower_bound must have found xCensNew <= xCensOld.front() but
383 // in that situation we should have exited before now
384 const double error2 = esOld[indAbove - 1];
385 double weight2 = xNew - xsOld[indAbove - 1];
386 if (weight2 < 1e-100) {
387 // the point is on an input point, all the weight is on this point ignore
388 // the other
389 return error2;
390 }
391 weight2 = 1 / weight2;
392
393 return (weight1 * error1 + weight2 * error2) / (weight1 + weight2);
394}
395
396} // namespace Mantid::Algorithms
#define DECLARE_ALGORITHM(classname)
Definition Algorithm.h:538
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.
Helper class for reporting progress from algorithms.
Definition Progress.h:25
A property class for workspaces.
void outputYandEValues(const API::MatrixWorkspace_const_sptr &inputW, const HistogramData::BinEdges &XValues_new, const API::MatrixWorkspace_sptr &outputW)
Calls the interpolation function for each histogram in the workspace.
void exec() override
Executes the rebin algorithm.
HistogramData::Histogram noInterpolation(const HistogramData::Histogram &oldHistogram, const HistogramData::BinEdges &xNew) const
This can be used whenever the original spectrum is filled with only one value.
std::map< std::string, std::string > validateInputs() override
Overwrites the parent (Rebin) validateInputs Avoids all bother about binning modes and powers.
double estimateError(const HistogramData::Points &xsOld, const HistogramData::HistogramE &esOld, const double xNew) const
Estimates the error on each interpolated point by assuming it is similar to the errors in near by inp...
HistogramData::Histogram cubicInterpolation(const HistogramData::Histogram &oldHistogram, const HistogramData::BinEdges &xNew) const
Uses cubic splines to interpolate the mean rate of change of the integral over the inputed data bins ...
void init() override
Overwrites the parent (Rebin) init method These properties will be declares instead.
void propagateMasks(const API::MatrixWorkspace_const_sptr &inputWS, const API::MatrixWorkspace_sptr &outputWS, const int hist, const bool IgnoreBinErrors=false)
Takes the masks in the input workspace and apportions the weights into the new bins that overlap with...
Definition Rebin.cpp:437
static std::vector< double > rebinParamsFromInput(const std::vector< double > &inParams, const API::MatrixWorkspace &inputWS, Kernel::Logger &logger, const std::string &binModeName="Default")
Return the rebin parameters from a user input.
Definition Rebin.cpp:76
Support for a property that holds an array of values.
static std::vector< Y > getSplinedYValues(std::span< X const > newX, std::span< X const > x, std::span< Y const > y)
Definition Spline.h:84
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 information(const std::string &msg)
Logs at information level.
Definition Logger.cpp:136
void report()
Increments the loop counter by 1, then sends the progress notification on behalf of its algorithm.
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
Kernel::Logger g_log("DetermineSpinStateOrder")
void MANTID_KERNEL_DLL convertToBinCentre(const std::vector< double > &bin_edges, std::vector< double > &bin_centres)
Convert an array of bin boundaries to bin center values.
bool MANTID_KERNEL_DLL isConstantValue(const std::vector< double > &arra)
Assess if all the values in the vector are equal or if there are some different values.
int MANTID_KERNEL_DLL createAxisFromRebinParams(const std::vector< double > &params, std::vector< double > &xnew, const bool resize_xnew=true, const bool full_bins_only=false, const double xMinHint=std::nan(""), const double xMaxHint=std::nan(""), const bool useReverseLogarithmic=false, const double power=-1)
Creates a new output X array given a 'standard' set of rebinning parameters.
static void makeDistribution(const MatrixWorkspace_sptr &workspace, const bool forwards=true)
Divides the data in a workspace by the bin width to make it a distribution.
@ Input
An input workspace.
Definition Property.h:53
@ Output
An output workspace.
Definition Property.h:54