19#include "MantidHistogramData/Histogram.h"
24#include <boost/iterator/counting_iterator.hpp>
34Variances subtractMean(std::vector<double> &signal, std::vector<double> &
error) {
35 double mean = std::accumulate(signal.cbegin(), signal.cend(), 0.0);
36 double errorMeanSquared =
38 const auto n = signal.size();
39 mean /=
static_cast<double>(
n);
40 errorMeanSquared /=
static_cast<double>(
n *
n);
41 double variance = 0.0, errorVariance = 0.0;
42 auto itY = signal.begin();
43 auto itE =
error.begin();
44 for (; itY != signal.end(); ++itY, ++itE) {
46 (*itE) = (*itE) * (*itE) + errorMeanSquared;
47 const double t = (*itY) * (*itY);
49 errorVariance += 4.0 * t * (*itE);
51 return {variance, errorVariance};
61using namespace Kernel;
63using namespace DataObjects;
64using namespace HistogramData;
68 auto wsValidator = std::make_shared<CompositeValidator>();
76 "A 2D workspace with X values of d-spacing");
78 "The name of the output workspace");
80 auto mustBePositive = std::make_shared<BoundedValidator<int>>();
81 mustBePositive->setLower(0);
84 "The Workspace Index of the spectra to correlate all other "
88 "The workspace index of the first member of the range of "
89 "spectra to cross-correlate against.");
91 " The workspace index of the last member of the range of "
92 "spectra to cross-correlate against.");
94 declareProperty(
"XMin", 0.0,
"The starting point of the region to be cross correlated.");
95 declareProperty(
"XMax", 0.0,
"The ending point of the region to be cross correlated.");
106 double maxDSpaceShift =
getProperty(
"MaxDSpaceShift");
107 int referenceSpectra =
getProperty(
"ReferenceSpectra");
110 const int wsIndexMin =
getProperty(
"WorkspaceIndexMin");
111 const int wsIndexMax =
getProperty(
"WorkspaceIndexMax");
113 const auto index_ref =
static_cast<size_t>(referenceSpectra);
115 if (wsIndexMin >= wsIndexMax)
116 throw std::runtime_error(
"Must specify WorkspaceIndexMin<WorkspaceIndexMax");
118 int numSpectra = 1 + wsIndexMax - wsIndexMin;
120 std::vector<size_t> indexes(boost::make_counting_iterator(wsIndexMin), boost::make_counting_iterator(wsIndexMax + 1));
122 if (numSpectra == 0) {
123 std::ostringstream message;
124 message <<
"No spectra in range between" << wsIndexMin <<
" and " << wsIndexMax;
125 throw std::runtime_error(message.str());
128 g_log.
information() <<
"There are " << numSpectra <<
" spectra in the range\n";
135 auto &referenceSpectraE = inputWS->e(index_ref);
136 auto &referenceSpectraX = inputWS->x(index_ref);
137 auto &referenceSpectraY = inputWS->y(index_ref);
139 using std::placeholders::_1;
141 std::find_if(referenceSpectraX.cbegin(), referenceSpectraX.cend(), std::bind(std::greater<double>(), _1, xmin));
142 if (rangeStart == referenceSpectraX.cend())
143 throw std::runtime_error(
"No data above XMin");
144 auto rangeEnd = std::find_if(rangeStart, referenceSpectraX.cend(), std::bind(std::greater<double>(), _1, xmax));
145 if (rangeStart == rangeEnd)
146 throw std::runtime_error(
"Range is not valid");
148 MantidVec::difference_type rangeStartCorrection = std::distance(referenceSpectraX.cbegin(), rangeStart);
149 MantidVec::difference_type rangeEndCorrection = std::distance(referenceSpectraX.cbegin(), rangeEnd);
151 const std::vector<double> referenceXVector(rangeStart, rangeEnd);
152 std::vector<double> referenceYVector(referenceSpectraY.cbegin() + rangeStartCorrection,
153 referenceSpectraY.cbegin() + (rangeEndCorrection - 1));
154 std::vector<double> referenceEVector(referenceSpectraE.cbegin() + rangeStartCorrection,
155 referenceSpectraE.cbegin() + (rangeEndCorrection - 1));
157 g_log.
information() <<
"min max " << referenceXVector.front() <<
" " << referenceXVector.back() <<
'\n';
161 auto numReferenceY =
static_cast<int>(referenceYVector.size());
164 int shiftCorrection = 0;
166 if (xmax - xmin < maxDSpaceShift)
168 <<
") is larger than specified range of xmin(" << xmin <<
") to xmax(" << xmax
169 <<
"), please make it smaller or removed it entirely!"
173 const auto maxBins = std::max(0.0 + maxDSpaceShift * 2, 0.1) / inputWS->getDimension(0)->getBinWidth();
175 shiftCorrection = (int)std::max(0.0, abs((-numReferenceY + 2) - (numReferenceY - 2)) - maxBins) / 2;
178 const int numPoints = 2 * (numReferenceY - shiftCorrection) - 3;
180 throw std::runtime_error(
"Range is not valid");
184 const auto referenceVariance = subtractMean(referenceYVector, referenceEVector);
186 const double referenceNorm = 1.0 / sqrt(referenceVariance.y);
187 double referenceNormE = 0.5 * pow(referenceNorm, 3) * sqrt(referenceVariance.e);
190 bool isDistribution = inputWS->isDistribution();
192 auto &outX = out->mutableX(0);
193 for (
int i = 0; i < static_cast<int>(outX.size()); ++i) {
194 outX[i] =
static_cast<double>(i - (numReferenceY - shiftCorrection) + 2);
198 m_progress = std::make_unique<Progress>(
this, 0.0, 1.0, numSpectra);
200 for (
int currentSpecIndex = 0; currentSpecIndex < numSpectra; ++currentSpecIndex)
203 size_t wsIndex = indexes[currentSpecIndex];
205 out->getSpectrum(currentSpecIndex).copyInfoFrom(inputWS->getSpectrum(wsIndex));
206 out->setSharedX(currentSpecIndex, out->sharedX(0));
208 const auto &inputXVector = inputWS->x(wsIndex);
209 const auto &inputYVector = inputWS->y(wsIndex);
210 const auto &inputEVector = inputWS->e(wsIndex);
213 std::vector<double> tempY(numReferenceY);
214 std::vector<double> tempE(numReferenceY);
216 VectorHelper::rebin(inputXVector.rawData(), inputYVector.rawData(), inputEVector.rawData(), referenceXVector, tempY,
217 tempE, isDistribution);
218 const auto tempVar = subtractMean(tempY, tempE);
221 const double tempNorm = 1.0 / sqrt(tempVar.y);
222 const double tempNormE = 0.5 * pow(tempNorm, 3) * sqrt(tempVar.e);
223 const double normalisation = referenceNorm * tempNorm;
224 const double normalisationE2 = pow((referenceNorm * tempNormE), 2) + pow((tempNorm * referenceNormE), 2);
226 auto &outY = out->mutableY(currentSpecIndex);
227 auto &outE = out->mutableE(currentSpecIndex);
229 for (
int dataIndex = -numReferenceY + 2 + shiftCorrection; dataIndex <= numReferenceY - 2 - shiftCorrection;
231 const int dataIndexP = abs(dataIndex);
232 double val = 0, err2 = 0,
x,
y, xE, yE;
233 for (
int j = numReferenceY - 1 - dataIndexP; j >= 0; --j) {
234 if (dataIndex >= 0) {
235 x = referenceYVector[j];
236 y = tempY[j + dataIndexP];
237 xE = referenceEVector[j];
238 yE = tempE[j + dataIndexP];
241 y = referenceYVector[j + dataIndexP];
243 yE = referenceEVector[j + dataIndexP];
246 err2 +=
x *
x * yE +
y *
y * xE;
248 outY[dataIndex + numReferenceY - shiftCorrection - 2] = (val * normalisation);
249 outE[dataIndex + numReferenceY - shiftCorrection - 2] =
250 sqrt(val * val * normalisationE2 + normalisation * normalisation * err2);
259 Unit_sptr unit = out->getAxis(0)->unit();
260 std::shared_ptr<Units::Label> label = std::dynamic_pointer_cast<Units::Label>(unit);
261 label->setLabel(
"Bins of Shift",
"\\mathbb{Z}");
#define DECLARE_ALGORITHM(classname)
#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 PARALLEL_FOR_IF(condition)
Empty definitions - to enable set your complier to enable openMP.
#define PARALLEL_CHECK_INTERRUPT_REGION
Adds a check after a Parallel region to see if it was interupted.
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.
A validator which checks that a workspace contains histogram data (the default) or point data as requ...
A validator which checks that a workspace contains raw counts in its bins.
A property class for workspaces.
A validator which checks that the unit of the workspace referred to by a WorkspaceProperty is the exp...
void init() override
Initialisation code.
std::unique_ptr< API::Progress > m_progress
Progress reporting.
void exec() override
Execution code.
IPropertyManager * setProperty(const std::string &name, const T &value)
Templated method to set the value of a PropertyWithValue.
void warning(const std::string &msg)
Logs at warning level.
void information(const std::string &msg)
Logs at information level.
static T & Instance()
Return a reference to the Singleton instance, creating it if it does not already exist Creation is do...
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
void MANTID_KERNEL_DLL rebin(const std::vector< double > &xold, const std::vector< double > &yold, const std::vector< double > &eold, const std::vector< double > &xnew, std::vector< double > &ynew, std::vector< double > &enew, bool distribution, bool addition=false)
Rebins data according to a new output X array.
std::shared_ptr< Unit > Unit_sptr
Shared pointer to the Unit base class.
std::enable_if< std::is_pointer< Arg >::value, bool >::type threadSafe(Arg workspace)
Thread-safety check Checks the workspace to ensure it is suitable for multithreaded access.
constexpr double EMPTY_DBL() noexcept
Returns what we consider an "empty" double within a property.
std::string to_string(const wide_integer< Bits, Signed > &n)
@ Input
An input workspace.
@ Output
An output workspace.
Functor to accumulate a sum of squares.