24using namespace Kernel;
26using namespace Algorithms::PeakParameterHelper;
29using namespace DataObjects;
38 "A 2D workspace with X values of d-spacing");
40 auto mustBePositive = std::make_shared<BoundedValidator<double>>();
41 mustBePositive->setLower(0);
43 declareProperty(
"Step", 0.001, mustBePositive,
"Step size used to bin d-spacing data");
44 declareProperty(
"DReference", 2.0, mustBePositive,
"Center of reference peak in d-space");
45 declareProperty(
"XMin", 0.0,
"Minimum of CrossCorrelation data to search for peak, usually negative");
46 declareProperty(
"XMax", 0.0,
"Maximum of CrossCorrelation data to search for peak, usually positive");
49 "Optional: The name of the output CalFile to save the "
50 "generated OffsetsWorkspace.");
52 "An output workspace containing the offsets.");
58 "Mask workspace (optional input / output workspace):"
59 " when specified, if the workspace already exists, any incoming masked detectors will be combined"
60 " with any additional outgoing masked detectors detected by the algorithm");
64 std::make_shared<StringListValidator>(FunctionFactory::Instance().getFunctionNames<IPeakFunction>()),
65 "The function type for fitting the peaks.");
67 "Whether to esimate FWHM of peak function when estimating fit parameters");
68 declareProperty(
"MaxOffset", 1.0,
"Maximum absolute value of offsets; default is 1");
71 std::vector<std::string> modes{
"Relative",
"Absolute",
"Signed"};
73 declareProperty(
"OffsetMode",
"Relative", std::make_shared<StringListValidator>(modes),
74 "Whether to calculate a relative, absolute, or signed offset");
76 "The known peak centre value from the NIST standard "
77 "information, this is only used in Absolute OffsetMode.");
81 std::map<std::string, std::string> result;
85 result[
"InputWorkspace"] =
"The InputWorkspace must be a MatrixWorkspace.";
89 const auto unit = inputWS->getAxis(0)->unit()->caption();
90 if (unit !=
"Bins of Shift" && unit !=
"d-Spacing") {
91 const auto unitErrorMsg =
"GetDetectorOffsets only supports input workspaces with units 'Bins of Shift' or "
92 "'d-Spacing', your unit was : " +
94 result[
"InputWorkspace"] = unitErrorMsg;
100 const bool maskContainsAllDetIDs = maskWS->containsDetIDs(inputWS->getInstrument()->getDetectorIDs(
true));
102 if (!maskContainsAllDetIDs) {
103 result[
"MaskWorkspace"] =
104 "incoming mask workspace must be fully specified, containing an entry for each detid in the input workspace";
105 }
else if (maskWS->getInstrument()->getNumberDetectors(
true) !=
106 inputWS->getInstrument()->getNumberDetectors(
true)) {
107 result[
"MaskWorkspace"] =
"incoming mask workspace must have the same instrument as the input workspace";
108 }
else if (maskWS->getNumberHistograms() != inputWS->getInstrument()->getNumberDetectors(
true)) {
109 result[
"MaskWorkspace"] =
"incoming mask workspace must have one spectrum per detector";
126 throw std::runtime_error(
"Must specify m_Xmin<m_Xmax");
133 if (mode_str ==
"Absolute") {
137 else if (mode_str ==
"Relative") {
141 else if (mode_str ==
"Signed") {
147 int64_t nspec =
static_cast<int64_t
>(
inputW->getNumberHistograms());
150 auto outputW = std::make_shared<OffsetsWorkspace>(
inputW->getInstrument());
151 size_t N_d = outputW->getNumberHistograms();
152 for (
size_t di = 0; di < N_d; ++di)
153 outputW->mutableY(di) = 0.0;
161 g_log.
debug() <<
"[GetDetectorOffsets]: CREATING new MaskWorkspace.\n";
163 maskWS = std::make_shared<MaskWorkspace>(
inputW->getInstrument());
165 g_log.
debug() <<
"[GetDetectorOffsets]: Using EXISTING MaskWorkspace.\n";
168 maskWS->combineFromDetectorMasks(
inputW->detectorInfo());
171 Progress prog(
this, 0.0, 1.0, nspec);
173 for (int64_t wi = 0; wi < nspec; ++wi) {
176 const auto &dets =
inputW->getSpectrum(
static_cast<size_t>(wi)).getDetectorIDs();
179 if (maskWS->isMasked(dets))
183 double offset =
fitSpectra(
static_cast<size_t>(wi));
184 bool spectrumIsMasked =
false;
186 g_log.
debug() <<
"[GetDetectorOffsets]: fit failure: offset: " << std::abs(offset)
187 <<
" is greater than maximum allowed: " <<
m_maxOffset <<
".\n";
188 spectrumIsMasked =
true;
195 for (
const auto &det : dets) {
196 if (spectrumIsMasked) {
197 maskWS->setMasked(det,
true);
201 if (!maskWS->isMasked(det))
202 outputW->setValue(det, offset);
210 maskWS->combineToDetectorMasks();
211 maskWS->combineToDetectorMasks(outputW->mutableDetectorInfo());
215 trace <<
"[GetDetectorOffsets]: Computed offsets:\n" << std::endl;
216 for (
size_t ns = 0; ns <
inputW->getNumberHistograms(); ++ns) {
217 const auto dets =
inputW->getSpectrum(ns).getDetectorIDs();
218 for (
const auto &det : dets)
219 trace <<
" " << outputW->getValue(det) << (maskWS->isMasked(det) ?
"*" :
"") <<
"\n";
230 std::string filename =
getProperty(
"GroupingFileName");
231 if (!filename.empty()) {
234 childAlg->setProperty(
"OffsetsWorkspace", outputW);
235 childAlg->setProperty(
"MaskWorkspace", maskWS);
236 childAlg->setPropertyValue(
"Filename", filename);
237 childAlg->executeAsChildAlg();
248 constexpr double FIT_FAILURE = DBL_MAX;
251 const auto &yValues =
inputW->y(wksp_index);
252 auto it = std::max_element(yValues.cbegin(), yValues.cend());
255 double peakHeight = *it;
256 const double peakLoc =
inputW->x(wksp_index)[
static_cast<size_t>(it - yValues.cbegin())];
260 if (std::isnan(peakHeight))
261 return (FIT_FAILURE);
266 const auto &histogram =
inputW->histogram(wksp_index);
267 const auto &vector_x = histogram.points();
271 if (start_index != stop_index) {
273 auto bkgdFunction = std::dynamic_pointer_cast<IBackgroundFunction>(fun_ptr->getFunction(0));
274 auto peakFunction = std::dynamic_pointer_cast<IPeakFunction>(fun_ptr->getFunction(1));
275 int result =
estimatePeakParameters(histogram, std::pair<size_t, size_t>(start_index, stop_index), peakFunction,
277 if (result != PeakFitResult::GOOD) {
279 <<
" bad result for estimating peak parameters, using default peak height and loc\n";
283 <<
" range size is zero when estimating peak parameters, using default peak height and loc\n";
295 fit_alg->setProperty(
"Function", fun_ptr);
297 fit_alg->setProperty(
"InputWorkspace",
inputW);
298 fit_alg->setProperty<
int>(
"WorkspaceIndex",
299 static_cast<int>(wksp_index));
300 fit_alg->setProperty(
"StartX",
m_Xmin);
301 fit_alg->setProperty(
"EndX",
m_Xmax);
302 fit_alg->setProperty(
"MaxIterations", 100);
304 fit_alg->executeAsChildAlg();
305 std::string fitStatus = fit_alg->getProperty(
"OutputStatus");
307 if (fitStatus !=
"success") {
308 g_log.
debug() <<
"[GetDetectorOffsets]: Fit algorithm failure: " << fitStatus <<
"\n";
309 return (FIT_FAILURE);
315 double offset = function->getParameter(3);
351 peak->setHeight(peakHeight);
352 peak->setCentre(peakLoc);
353 const double sigma(10.0);
354 peak->setFwhm(2.0 * std::sqrt(2.0 * M_LN2) *
sigma);
357 fitFunc->addFunction(background);
358 fitFunc->addFunction(peak);
360 return std::shared_ptr<IFunction>(fitFunc);
#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_CRITICAL(name)
#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.
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.
void progress(double p, const std::string &msg="", double estimatedTime=0.0, int progressPrecision=0)
Sends ProgressNotification.
bool isDefault(const std::string &name) const
A composite function is a function containing other functions.
@ OptionalSave
to specify a file to write to but an empty string is
The FunctionFactory class is in charge of the creation of concrete instances of fitting functions.
std::shared_ptr< IFunction > createFunction(const std::string &type) const
Creates an instance of a function.
Helper class for reporting progress from algorithms.
A property class for workspaces.
API::MatrixWorkspace_sptr inputW
A pointer to the input workspace.
double m_Xmax
The end of the X range for fitting.
API::IFunction_sptr createFunction(const double peakHeight, const double peakLoc)
Create a function string from the given parameters and the algorithm inputs.
double m_dreference
The expected peak position in d-spacing (?)
std::map< std::string, std::string > validateInputs() override
Method checking errors on ALL the inputs, before execution.
double fitSpectra(const size_t wksp_index)
Call Gaussian as a Child Algorithm to fit the peak in a spectrum.
void init() override
Initialisation method.
bool m_estimateFWHM
Flag to estimate fwhm fit parameter.
void exec() override
Executes the algorithm.
double m_Xmin
The start of the X range for fitting.
double m_dideal
The known peak centre value from the NIST standard.
double m_maxOffset
The maximum absolute value of offsets.
Exception for when an item is not found in a collection.
IPropertyManager * setProperty(const std::string &name, const T &value)
Templated method to set the value of a PropertyWithValue.
The Logger class is in charge of the publishing messages from the framework through various channels.
void debug(const std::string &msg)
Logs at debug level.
void notice(const std::string &msg)
Logs at notice level.
void error(const std::string &msg)
Logs at error level.
std::ostream & getLogStream(const Priority &priority)
gets the correct log stream for a priority
int getLevel() const
Returns the Logger's log level.
void report()
Increments the loop counter by 1, then sends the progress notification on behalf of its algorithm.
The concrete, templated class for properties.
std::shared_ptr< IAlgorithm > IAlgorithm_sptr
shared pointer to Mantid::API::IAlgorithm
std::shared_ptr< const MatrixWorkspace > MatrixWorkspace_const_sptr
shared pointer to the matrix workspace base class (const version)
std::shared_ptr< IFunction > IFunction_sptr
shared pointer to the function base class
MANTID_ALGORITHMS_DLL size_t findXIndex(const vector_like &vecx, const double x, const size_t startindex=0)
Get an index of a value in a sorted vector.
MANTID_ALGORITHMS_DLL int estimatePeakParameters(const HistogramData::Histogram &histogram, const std::pair< size_t, size_t > &peak_window, const API::IPeakFunction_sptr &peakfunction, const API::IBackgroundFunction_sptr &bkgdfunction, bool observe_peak_width, const EstimatePeakWidth peakWidthEstimateApproach, const double peakWidthPercentage, const double minPeakHeight)
Estimate peak parameters by 'observation'.
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::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.
@ Input
An input workspace.
@ Output
An output workspace.