Mantid
Loading...
Searching...
No Matches
IntegratePeaksMD2.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 +
9#include "MantidAPI/Column.h"
16#include "MantidAPI/Progress.h"
17#include "MantidAPI/Run.h"
18#include "MantidAPI/TableRow.h"
19#include "MantidAPI/TextAxis.h"
31#include "MantidHistogramData/LinearGenerator.h"
37#include "MantidKernel/Utils.h"
41
42#include "boost/math/distributions.hpp"
43
44#include <cmath>
45#include <fstream>
46#include <gsl/gsl_integration.h>
47
48namespace Mantid::MDAlgorithms {
49
50// Register the algorithm into the AlgorithmFactory
51DECLARE_ALGORITHM(IntegratePeaksMD2)
52
53using namespace Mantid::HistogramData;
54using namespace Mantid::Kernel;
55using namespace Mantid::API;
56using namespace Mantid::DataObjects;
57using namespace Mantid::DataObjects;
58using namespace Mantid::Geometry;
59
63 declareProperty(std::make_unique<WorkspaceProperty<IMDEventWorkspace>>("InputWorkspace", "", Direction::Input),
64 "An input MDEventWorkspace.");
65
66 auto radiiValidator = std::make_shared<ArrayBoundedValidator<double>>();
67 radiiValidator->setLower(0.0);
68 radiiValidator->setLowerExclusive(true);
69 declareProperty(std::make_unique<ArrayProperty<double>>("PeakRadius", std::vector<double>({1.0}), radiiValidator,
71 "Fixed radius around each peak position in which to integrate, or the "
72 "semi-axis lengths (a,b,c) describing an ellipsoid shape used for "
73 "integration (in the same units as the workspace).");
74
75 radiiValidator->setLowerExclusive(false);
76 declareProperty(std::make_unique<ArrayProperty<double>>("BackgroundInnerRadius", std::vector<double>({0.0}),
77 radiiValidator, Direction::Input),
78 "Inner radius, or three values for semi-axis lengths (a,b,c) of the "
79 "ellipsoid shape, used to evaluate the background of the peak.\n"
80 "If smaller than PeakRadius, then we assume BackgroundInnerRadius = "
81 "PeakRadius.");
82
83 declareProperty(std::make_unique<ArrayProperty<double>>("BackgroundOuterRadius", std::vector<double>({0.0}),
84 radiiValidator, Direction::Input),
85 "Outer radius, or three values for semi-axis lengths (a,b,c) of the "
86 "ellipsoid shape, to use to evaluate the background of the peak.\n"
87 "The signal density around the peak (BackgroundInnerRadius < r < "
88 "BackgroundOuterRadius) is used to estimate the background under the "
89 "peak.\n"
90 "If smaller than PeakRadius, no background measurement is done.");
91
92 declareProperty(std::make_unique<WorkspaceProperty<IPeaksWorkspace>>("PeaksWorkspace", "", Direction::Input),
93 "A PeaksWorkspace containing the peaks to integrate.");
94
95 declareProperty(std::make_unique<WorkspaceProperty<IPeaksWorkspace>>("OutputWorkspace", "", Direction::Output),
96 "The output PeaksWorkspace will be a copy of the input PeaksWorkspace "
97 "with the peaks' integrated intensities.");
98
99 declareProperty("ReplaceIntensity", true,
100 "Always replace intensity in PeaksWorkspacem (default).\n"
101 "If false, then do not replace intensity if calculated value "
102 "is 0 (used for SNSSingleCrystalReduction)");
103
104 declareProperty("IntegrateIfOnEdge", true,
105 "Only warning if all of peak outer radius is not on detector (default).\n"
106 "If false, do not integrate if the outer radius is not on a detector.");
107
108 declareProperty("AdaptiveQBackground", false,
109 "Default is false. If true, "
110 "BackgroundOuterRadius + AdaptiveQMultiplier * **|Q|** and "
111 "BackgroundInnerRadius + AdaptiveQMultiplier * **|Q|**");
112
113 declareProperty("Ellipsoid", false, "Default is sphere.");
114
115 declareProperty("FixQAxis", false, "Fix one axis of ellipsoid to be along direction of Q.");
116
117 declareProperty("Cylinder", false, "Default is sphere. Use next five parameters for cylinder.");
118
119 declareProperty(std::make_unique<PropertyWithValue<double>>("CylinderLength", 0.0, Direction::Input),
120 "Length of cylinder in which to integrate (in the same units as the "
121 "workspace).");
122
123 declareProperty(std::make_unique<PropertyWithValue<double>>("PercentBackground", 0.0, Direction::Input),
124 "Percent of CylinderLength that is background (20 is 20%)");
125
126 std::vector<std::string> peakNames = FunctionFactory::Instance().getFunctionNames<IPeakFunction>();
127 peakNames.emplace_back("NoFit");
128 declareProperty("ProfileFunction", "Gaussian", std::make_shared<StringListValidator>(peakNames),
129 "Fitting function for profile that is used only with "
130 "Cylinder integration.");
131
132 std::vector<std::string> integrationOptions(2);
133 integrationOptions[0] = "Sum";
134 integrationOptions[1] = "GaussianQuadrature";
135 auto integrationvalidator = std::make_shared<StringListValidator>(integrationOptions);
136 declareProperty("IntegrationOption", "GaussianQuadrature", integrationvalidator,
137 "Integration method for calculating intensity "
138 "used only with Cylinder integration.");
139
140 declareProperty(std::make_unique<FileProperty>("ProfilesFile", "", FileProperty::OptionalSave,
141 std::vector<std::string>(1, "profiles")),
142 "Save (Optionally) as Isaw peaks file with profiles included");
143
144 declareProperty("AdaptiveQMultiplier", 0.0,
145 "PeakRadius + AdaptiveQMultiplier * **|Q|** "
146 "so each peak has a "
147 "different integration radius. Q includes the 2*pi factor.");
148
149 declareProperty("CorrectIfOnEdge", false,
150 "Only warning if all of peak outer radius is not on detector (default).\n"
151 "If false, correct for volume off edge for both background and "
152 "intensity (the peak is assumed uniform Gaussian so this only applies "
153 "to spherical integration).");
154
155 declareProperty("UseOnePercentBackgroundCorrection", true,
156 "If this options is enabled, then the top 1% of the "
157 "background will be removed"
158 "before the background subtraction.");
159
160 // continued ellipsoid args
161 declareProperty("FixMajorAxisLength", true,
162 "This option is ignored if all peak radii are specified. "
163 "Otherwise, if True the ellipsoid radidi (proportional to "
164 "the sqrt of the eigenvalues of the covariance matrix) are "
165 "scaled such that the major axis radius is equal to the "
166 "PeakRadius. If False then the ellipsoid radii are set to "
167 "3 times the sqrt of the eigenvalues of the covariance "
168 "matrix");
169
170 declareProperty("UseCentroid", false,
171 "Perform integration on estimated centroid not peak position "
172 "(ignored if all three peak radii are specified).");
173
174 auto maxIterValidator = std::make_shared<BoundedValidator<int>>();
175 maxIterValidator->setLower(1);
176 declareProperty("MaxIterations", 1, maxIterValidator,
177 "Number of iterations in covariance estimation (ignored if all "
178 "peak radii are specified). 2-3 should be sufficient.");
179
181 "MaskEdgeTubes", true,
182 "Mask tubes on the edge of all banks in the PeaksWorkspace instrument (note the edge pixels at top/bottom of all "
183 "tubes will always be masked even if this property is False). Note the algorithm will treat "
184 "any masked pixels as edges (including pixels already masked prior to the execution of this algorithm) - this "
185 "means a custom mask can be applied to the PeaksWorkspace before integration.");
186
187 // Group Properties
188 std::string general_grp = "General Inputs";
189 std::string cylin_grp = "Cylindrical Integration";
190 std::string ellip_grp = "Ellipsoid Integration";
191
192 setPropertyGroup("InputWorkspace", general_grp);
193 setPropertyGroup("PeakRadius", general_grp);
194 setPropertyGroup("BackgroundInnerRadius", general_grp);
195 setPropertyGroup("BackgroundOuterRadius", general_grp);
196 setPropertyGroup("PeaksWorkspace", general_grp);
197 setPropertyGroup("OutputWorkspace", general_grp);
198 setPropertyGroup("ReplaceIntensity", general_grp);
199 setPropertyGroup("IntegrateIfOnEdge", general_grp);
200 setPropertyGroup("AdaptiveQBackground", general_grp);
201
202 setPropertyGroup("Ellipsoid", ellip_grp);
203 setPropertyGroup("FixQAxis", ellip_grp);
204
205 setPropertyGroup("Cylinder", cylin_grp);
206 setPropertyGroup("CylinderLength", cylin_grp);
207 setPropertyGroup("PercentBackground", cylin_grp);
208 setPropertyGroup("ProfileFunction", cylin_grp);
209 setPropertyGroup("IntegrationOption", cylin_grp);
210 setPropertyGroup("ProfilesFile", cylin_grp);
211
212 setPropertyGroup("AdaptiveQMultiplier", general_grp);
213 setPropertyGroup("CorrectIfOnEdge", general_grp);
214 setPropertyGroup("UseOnePercentBackgroundCorrection", general_grp);
215
216 setPropertyGroup("FixMajorAxisLength", ellip_grp);
217 setPropertyGroup("UseCentroid", ellip_grp);
218 setPropertyGroup("MaxIterations", ellip_grp);
219
220 setPropertyGroup("MaskEdgeTubes", general_grp);
221
222 // SetValue when another property value changes
223 setPropertySettings("Ellipsoid", std::make_unique<SetValueWhenProperty>(
224 "Cylinder", [](std::string ellipsoid, const std::string &cylinder) {
225 // Set Ellipsoid to 0, if user has set Cylinder to 1
226 if (ellipsoid == "1" && cylinder == "1") {
227 return std::string{"0"};
228 } else {
229 return ellipsoid;
230 };
231 }));
232 setPropertySettings("Cylinder", std::make_unique<SetValueWhenProperty>(
233 "Ellipsoid", [](std::string cylinder, const std::string &ellipsoid) {
234 // Set Cylinder to 0, if user has set Ellipsoid to 1
235 if (cylinder == "1" && ellipsoid == "1") {
236 return std::string{"0"};
237 } else {
238 return cylinder;
239 };
240 }));
241
242 // Set these Properties as visible only when Cylinder = 1
243 setPropertySettings("CylinderLength", std::make_unique<VisibleWhenProperty>("Cylinder", IS_EQUAL_TO, "1"));
244 setPropertySettings("PercentBackground", std::make_unique<VisibleWhenProperty>("Cylinder", IS_EQUAL_TO, "1"));
245 setPropertySettings("ProfileFunction", std::make_unique<VisibleWhenProperty>("Cylinder", IS_EQUAL_TO, "1"));
246 setPropertySettings("IntegrationOption", std::make_unique<VisibleWhenProperty>("Cylinder", IS_EQUAL_TO, "1"));
247 setPropertySettings("ProfilesFile", std::make_unique<VisibleWhenProperty>("Cylinder", IS_EQUAL_TO, "1"));
248
249 // Set these Properties as visible only when Ellipsoid = 1
250 setPropertySettings("FixQAxis", std::make_unique<VisibleWhenProperty>("Ellipsoid", IS_EQUAL_TO, "1"));
251 setPropertySettings("FixMajorAxisLength", std::make_unique<VisibleWhenProperty>("Ellipsoid", IS_EQUAL_TO, "1"));
252 setPropertySettings("UseCentroid", std::make_unique<VisibleWhenProperty>("Ellipsoid", IS_EQUAL_TO, "1"));
253 setPropertySettings("MaxIterations", std::make_unique<VisibleWhenProperty>("Ellipsoid", IS_EQUAL_TO, "1"));
254
255 // Disable / greyed out these Properties based on the value of another
256 setPropertySettings("CorrectIfOnEdge",
257 std::make_unique<Kernel::EnabledWhenProperty>("IntegrateIfOnEdge", IS_EQUAL_TO, "1"));
258}
259
260std::map<std::string, std::string> IntegratePeaksMD2::validateInputs() {
261 std::map<std::string, std::string> result;
262
263 std::vector<double> PeakRadius = getProperty("PeakRadius");
264 std::vector<double> BackgroundInnerRadius = getProperty("BackgroundInnerRadius");
265 std::vector<double> BackgroundOuterRadius = getProperty("BackgroundOuterRadius");
266 bool ellipsoid = getProperty("Ellipsoid");
267 bool cylinder = getProperty("Cylinder");
268 std::shared_ptr<IMDEventWorkspace> inputWS = getProperty("InputWorkspace");
269
270 if (PeakRadius.size() != 1 && PeakRadius.size() != 3) {
271 std::stringstream errmsg;
272 errmsg << "Only one or three values should be specified";
273 result["PeakRadius"] = errmsg.str();
274 }
275
276 if (!ellipsoid && PeakRadius.size() != 1) {
277 std::stringstream errmsg;
278 errmsg << "One value must be specified when Ellipsoid is false";
279 result["PeakRadius"] = errmsg.str();
280 }
281
282 if (BackgroundInnerRadius.size() != 1 && BackgroundInnerRadius.size() != 3) {
283 std::stringstream errmsg;
284 errmsg << "Only one or three values should be specified";
285 result["BackgroundInnerRadius"] = errmsg.str();
286 }
287
288 if (!ellipsoid && BackgroundInnerRadius.size() != 1) {
289 std::stringstream errmsg;
290 errmsg << "One value must be specified when Ellipsoid is false";
291 result["BackgroundInnerRadius"] = errmsg.str();
292 }
293
294 if (BackgroundOuterRadius.size() != 1 && BackgroundOuterRadius.size() != 3) {
295 std::stringstream errmsg;
296 errmsg << "Only one or three values should be specified";
297 result["BackgroundOuterRadius"] = errmsg.str();
298 }
299
300 if (!ellipsoid && BackgroundOuterRadius.size() != 1) {
301 std::stringstream errmsg;
302 errmsg << "One value must be specified when Ellipsoid is false";
303 result["BackgroundOuterRadius"] = errmsg.str();
304 }
305
306 if (ellipsoid && cylinder) {
307 std::stringstream errmsg;
308 errmsg << "Ellipsoid and Cylinder cannot both be true";
309 result["Ellipsoid"] = errmsg.str();
310 result["Cylinder"] = errmsg.str();
311 }
312
313 if (inputWS->getSpecialCoordinateSystem() == Mantid::Kernel::None) {
314 result["InputWorkspace"] = "Must have Special Coordinate System of QLab, QSample or HKL";
315 }
316
317 return result;
318}
319
320//----------------------------------------------------------------------------------------------
325template <typename MDE, size_t nd> void IntegratePeaksMD2::integrate(typename MDEventWorkspace<MDE, nd>::sptr ws) {
326 if (nd != 3)
327 throw std::invalid_argument("For now, we expect the input MDEventWorkspace "
328 "to have 3 dimensions only.");
329
331 IPeaksWorkspace_sptr inPeakWS = getProperty("PeaksWorkspace");
332
334 IPeaksWorkspace_sptr peakWS = getProperty("OutputWorkspace");
335 if (peakWS != inPeakWS)
336 peakWS = inPeakWS->clone();
337 // This only fails in the unit tests which say that MaskBTP is not registered
338 bool maskTubes = getProperty("MaskEdgeTubes");
339 try {
340 PeaksWorkspace_sptr p = std::dynamic_pointer_cast<PeaksWorkspace>(inPeakWS);
341 if (p) {
342 if (maskTubes) {
343 runMaskDetectors(p, "Tube", "edges");
344 }
345 runMaskDetectors(p, "Pixel", "edges");
346 }
347 } catch (...) {
348 g_log.error("Can't execute MaskBTP algorithm for this instrument to set "
349 "edge for IntegrateIfOnEdge option");
350 }
351
352 calculateE1(inPeakWS->detectorInfo()); // fill E1Vec for use in detectorQ
354
356 std::vector<double> PeakRadius = getProperty("PeakRadius");
358 std::vector<double> BackgroundOuterRadius = getProperty("BackgroundOuterRadius");
360 std::vector<double> BackgroundInnerRadius = getProperty("BackgroundInnerRadius");
362 bool useOnePercentBackgroundCorrection = getProperty("UseOnePercentBackgroundCorrection");
363
364 bool manualEllip = false;
365 if (PeakRadius.size() > 1) {
366 manualEllip = true;
367 // make sure the background radii are 3 values (they default to 1)
368 if (BackgroundInnerRadius.size() == 1)
369 BackgroundInnerRadius.resize(3, BackgroundInnerRadius[0]);
370 if (BackgroundOuterRadius.size() == 1)
371 BackgroundOuterRadius.resize(3, BackgroundOuterRadius[0]);
372 }
373
374 double minInnerRadius = PeakRadius[0];
375 for (size_t r = 0; r < BackgroundInnerRadius.size(); r++) {
376 if (manualEllip) {
377 minInnerRadius = PeakRadius[r];
378 }
379 if (BackgroundInnerRadius[r] < minInnerRadius)
380 BackgroundInnerRadius[r] = minInnerRadius;
381 }
382 // Ellipsoid
383 bool isEllipse = getProperty("Ellipsoid");
384 bool qAxisIsFixed = getProperty("FixQAxis");
385 bool majorAxisLengthFixed = getProperty("FixMajorAxisLength");
386 bool useCentroid = getProperty("UseCentroid");
387 int maxCovarIter = getProperty("MaxIterations");
389 double cylinderLength = getProperty("CylinderLength");
390 Workspace2D_sptr wsProfile2D, wsFit2D, wsDiff2D;
391 size_t numSteps = 0;
392 bool cylinderBool = getProperty("Cylinder");
393 bool adaptiveQBackground = getProperty("AdaptiveQBackground");
394 double adaptiveQMultiplier = getProperty("AdaptiveQMultiplier");
395 double adaptiveQBackgroundMultiplier = 0.0;
396 if (adaptiveQBackground)
397 adaptiveQBackgroundMultiplier = adaptiveQMultiplier;
398 std::vector<double> PeakRadiusVector(peakWS->getNumberPeaks(), PeakRadius[0]);
399 std::vector<double> BackgroundInnerRadiusVector(peakWS->getNumberPeaks(), BackgroundInnerRadius[0]);
400 std::vector<double> BackgroundOuterRadiusVector(peakWS->getNumberPeaks(), BackgroundOuterRadius[0]);
401 if (cylinderBool) {
402 numSteps = 100;
403 size_t histogramNumber = peakWS->getNumberPeaks();
404 Workspace_sptr wsProfile = WorkspaceFactory::Instance().create("Workspace2D", histogramNumber, numSteps, numSteps);
405 wsProfile2D = std::dynamic_pointer_cast<Workspace2D>(wsProfile);
406 AnalysisDataService::Instance().addOrReplace("ProfilesData", wsProfile2D);
407 Workspace_sptr wsFit = WorkspaceFactory::Instance().create("Workspace2D", histogramNumber, numSteps, numSteps);
408 wsFit2D = std::dynamic_pointer_cast<Workspace2D>(wsFit);
409 AnalysisDataService::Instance().addOrReplace("ProfilesFit", wsFit2D);
410 Workspace_sptr wsDiff = WorkspaceFactory::Instance().create("Workspace2D", histogramNumber, numSteps, numSteps);
411 wsDiff2D = std::dynamic_pointer_cast<Workspace2D>(wsDiff);
412 AnalysisDataService::Instance().addOrReplace("ProfilesFitDiff", wsDiff2D);
413 auto newAxis1 = std::make_unique<TextAxis>(peakWS->getNumberPeaks());
414 auto newAxis2 = std::make_unique<TextAxis>(peakWS->getNumberPeaks());
415 auto newAxis3 = std::make_unique<TextAxis>(peakWS->getNumberPeaks());
416 auto newAxis1Raw = newAxis1.get();
417 auto newAxis2Raw = newAxis2.get();
418 auto newAxis3Raw = newAxis3.get();
419 wsProfile2D->replaceAxis(1, std::move(newAxis1));
420 wsFit2D->replaceAxis(1, std::move(newAxis2));
421 wsDiff2D->replaceAxis(1, std::move(newAxis3));
422 for (int i = 0; i < peakWS->getNumberPeaks(); ++i) {
423 // Get a direct ref to that peak.
424 IPeak &p = peakWS->getPeak(i);
425 std::ostringstream label;
426 label << Utils::round(p.getH()) << "_" << Utils::round(p.getK()) << "_" << Utils::round(p.getL()) << "_"
427 << p.getRunNumber();
428 newAxis1Raw->setLabel(i, label.str());
429 newAxis2Raw->setLabel(i, label.str());
430 newAxis3Raw->setLabel(i, label.str());
431 }
432 }
433 double percentBackground = getProperty("PercentBackground");
434 size_t peakMin = 0;
435 size_t peakMax = numSteps;
436 double ratio = 0.0;
437 if (cylinderBool) {
438 peakMin = static_cast<size_t>(static_cast<double>(numSteps) * percentBackground / 100.);
439 peakMax = numSteps - peakMin - 1;
440 size_t numPeakCh = peakMax - peakMin + 1; // number of peak channels
441 size_t numBkgCh = numSteps - numPeakCh; // number of background channels
442 ratio = static_cast<double>(numPeakCh) / static_cast<double>(numBkgCh);
443 }
445 bool replaceIntensity = getProperty("ReplaceIntensity");
446 bool integrateEdge = getProperty("IntegrateIfOnEdge");
447 bool correctEdge = getProperty("CorrectIfOnEdge");
448
449 std::string profileFunction = getProperty("ProfileFunction");
450 std::string integrationOption = getProperty("IntegrationOption");
451 std::ofstream out;
452 if (cylinderBool && profileFunction != "NoFit") {
453 std::string outFile = getProperty("InputWorkspace");
454 outFile.append(profileFunction);
455 outFile.append(".dat");
456 std::string save_path = ConfigService::Instance().getString("defaultsave.directory");
457 outFile = save_path + outFile;
458 out.open(outFile.c_str(), std::ofstream::out);
459 }
460 // volume of Background sphere with inner volume subtracted
461 double volumeBkg = 4.0 / 3.0 * M_PI * (std::pow(BackgroundOuterRadius[0], 3) - std::pow(BackgroundOuterRadius[0], 3));
462 // volume of PeakRadius sphere
463 double volumeRadius = 4.0 / 3.0 * M_PI * std::pow(PeakRadius[0], 3);
464
465 // Initialize progress reporting
466 int nPeaks = peakWS->getNumberPeaks();
467 Progress progress(this, 0., 1., nPeaks);
468 int zeroHKLCount = 0;
469 bool doParallel = cylinderBool ? false : Kernel::threadSafe(*ws, *peakWS);
470 PARALLEL_SET_CONFIG_THREADS PRAGMA_OMP(parallel for if(doParallel) reduction(+:zeroHKLCount))
471 for (int i = 0; i < nPeaks; ++i) {
473 progress.report();
474
475 // Get a direct ref to that peak.
476 IPeak &p = peakWS->getPeak(i);
477
478 // Get the peak center as a position in the dimensions of the workspace
479 bool missingIndex = false; // True if in HKL mode and pos is 0,0,0
480 V3D pos;
481 if (CoordinatesToUse == Mantid::Kernel::QLab) //"Q (lab frame)"
482 pos = p.getQLabFrame();
483 else if (CoordinatesToUse == Mantid::Kernel::QSample) //"Q (sample frame)"
484 pos = p.getQSampleFrame();
485 else if (CoordinatesToUse == Mantid::Kernel::HKL) { //"HKL"
486 pos = p.getHKL();
487 if (pos.X() == 0 && pos.Y() == 0 && pos.Z() == 0) {
488 ++zeroHKLCount;
489 missingIndex = true;
490 }
491 } else {
492 throw std::runtime_error("Workspace does not have a coordinate system set");
493 }
494 // Do not integrate if sphere is off edge of detector
495
496 const double edgeDist = calculateDistanceToEdge(p.getQLabFrame());
497 if (edgeDist < std::max(BackgroundOuterRadius[0], PeakRadius[0])) {
498 g_log.warning() << "Warning: sphere/cylinder for integration is off edge "
499 "of detector for peak "
500 << i << "; radius of edge = " << edgeDist << '\n';
501 if (!integrateEdge) {
502 if (replaceIntensity) {
503 p.setIntensity(0.0);
504 p.setSigmaIntensity(0.0);
505 }
506 continue;
507 }
508 }
509
510 // Build the sphere transformation
511 bool dimensionsUsed[nd];
512 coord_t center[nd];
513 for (size_t d = 0; d < nd; ++d) {
514 dimensionsUsed[d] = true; // Use all dimensions
515 center[d] = static_cast<coord_t>(pos[d]);
516 }
517 signal_t signal = 0;
518 signal_t errorSquared = 0;
519 signal_t bgSignal = 0;
520 signal_t bgErrorSquared = 0;
521 double background_total = 0.0;
522 if (!cylinderBool) {
523 // modulus of Q
524 coord_t lenQpeak = 0.0;
525 if (adaptiveQMultiplier != 0.0) {
526 lenQpeak = 0.0;
527 for (size_t d = 0; d < nd; ++d) {
528 lenQpeak += center[d] * center[d];
529 }
530 lenQpeak = std::sqrt(lenQpeak);
531 }
532 double adaptiveRadius = adaptiveQMultiplier * lenQpeak + *std::max_element(PeakRadius.begin(), PeakRadius.end());
533 if (adaptiveRadius <= 0.0) {
534 g_log.error() << "Error: Radius for integration sphere of peak " << i << " is negative = " << adaptiveRadius
535 << '\n';
536 adaptiveRadius = 0.;
537 p.setIntensity(0.0);
538 p.setSigmaIntensity(0.0);
539 PeakRadiusVector[i] = 0.0;
540 BackgroundInnerRadiusVector[i] = 0.0;
541 BackgroundOuterRadiusVector[i] = 0.0;
542 continue;
543 }
544 PeakRadiusVector[i] = adaptiveRadius;
545 BackgroundInnerRadiusVector[i] = adaptiveQBackgroundMultiplier * lenQpeak +
546 *std::max_element(BackgroundInnerRadius.begin(), BackgroundInnerRadius.end());
547 BackgroundOuterRadiusVector[i] = adaptiveQBackgroundMultiplier * lenQpeak +
548 *std::max_element(BackgroundOuterRadius.begin(), BackgroundOuterRadius.end());
549 // define the radius squared for a sphere intially
550 CoordTransformDistance getRadiusSq(nd, center, dimensionsUsed);
551 // set spherical shape
552 PeakShape *sphereShape =
553 new PeakShapeSpherical(PeakRadiusVector[i], BackgroundInnerRadiusVector[i], BackgroundOuterRadiusVector[i],
554 CoordinatesToUse, this->name(), this->version());
555 p.setPeakShape(sphereShape);
556 const double scaleFactor = pow(PeakRadiusVector[i], 3) /
557 (pow(BackgroundOuterRadiusVector[i], 3) - pow(BackgroundInnerRadiusVector[i], 3));
558 // Integrate spherical background shell if specified
559 if (BackgroundOuterRadius[0] > PeakRadius[0]) {
560 // Get the total signal inside background shell
561 ws->getBox()->integrateSphere(
562 getRadiusSq, static_cast<coord_t>(pow(BackgroundOuterRadiusVector[i], 2)), bgSignal, bgErrorSquared,
563 static_cast<coord_t>(pow(BackgroundInnerRadiusVector[i], 2)), useOnePercentBackgroundCorrection);
564 // correct bg signal by Vpeak/Vshell (same for sphere and ellipse)
565 bgSignal *= scaleFactor;
566 bgErrorSquared *= scaleFactor * scaleFactor;
567 }
568 // if ellipsoid find covariance and centroid in spherical region
569 // using one-pass algorithm from https://doi.org/10.1145/359146.359153
570 bool integrateAsEllipse = isEllipse;
571 if (isEllipse && missingIndex) {
572 integrateAsEllipse = false;
573 g_log.warning() << "Integrating un-indexed peak. Falling back to sphere peak shape to avoid numeric issues\n";
574 }
575
576 if (integrateAsEllipse) {
577 // flat bg to subtract
578 const auto bgDensity = bgSignal / (4 * M_PI * pow(PeakRadiusVector[i], 3) / 3);
579 std::array<V3D, 3> eigenvects;
580 std::array<double, 3> eigenvals;
581 V3D translation(0.0, 0.0, 0.0); // translation from peak pos to centroid
582 if (PeakRadius.size() == 1) {
583 V3D mean(0.0, 0.0, 0.0); // vector to hold centroid
584 findEllipsoid<MDE, nd>(ws, getRadiusSq, pos, static_cast<coord_t>(pow(PeakRadiusVector[i], 2)), qAxisIsFixed,
585 useCentroid, bgDensity, eigenvects, eigenvals, mean, maxCovarIter);
586 if (!majorAxisLengthFixed) {
587 // replace radius for this peak with 3*stdev along major axis
588 auto max_stdev = sqrt(*std::max_element(eigenvals.begin(), eigenvals.end()));
589 BackgroundOuterRadiusVector[i] = 3 * max_stdev * (BackgroundOuterRadiusVector[i] / PeakRadiusVector[i]);
590 BackgroundInnerRadiusVector[i] = 3 * max_stdev * (BackgroundInnerRadiusVector[i] / PeakRadiusVector[i]);
591 PeakRadiusVector[i] = 3 * max_stdev;
592 }
593 if (useCentroid) {
594 // calculate translation to apply when drawing
595 translation = mean - pos;
596 // update integration center with mean
597 for (size_t d = 0; d < 3; ++d) {
598 center[d] = static_cast<coord_t>(mean[d]);
599 }
600 }
601 } else {
602 // Use the manually specified radii instead of finding them via
603 // findEllipsoid
604 std::transform(PeakRadius.cbegin(), PeakRadius.cend(), eigenvals.begin(),
605 [](double r) { return std::pow(r, 2.0); });
606 eigenvects[0] = V3D(1.0, 0.0, 0.0);
607 eigenvects[1] = V3D(0.0, 1.0, 0.0);
608 eigenvects[2] = V3D(0.0, 0.0, 1.0);
609 }
610 // transform ellispoid onto sphere of radius = R
611 getRadiusSq = CoordTransformDistance(nd, center, dimensionsUsed, 1, /* outD */
612 eigenvects, eigenvals);
613 // Integrate ellipsoid background shell if specified
614 if (PeakRadius.size() == 1) {
615 if (BackgroundOuterRadius[0] > PeakRadius[0]) {
616 // Get the total signal inside "BackgroundOuterRadius"
617 bgSignal = 0;
618 bgErrorSquared = 0;
619 ws->getBox()->integrateSphere(
620 getRadiusSq, static_cast<coord_t>(pow(BackgroundOuterRadiusVector[i], 2)), bgSignal, bgErrorSquared,
621 static_cast<coord_t>(pow(BackgroundInnerRadiusVector[i], 2)), useOnePercentBackgroundCorrection);
622 // correct bg signal by Vpeak/Vshell (same as previously
623 // calculated for sphere)
624 bgSignal *= scaleFactor;
625 bgErrorSquared *= scaleFactor * scaleFactor;
626 }
627 // set peak shape
628 // get radii in same proprtion as eigenvalues
629 auto max_stdev = pow(*std::max_element(eigenvals.begin(), eigenvals.end()), 0.5);
630 PeakEllipsoidExtent peakRadii{0., 0., 0.};
631 PeakEllipsoidExtent backgroundInnerRadii{0., 0., 0.};
632 PeakEllipsoidExtent backgroundOuterRadii{0., 0., 0.};
633 for (size_t irad = 0; irad < peakRadii.size(); irad++) {
634 auto scale = pow(eigenvals[irad], 0.5) / max_stdev;
635 peakRadii[irad] = PeakRadiusVector[i] * scale;
636 backgroundInnerRadii[irad] = BackgroundInnerRadiusVector[i] * scale;
637 backgroundOuterRadii[irad] = BackgroundOuterRadiusVector[i] * scale;
638 }
639 PeakShape *ellipsoidShape =
640 new PeakShapeEllipsoid(eigenvects, peakRadii, backgroundInnerRadii, backgroundOuterRadii,
641 CoordinatesToUse, this->name(), this->version(), translation);
642 p.setPeakShape(ellipsoidShape);
643 } else {
644 // Use the manually specified radii instead of finding them via
645 // findEllipsoid
646 std::array<double, 3> eigenvals_background_inner;
647 std::array<double, 3> eigenvals_background_outer;
648 std::transform(BackgroundInnerRadius.cbegin(), BackgroundInnerRadius.cend(),
649 eigenvals_background_inner.begin(), [](double r) { return std::pow(r, 2.0); });
650 std::transform(BackgroundOuterRadius.cbegin(), BackgroundOuterRadius.cend(),
651 eigenvals_background_outer.begin(), [](double r) { return std::pow(r, 2.0); });
652
653 if (BackgroundOuterRadiusVector[0] > PeakRadiusVector[0]) {
654 // transform ellispoid onto sphere of radius = R
655 auto getRadiusSqInner = CoordTransformDistance(nd, center, dimensionsUsed, 1, /* outD */
656 eigenvects, eigenvals_background_inner);
657 auto getRadiusSqOuter = CoordTransformDistance(nd, center, dimensionsUsed, 1, /* outD */
658 eigenvects, eigenvals_background_outer);
659 // Get the total signal inside "BackgroundOuterRadius"
660 bgSignal = 0;
661 bgErrorSquared = 0;
662 signal_t bgSignalInner = 0;
663 signal_t bgSignalOuter = 0;
664 signal_t bgErrorSquaredInner = 0;
665 signal_t bgErrorSquaredOuter = 0;
666 ws->getBox()->integrateSphere(getRadiusSqInner,
667 static_cast<coord_t>(pow(BackgroundInnerRadiusVector[i], 2)), bgSignalInner,
668 bgErrorSquaredInner, 0.0, useOnePercentBackgroundCorrection);
669 ws->getBox()->integrateSphere(getRadiusSqOuter,
670 static_cast<coord_t>(pow(BackgroundOuterRadiusVector[i], 2)), bgSignalOuter,
671 bgErrorSquaredOuter, 0.0, useOnePercentBackgroundCorrection);
672 // correct bg signal by Vpeak/Vshell (same as previously
673 // calculated for sphere)
674 bgSignal = bgSignalOuter - bgSignalInner;
675 bgErrorSquared = bgErrorSquaredInner + bgErrorSquaredOuter;
676 g_log.debug() << "unscaled background signal from ellipsoid integration = " << bgSignal << '\n';
677 const double scale = (PeakRadius[0] * PeakRadius[1] * PeakRadius[2]) /
678 (BackgroundOuterRadius[0] * BackgroundOuterRadius[1] * BackgroundOuterRadius[2] -
679 BackgroundInnerRadius[0] * BackgroundInnerRadius[1] * BackgroundInnerRadius[2]);
680 bgSignal *= scale;
681 bgErrorSquared *= pow(scale, 2);
682 }
683 // set peak shape
684 // get radii in same proprtion as eigenvalues
685 auto max_stdev = pow(*std::max_element(eigenvals.begin(), eigenvals.end()), 0.5);
686 auto max_stdev_inner =
687 pow(*std::max_element(eigenvals_background_inner.begin(), eigenvals_background_inner.end()), 0.5);
688 auto max_stdev_outer =
689 pow(*std::max_element(eigenvals_background_outer.begin(), eigenvals_background_outer.end()), 0.5);
690 PeakEllipsoidExtent peakRadii{0., 0., 0.};
691 PeakEllipsoidExtent backgroundInnerRadii{0., 0., 0.};
692 PeakEllipsoidExtent backgroundOuterRadii{0., 0., 0.};
693 for (size_t irad = 0; irad < peakRadii.size(); irad++) {
694 peakRadii[irad] = PeakRadiusVector[i] * pow(eigenvals[irad], 0.5) / max_stdev;
695 backgroundInnerRadii[irad] =
696 BackgroundInnerRadiusVector[i] * pow(eigenvals_background_inner[irad], 0.5) / max_stdev_inner;
697 backgroundOuterRadii[irad] =
698 BackgroundOuterRadiusVector[i] * pow(eigenvals_background_outer[irad], 0.5) / max_stdev_outer;
699 }
700 PeakShape *ellipsoidShape =
701 new PeakShapeEllipsoid(eigenvects, peakRadii, backgroundInnerRadii, backgroundOuterRadii,
702 CoordinatesToUse, this->name(), this->version());
703 p.setPeakShape(ellipsoidShape);
704 }
705 }
706 ws->getBox()->integrateSphere(getRadiusSq, static_cast<coord_t>(PeakRadiusVector[i] * PeakRadiusVector[i]),
707 signal, errorSquared, 0.0 /* innerRadiusSquared */,
708 useOnePercentBackgroundCorrection);
709 //
710 } else {
711 CoordTransformDistance cylinder(nd, center, dimensionsUsed, 2);
712
713 // Perform the integration into whatever box is contained within.
714 Counts signal_fit(numSteps);
715 signal_fit = 0;
716
717 ws->getBox()->integrateCylinder(cylinder, static_cast<coord_t>(PeakRadius[0]),
718 static_cast<coord_t>(cylinderLength), signal, errorSquared,
719 signal_fit.mutableRawData());
720
721 // Integrate around the background radius
722 if (BackgroundOuterRadius[0] > PeakRadius[0]) {
723 // Get the total signal inside "BackgroundOuterRadius"
724 signal_fit = 0;
725
726 ws->getBox()->integrateCylinder(cylinder, static_cast<coord_t>(BackgroundOuterRadius[0]),
727 static_cast<coord_t>(cylinderLength), bgSignal, bgErrorSquared,
728 signal_fit.mutableRawData());
729
730 Points points(signal_fit.size(), LinearGenerator(0, 1));
731 wsProfile2D->setHistogram(i, points, signal_fit);
732
733 // Evaluate the signal inside "BackgroundInnerRadius"
734 signal_t interiorSignal = 0;
735 signal_t interiorErrorSquared = 0;
736
737 // Integrate this 3rd radius, if needed
738 if (BackgroundInnerRadius[0] != PeakRadius[0]) {
739 ws->getBox()->integrateCylinder(cylinder, static_cast<coord_t>(BackgroundInnerRadius[0]),
740 static_cast<coord_t>(cylinderLength), interiorSignal, interiorErrorSquared,
741 signal_fit.mutableRawData());
742 } else {
743 // PeakRadius == BackgroundInnerRadius, so use the previous
744 // value
745 interiorSignal = signal;
746 interiorErrorSquared = errorSquared;
747 }
748 // Subtract the peak part to get the intensity in the shell
749 // (BackgroundInnerRadius < r < BackgroundOuterRadius)
750 bgSignal -= interiorSignal;
751 // We can subtract the error (instead of adding) because the two
752 // values are 100% dependent; this is the same as integrating a
753 // shell.
754 bgErrorSquared -= interiorErrorSquared;
755 // Relative volume of peak vs the BackgroundOuterRadius cylinder
756 const double radiusRatio = (PeakRadius[0] / BackgroundOuterRadius[0]);
757 const double peakVolume = radiusRatio * radiusRatio * cylinderLength;
758
759 // Relative volume of the interior of the shell vs overall
760 // background
761 const double interiorRatio = (BackgroundInnerRadius[0] / BackgroundOuterRadius[0]);
762 // Volume of the bg shell, relative to the volume of the
763 // BackgroundOuterRadius cylinder
764 const double bgVolume = 1.0 - interiorRatio * interiorRatio * cylinderLength;
765
766 // Finally, you will multiply the bg intensity by this to get the
767 // estimated background under the peak volume
768 const double scaleFactor = peakVolume / bgVolume;
769 bgSignal *= scaleFactor;
770 bgErrorSquared *= scaleFactor * scaleFactor;
771 } else {
772 Points points(signal_fit.size(), LinearGenerator(0, 1));
773 wsProfile2D->setHistogram(i, points, signal_fit);
774 }
775
776 if (profileFunction == "NoFit") {
777 signal = 0.;
778 for (size_t j = 0; j < numSteps; j++) {
779 if (j < peakMin || j > peakMax)
780 background_total = background_total + wsProfile2D->y(i)[j];
781 else
782 signal = signal + wsProfile2D->y(i)[j];
783 }
784 errorSquared = std::fabs(signal);
785 } else {
786
787 auto fitAlgorithm = createChildAlgorithm("Fit", -1, -1, false);
788 // fitAlgorithm->setProperty("CreateOutput", true);
789 // fitAlgorithm->setProperty("Output", "FitPeaks1D");
790 std::string myFunc = std::string("name=LinearBackground;name=") + profileFunction;
791 auto maxPeak = std::max_element(signal_fit.begin(), signal_fit.end());
792
793 std::ostringstream strs;
794 strs << maxPeak[0];
795 std::string strMax = strs.str();
796 if (profileFunction == "Gaussian") {
797 myFunc += ", PeakCentre=50, Height=" + strMax;
798 fitAlgorithm->setProperty("Constraints", "40<f1.PeakCentre<60");
799 } else if (profileFunction == "BackToBackExponential" || profileFunction == "IkedaCarpenterPV") {
800 myFunc += ", X0=50, I=" + strMax;
801 fitAlgorithm->setProperty("Constraints", "40<f1.X0<60");
802 }
803 fitAlgorithm->setProperty("CalcErrors", true);
804 fitAlgorithm->setProperty("Function", myFunc);
805 fitAlgorithm->setProperty("InputWorkspace", wsProfile2D);
806 fitAlgorithm->setProperty("WorkspaceIndex", static_cast<int>(i));
807 try {
808 fitAlgorithm->executeAsChildAlg();
809 } catch (...) {
810 g_log.error("Can't execute Fit algorithm");
811 continue;
812 }
813
814 IFunction_sptr ifun = fitAlgorithm->getProperty("Function");
815 if (i == 0) {
816 out << std::setw(20) << "spectrum"
817 << " ";
818 for (size_t j = 0; j < ifun->nParams(); ++j)
819 out << std::setw(20) << ifun->parameterName(j) << " ";
820 out << std::setw(20) << "chi2"
821 << " ";
822 out << "\n";
823 }
824 out << std::setw(20) << i << " ";
825 for (size_t j = 0; j < ifun->nParams(); ++j) {
826 out << std::setw(20) << std::fixed << std::setprecision(10) << ifun->getParameter(j) << " ";
827 }
828 double chi2 = fitAlgorithm->getProperty("OutputChi2overDoF");
829 out << std::setw(20) << std::fixed << std::setprecision(10) << chi2 << "\n";
830
831 std::shared_ptr<const CompositeFunction> fun = std::dynamic_pointer_cast<const CompositeFunction>(ifun);
832
833 const auto &x = wsProfile2D->x(i);
834 wsFit2D->setSharedX(i, wsProfile2D->sharedX(i));
835 wsDiff2D->setSharedX(i, wsProfile2D->sharedX(i));
836
837 FunctionDomain1DVector domain(x.rawData());
838 FunctionValues yy(domain);
839 fun->function(domain, yy);
840 auto funcValues = yy.toVector();
841
842 wsFit2D->mutableY(i) = funcValues;
843 wsDiff2D->setSharedY(i, wsProfile2D->sharedY(i));
844 wsDiff2D->mutableY(i) -= wsFit2D->y(i);
845
846 // Calculate intensity
847 signal = 0.0;
848 if (integrationOption == "Sum") {
849 for (size_t j = peakMin; j <= peakMax; j++)
850 if (std::isfinite(yy[j]))
851 signal += yy[j];
852 } else {
853 gsl_integration_workspace *w = gsl_integration_workspace_alloc(1000);
854
855 double error;
856
857 gsl_function F;
858 F.function = &Mantid::MDAlgorithms::f_eval2;
859 F.params = &fun;
860
861 gsl_integration_qags(&F, x[peakMin], x[peakMax], 0, 1e-7, 1000, w, &signal, &error);
862
863 gsl_integration_workspace_free(w);
864 }
865 errorSquared = std::fabs(signal);
866 // Get background counts
867 for (size_t j = 0; j < numSteps; j++) {
868 double background = ifun->getParameter(0) + ifun->getParameter(1) * x[j];
869 if (j < peakMin || j > peakMax)
870 background_total = background_total + background;
871 }
872 }
873 }
874 checkOverlap(i, peakWS, CoordinatesToUse, 2.0 * std::max(PeakRadiusVector[i], BackgroundOuterRadiusVector[i]));
875 // Save it back in the peak object.
876 if (signal != 0. || replaceIntensity) {
877 double edgeMultiplier = 1.0;
878 double peakMultiplier = 1.0;
879 if (correctEdge) {
880 if (edgeDist < BackgroundOuterRadius[0]) {
881 double e1 = BackgroundOuterRadius[0] - edgeDist;
882 // volume of cap of sphere with h = edge
883 double f1 = M_PI * std::pow(e1, 2) / 3 * (3 * BackgroundOuterRadius[0] - e1);
884 edgeMultiplier = volumeBkg / (volumeBkg - f1);
885 }
886 if (edgeDist < PeakRadius[0]) {
887 double sigma = PeakRadius[0] / 3.0;
888 // assume gaussian peak
889 double e1 = std::exp(-std::pow(edgeDist, 2) / (2 * sigma * sigma)) * PeakRadius[0];
890 // volume of cap of sphere with h = edge
891 double f1 = M_PI * std::pow(e1, 2) / 3 * (3 * PeakRadius[0] - e1);
892 peakMultiplier = volumeRadius / (volumeRadius - f1);
893 }
894 }
895
896 p.setIntensity(peakMultiplier * signal - edgeMultiplier * (ratio * background_total + bgSignal));
897 p.setSigmaIntensity(sqrt(peakMultiplier * errorSquared +
898 edgeMultiplier * (ratio * ratio * std::fabs(background_total) + bgErrorSquared)));
899 }
900
901 g_log.information() << "Peak " << i << " at " << pos << ": signal " << signal << " (sig^2 " << errorSquared
902 << "), with background " << bgSignal + ratio * background_total << " (sig^2 "
903 << bgErrorSquared + ratio * ratio * std::fabs(background_total) << ") subtracted.\n";
905 }
907
908 if (CoordinatesToUse == Kernel::HKL && zeroHKLCount > 0) {
909 if (zeroHKLCount == nPeaks) {
910 g_log.error() << "No peaks are indexed, integrated at HKL=0,0,0. You might need to run IndexPeaks.\n";
911 } else {
912 g_log.warning() << zeroHKLCount << " of " << nPeaks << " peaks are not indexed, integrated at HKL=0,0,0\n";
913 }
914 }
915
916 // This flag is used by the PeaksWorkspace to evaluate whether it has
917 // been integrated.
918 peakWS->mutableRun().addProperty("PeaksIntegrated", 1, true);
919 // These flags are specific to the algorithm.
920 peakWS->mutableRun().addProperty("PeakRadius", PeakRadiusVector, true);
921 peakWS->mutableRun().addProperty("BackgroundInnerRadius", BackgroundInnerRadiusVector, true);
922 peakWS->mutableRun().addProperty("BackgroundOuterRadius", BackgroundOuterRadiusVector, true);
923
924 // save profiles in peaks file
925 const std::string outfile = getProperty("ProfilesFile");
926 if (outfile.length() > 0) {
927 IAlgorithm_sptr alg;
928 try {
929 alg = createChildAlgorithm("SaveIsawPeaks", -1, -1, false);
930 } catch (Exception::NotFoundError &) {
931 g_log.error("Can't locate SaveIsawPeaks algorithm");
932 throw;
933 }
934 alg->setProperty("InputWorkspace", peakWS);
935 alg->setProperty("ProfileWorkspace", wsProfile2D);
936 alg->setPropertyValue("Filename", outfile);
937 alg->execute();
938 }
939 // Save the output
940 setProperty("OutputWorkspace", peakWS);
941}
942
960template <typename MDE, size_t nd>
962 const CoordTransform &getRadiusSq, const V3D &pos, const coord_t &radiusSquared,
963 const bool &qAxisIsFixed, const bool &useCentroid, const double &bgDensity,
964 std::array<V3D, 3> &eigenvects, std::array<double, 3> &eigenvals, V3D &mean,
965 const int maxIter) {
966
967 // get leaf-only iterators over all boxes in ws
968 auto function = std::make_unique<Geometry::MDAlgorithms::MDBoxMaskFunction>(pos, radiusSquared);
969 MDBoxBase<MDE, nd> *baseBox = ws->getBox();
970 MDBoxIterator<MDE, nd> MDiter(baseBox, 1000, true, function.get());
971
972 // get initial vector of events inside sphere
973 std::vector<std::pair<V3D, double>> peak_events;
974
975 do {
976 auto *box = dynamic_cast<MDBox<MDE, nd> *>(MDiter.getBox());
977 if (box && !box->getIsMasked()) {
978
979 // simple check whether box is defintely not contained
980 coord_t boxCenter[nd];
981 box->getCenter(boxCenter);
982 V3D displacement; // vector between peak pos and box center
983 coord_t rboxSq = 0; // dist from center to vertex sq
984 for (size_t d = 0; d < nd; ++d) {
985 auto dim = box->getExtents(d);
986 rboxSq += static_cast<coord_t>(0.25 * dim.getSize() * dim.getSize());
987 displacement[d] = pos[d] - static_cast<double>(boxCenter[d]);
988 }
989
990 if (displacement.norm() < static_cast<double>(sqrt(rboxSq)) + static_cast<double>(sqrt(radiusSquared))) {
991 // box MIGHT intersect peak spherical region so go through events
992 const std::vector<MDE> &events = box->getConstEvents();
993 auto bg = bgDensity / (static_cast<double>(events.size()) * (box->getInverseVolume()));
994 // For each event
995 for (const auto &evnt : events) {
996
997 coord_t center_array[nd];
998 for (size_t d = 0; d < nd; ++d) {
999 center_array[d] = evnt.getCenter(d);
1000 }
1001 coord_t out[1];
1002 auto *cen_ptr = center_array; // pointer to first element
1003 getRadiusSq.apply(cen_ptr, out);
1004
1005 if (evnt.getSignal() > bg && out[0] < radiusSquared) {
1006 // need in V3D for matrix maths later
1007 V3D center;
1008 for (size_t d = 0; d < nd; ++d) {
1009 center[d] = static_cast<double>(center_array[d]);
1010 }
1011 peak_events.emplace_back(center, evnt.getSignal() - bg);
1012 }
1013 }
1014 }
1015 }
1016 box->releaseEvents();
1017 } while (MDiter.next());
1018 calcCovar(peak_events, pos, radiusSquared, qAxisIsFixed, useCentroid, eigenvects, eigenvals, mean, maxIter);
1019}
1020
1036void IntegratePeaksMD2::calcCovar(const std::vector<std::pair<V3D, double>> &peak_events, const V3D &pos,
1037 const coord_t &radiusSquared, const bool &qAxisIsFixed, const bool &useCentroid,
1038 std::array<V3D, 3> &eigenvects, std::array<double, 3> &eigenvals, V3D &mean,
1039 const int maxIter) {
1040
1041 constexpr size_t num_dimensions = 3;
1042
1043 // to calc threshold mdsq to exclude events over 3 stdevs away
1044 boost::math::chi_squared chisq(static_cast<double>(num_dimensions));
1045 auto mdsq_max = boost::math::quantile(chisq, 0.997);
1046 Matrix<double> invCov; // required to calc mdsq
1047 double prev_cov_det = DBL_MAX;
1048
1049 // initialise mean with pos
1050 mean = pos;
1051 Matrix<double> Pinv(num_dimensions, num_dimensions);
1052 if (qAxisIsFixed) {
1053 // transformation from Qlab to Qhat, vhat and uhat,
1054 getPinv(pos, Pinv);
1055 mean = Pinv * mean;
1056 }
1057 Matrix<double> cov_mat(num_dimensions, num_dimensions);
1058
1059 for (int nIter = 0; nIter < maxIter; nIter++) {
1060
1061 // reset on each loop
1062 cov_mat.zeroMatrix();
1063 double w_sum = 0; // sum of weights
1064 size_t nmasked = 0; // num masked events outside 3stdevs
1065 auto prev_pos = mean;
1066
1067 for (size_t ievent = 0; ievent < peak_events.size(); ievent++) {
1068
1069 const auto event = peak_events[ievent];
1070 auto center = event.first;
1071 if (qAxisIsFixed) {
1072 // transform coords to Q, uhat, vhat basis
1073 center = Pinv * center;
1074 }
1075
1076 bool useEvent = true;
1077 if (nIter > 0) {
1078 // check if point within 3 stdevs of mean (in MD frame)
1079 // prev_pos is the mean if useCentroid
1080 const auto displ = center - prev_pos;
1081 auto mdsq = displ.scalar_prod(invCov * displ);
1082 if (mdsq > mdsq_max) {
1083 // exclude points outside 3 stdevs
1084 useEvent = false;
1085 nmasked += 1;
1086 }
1087 }
1088
1089 // if prev cov_mat chec
1090 if (useEvent) {
1091 const auto signal = event.second;
1092 w_sum += signal;
1093
1094 if (useCentroid) {
1095 // update mean
1096 mean += (center - mean) * (signal / w_sum);
1097 }
1098
1099 // weight for variance
1100 auto wi = signal * (w_sum - signal) / w_sum;
1101 size_t istart = 0;
1102 if (qAxisIsFixed) {
1103 // variance along Q (skipped in next nested loops below)
1104 cov_mat[0][0] += wi * pow((center[0] - mean[0]), 2);
1105 istart = 1;
1106 }
1107 for (size_t row = istart; row < cov_mat.numRows(); ++row) {
1108 for (size_t col = istart; col < cov_mat.numRows(); ++col) {
1109 // symmeteric matrix
1110 if (row <= col) {
1111 auto cov = wi * (center[row] - mean[row]) * (center[col] - mean[col]);
1112 if (row == col) {
1113 cov_mat[row][col] += cov;
1114 } else {
1115 cov_mat[row][col] += cov;
1116 cov_mat[col][row] += cov;
1117 }
1118 }
1119 }
1120 }
1121 }
1122 }
1123 // normalise the covariance matrix
1124 cov_mat /= w_sum; // normalise by sum of weights
1125
1126 // check if another iteration is required
1127 bool anyMasked = (nIter > 0) ? (nmasked > 0) : true;
1128 // check if ellipsoid volume greater than sphere
1129 auto cov_det = cov_mat.determinant();
1130 bool isEllipVolGreater = cov_det > pow(static_cast<double>(radiusSquared / 9), 3);
1131 // check for convergence of variances
1132 bool isConverged = (cov_det > 0.95 * prev_cov_det);
1133
1134 if (!anyMasked || isEllipVolGreater || isConverged) {
1135 break;
1136 } else {
1137 prev_cov_det = cov_det;
1138 // required to eval Mahalanobis distance
1139 invCov = Matrix<double>(cov_mat);
1140 invCov.Invert();
1141 }
1142 }
1143
1144 if (qAxisIsFixed) {
1145 // transform back to MD basis
1146 Matrix<double> P(Pinv);
1147 P.Transpose();
1148 mean = P * mean;
1149 cov_mat = P * cov_mat * Pinv;
1150 }
1151 Matrix<double> evecs; // hold eigenvectors
1152 Matrix<double> evals; // hold eigenvals in diag
1153 cov_mat.Diagonalise(evecs, evals);
1154
1155 auto min_eval = evals[0][0];
1156 for (size_t d = 1; d < num_dimensions; ++d) {
1157 min_eval = std::min(min_eval, evals[d][d]);
1158 }
1159 if (min_eval > static_cast<double>(radiusSquared / 9)) {
1160 // haven't found good covar - set to spherical region
1161 evals.identityMatrix();
1162 evals = evals * (static_cast<double>(radiusSquared) / 9);
1163 g_log.warning() << "Covariance of peak at ";
1164 pos.printSelf(g_log.warning());
1165 g_log.warning() << " is not well constrained, it has been set to spherical" << std::endl;
1166 }
1167
1168 // convert to vectors for output
1169 std::copy_n(evals.Diagonal().cbegin(), num_dimensions, eigenvals.begin());
1170 // set min eigenval to be small but non-zero (1e-6)
1171 // when no discernible peak above background
1172 std::replace_if(eigenvals.begin(), eigenvals.end(), [&](auto x) { return x < 1e-6; }, 1e-6);
1173
1174 // populate V3D vector of eigenvects (needed for ellipsoid shape)
1175 for (size_t ivect = 0; ivect < num_dimensions; ++ivect) {
1176 eigenvects[ivect] = V3D(evecs[0][ivect], evecs[1][ivect], evecs[2][ivect]);
1177 }
1178}
1179
1190 // loop over 3 mutually-orthogonal vectors until get one with
1191 // a component perp to Q (within tolerance)
1192 double dotprod = 1;
1193 size_t ii = 0;
1194 V3D qhat = q / q.norm();
1195 V3D tmp;
1196 do {
1197 tmp = V3D(0, 0, 0); // reset u
1198 tmp[ii] = 1.0;
1199 dotprod = qhat.scalar_prod(tmp);
1200 ii++;
1201 } while (abs(dotprod) > 1.0 - 1e-6);
1202 // populate Pinv with basis vector rows
1203 Pinv.setRow(0, qhat);
1204 tmp = qhat.cross_prod(tmp);
1205 Pinv.setRow(1, tmp / tmp.norm());
1206 tmp = qhat.cross_prod(tmp);
1207 Pinv.setRow(2, tmp / tmp.norm());
1208}
1209
1210/*
1211 * Define edges for each instrument by masking. For CORELLI, tubes 1 and
1212 *16, and pixels 0 and 255. Get Q in the lab frame for every peak, call it
1213 *C For every point on the edge, the trajectory in reciprocal space is a
1214 *straight line, going through O=V3D(0,0,0). Calculate a point at a fixed
1215 *momentum, say k=1. Q in the lab frame
1216 *E=V3D(-k*sin(tt)*cos(ph),-k*sin(tt)*sin(ph),k-k*cos(ph)).
1217 * Normalize E to 1: E=E*(1./E.norm())
1218 *
1219 * @param inst: instrument
1220 */
1222 for (size_t i = 0; i < detectorInfo.size(); ++i) {
1223 if (detectorInfo.isMonitor(i))
1224 continue; // skip monitor
1225 if (!detectorInfo.isMasked(i))
1226 continue; // edge is masked so don't check if not masked
1227 const auto &det = detectorInfo.detector(i);
1228 double tt1 = det.getTwoTheta(V3D(0, 0, 0), V3D(0, 0, 1)); // two theta
1229 double ph1 = det.getPhi(); // phi
1230 V3D E1 = V3D(-std::sin(tt1) * std::cos(ph1), -std::sin(tt1) * std::sin(ph1),
1231 1. - std::cos(tt1)); // end of trajectory
1232 E1 = E1 * (1. / E1.norm()); // normalize
1233 E1Vec.emplace_back(E1);
1234 }
1235}
1236
1247 double edgeDist = DBL_MAX;
1248 for (auto &E1 : E1Vec) {
1249 V3D distv = QLabFrame - E1 * (QLabFrame.scalar_prod(E1)); // distance to the
1250 // trajectory as a
1251 // vector
1252 edgeDist = std::min(edgeDist, distv.norm()); // want smallest dist to peak
1253 }
1254 return edgeDist;
1255}
1256
1258 const std::string &property, const std::string &values) {
1259 // For CORELLI do not count as edge if next to another detector bank
1260 if (property == "Tube" && peakWS->getInstrument()->getName() == "CORELLI") {
1261 auto alg = createChildAlgorithm("MaskBTP");
1262 alg->setProperty<Workspace_sptr>("Workspace", peakWS);
1263 alg->setProperty("Bank", "1,7,12,17,22,27,30,59,63,69,74,79,84,89");
1264 alg->setProperty(property, "1");
1265 if (!alg->execute())
1266 throw std::runtime_error("MaskDetectors Child Algorithm has not executed successfully");
1267 auto alg2 = createChildAlgorithm("MaskBTP");
1268 alg2->setProperty<Workspace_sptr>("Workspace", peakWS);
1269 alg2->setProperty("Bank", "6,11,16,21,26,29,58,62,68,73,78,83,88,91");
1270 alg2->setProperty(property, "16");
1271 if (!alg2->execute())
1272 throw std::runtime_error("MaskDetectors Child Algorithm has not executed successfully");
1273 } else {
1274 auto alg = createChildAlgorithm("MaskBTP");
1275 alg->setProperty<Workspace_sptr>("Workspace", peakWS);
1276 alg->setProperty(property, values);
1277 if (!alg->execute())
1278 throw std::runtime_error("MaskDetectors Child Algorithm has not executed successfully");
1279 }
1280}
1281
1283 Mantid::Kernel::SpecialCoordinateSystem CoordinatesToUse, double radius) {
1284 // Get a direct ref to that peak.
1285 const IPeak &p1 = peakWS->getPeak(i);
1286 V3D pos1;
1287 if (CoordinatesToUse == Kernel::QLab) //"Q (lab frame)"
1288 pos1 = p1.getQLabFrame();
1289 else if (CoordinatesToUse == Kernel::QSample) //"Q (sample frame)"
1290 pos1 = p1.getQSampleFrame();
1291 else if (CoordinatesToUse == Kernel::HKL) //"HKL"
1292 pos1 = p1.getHKL();
1293 for (int j = i + 1; j < peakWS->getNumberPeaks(); ++j) {
1294 // Get a direct ref to rest of peaks peak.
1295 const IPeak &p2 = peakWS->getPeak(j);
1296 V3D pos2;
1297 if (CoordinatesToUse == Kernel::QLab) //"Q (lab frame)"
1298 pos2 = p2.getQLabFrame();
1299 else if (CoordinatesToUse == Kernel::QSample) //"Q (sample frame)"
1300 pos2 = p2.getQSampleFrame();
1301 else if (CoordinatesToUse == Kernel::HKL) //"HKL"
1302 pos2 = p2.getHKL();
1303 if (pos1.distance(pos2) < radius) {
1304 g_log.warning() << " Warning: Peak integration spheres for peaks " << i << " and " << j
1305 << " overlap. Distance between peaks is " << pos1.distance(pos2) << '\n';
1306 }
1307 }
1308}
1309
1310//----------------------------------------------------------------------------------------------
1314 inWS = getProperty("InputWorkspace");
1315
1317}
1318
1319double f_eval2(double x, void *params) {
1320 std::shared_ptr<const API::CompositeFunction> fun =
1321 *reinterpret_cast<std::shared_ptr<const API::CompositeFunction> *>(params);
1322 FunctionDomain1DVector domain(x);
1323 FunctionValues yval(domain);
1324 fun->function(domain, yval);
1325 return yval[0];
1326}
1327
1328} // namespace Mantid::MDAlgorithms
#define DECLARE_ALGORITHM(classname)
Definition Algorithm.h:538
gsl_vector * tmp
double sigma
Definition GetAllEi.cpp:156
double error
#define CALL_MDEVENT_FUNCTION(funcname, workspace)
Macro that makes it possible to call a templated method for a MDEventWorkspace using a IMDEventWorksp...
#define PARALLEL_SET_CONFIG_THREADS
#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.
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.
Unique SingleValueParameter Declaration for InputNDimensions.
virtual void apply(const coord_t *inputVector, coord_t *outVector) const =0
@ OptionalSave
to specify a file to write to but an empty string is
Implements FunctionDomain1D with its own storage in form of a std::vector.
A class to store values calculated by a function.
const std::vector< double > & toVector() const
Return the calculated values as a vector.
An interface to a peak function, which extend the interface of IFunctionWithLocation by adding method...
Helper class for reporting progress from algorithms.
Definition Progress.h:25
A property class for workspaces.
Unique CoordCenterVectorParam type declaration for ndimensional coordinate centers.
Templated super-class of a multi-dimensional event "box".
Definition MDBoxBase.h:50
void getCenter(coord_t *const center) const override
Get the center of the box.
Definition MDBoxBase.h:210
void integrateSphere(Mantid::API::CoordTransform &radiusTransform, const coord_t radiusSquared, signal_t &signal, signal_t &errorSquared, const coord_t innerRadiusSquared=0.0, const bool useOnePercentBackgroundCorrection=true) const override=0
Sphere (peak) integration.
void integrateCylinder(Mantid::API::CoordTransform &radiusTransform, const coord_t radius, const coord_t length, signal_t &signal, signal_t &errorSquared, std::vector< signal_t > &signal_fit) const override=0
Cylinder (peak) integration.
MDBoxIterator: iterate through MDBoxBase hierarchy down to a given maximum depth.
bool next() override
Advance to the next cell.
MDBoxBase< MDE, nd > * getBox() const
Return a pointer to the current box pointed to by the iterator.
Templated class for a multi-dimensional event "box".
Definition MDBox.h:45
std::shared_ptr< MDEventWorkspace< MDE, nd > > sptr
Typedef for a shared pointer of this kind of event workspace.
Kernel::SpecialCoordinateSystem getSpecialCoordinateSystem() const override
Get the coordinate system.
PeakShapeEllipsoid : PeakShape representing a 3D ellipsoid.
PeakShapeSpherical : PeakShape for a spherical peak.
Geometry::DetectorInfo is an intermediate step towards a DetectorInfo that is part of Instrument-2....
bool isMasked(const size_t index) const
Returns true if the detector is masked.
const Geometry::IDetector & detector(const size_t index) const
Return a const reference to the detector with given index.
size_t size() const
Returns the size of the DetectorInfo, i.e., the number of detectors in the instrument.
bool isMonitor(const size_t index) const
Returns true if the detector is a monitor.
virtual double getTwoTheta(const Kernel::V3D &observer, const Kernel::V3D &axis) const =0
Gives the angle of this detector object with respect to an axis.
Structure describing a single-crystal peak.
Definition IPeak.h:26
virtual double getH() const =0
virtual void setIntensity(double m_Intensity)=0
virtual void setPeakShape(Mantid::Geometry::PeakShape *shape)=0
virtual double getK() const =0
virtual int getRunNumber() const =0
virtual Mantid::Kernel::V3D getQSampleFrame() const =0
virtual double getL() const =0
virtual void setSigmaIntensity(double m_SigmaIntensity)=0
virtual Mantid::Kernel::V3D getQLabFrame() const =0
virtual Mantid::Kernel::V3D getHKL() const =0
PeakShape : Abstract type to describes the shape of a peak.
Definition PeakShape.h:20
Support for a property that holds an array of values.
Exception for when an item is not found in a collection.
Definition Exception.h:145
IPropertyManager * setProperty(const std::string &name, const T &value)
Templated method to set the value of a PropertyWithValue.
void setPropertySettings(const std::string &name, std::unique_ptr< IPropertySettings const > settings)
Add a PropertySettings instance to the chain of settings for a given property.
void setPropertyGroup(const std::string &name, const std::string &group)
Set the group for a given property.
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
Numerical Matrix class.
Definition Matrix.h:42
T determinant() const
Calculate the determinant.
Definition Matrix.cpp:1048
T Invert()
LU inversion routine.
Definition Matrix.cpp:924
void zeroMatrix()
Set the matrix to zero.
Definition Matrix.cpp:646
void identityMatrix()
Makes the matrix an idenity matrix.
Definition Matrix.cpp:661
int Diagonalise(Matrix< T > &, Matrix< T > &) const
(only Symmetric matrix)
Definition Matrix.cpp:1319
void setRow(const size_t nRow, const std::vector< T > &newRow)
Definition Matrix.cpp:686
std::vector< T > Diagonal() const
Returns a vector of the diagonal.
Definition Matrix.cpp:1265
size_t numRows() const
Return the number of rows in the matrix.
Definition Matrix.h:144
Matrix< T > & Transpose()
Transpose the matrix.
Definition Matrix.cpp:793
The concrete, templated class for properties.
static T & Instance()
Return a reference to the Singleton instance, creating it if it does not already exist Creation is do...
Class for 3D vectors.
Definition V3D.h:34
double distance(const V3D &v) const noexcept
Calculates the distance between two vectors.
Definition V3D.h:293
constexpr double scalar_prod(const V3D &v) const noexcept
Calculates the cross product.
Definition V3D.h:280
void printSelf(std::ostream &) const
Prints a text representation of itself in format "[x,y,z]".
Definition V3D.cpp:349
constexpr double X() const noexcept
Get x.
Definition V3D.h:238
constexpr V3D cross_prod(const V3D &v) const noexcept
Cross product (this * argument)
Definition V3D.h:284
constexpr double Y() const noexcept
Get y.
Definition V3D.h:239
double norm() const noexcept
Definition V3D.h:269
constexpr double Z() const noexcept
Get z.
Definition V3D.h:240
std::vector< Kernel::V3D > E1Vec
save for all detector pixels
void runMaskDetectors(const Mantid::DataObjects::PeaksWorkspace_sptr &peakWS, const std::string &property, const std::string &values)
void calculateE1(const Geometry::DetectorInfo &detectorInfo)
Calculate if this Q is on a detector.
void getPinv(const Mantid::Kernel::V3D &q, Mantid::Kernel::Matrix< double > &Pinv)
Get the inverse of the matrix P.
void exec() override
Run the algorithm.
void findEllipsoid(const typename DataObjects::MDEventWorkspace< MDE, nd >::sptr ws, const Mantid::API::CoordTransform &getRadiusSq, const Mantid::Kernel::V3D &pos, const coord_t &radiusSquared, const bool &qAxisBool, const bool &useCentroid, const double &bgDensity, std::array< Mantid::Kernel::V3D, 3 > &eigenvects, std::array< double, 3 > &eigenvals, Mantid::Kernel::V3D &mean, const int maxIter=1)
Calculate the covariance matrix of a spherical region and store the eigenvectors and eigenvalues that...
void calcCovar(const std::vector< std::pair< Mantid::Kernel::V3D, double > > &peak_events, const Mantid::Kernel::V3D &pos, const coord_t &radiusSquared, const bool &qAxisIsFixed, const bool &useCentroid, std::array< Mantid::Kernel::V3D, 3 > &eigenvects, std::array< double, 3 > &eigenvals, Mantid::Kernel::V3D &mean, const int maxIter)
Calculate the covariance matrix of a spherical region and store the eigenvectors and eigenvalues that...
double calculateDistanceToEdge(const Mantid::Kernel::V3D &QLabFrame)
Calculate if this Q is on a detector The distance from C to OE is given by dv=C-E*(C....
int version() const override
Algorithm's version for identification.
void checkOverlap(int i, const Mantid::API::IPeaksWorkspace_sptr &peakWS, Mantid::Kernel::SpecialCoordinateSystem CoordinatesToUse, double radius)
Check if peaks overlap.
void init() override
Initialise the properties.
void integrate(typename DataObjects::MDEventWorkspace< MDE, nd >::sptr ws)
Integrate the peaks of the workspace using parameters saved in the algorithm class.
Mantid::API::IMDEventWorkspace_sptr inWS
Input MDEventWorkspace.
const std::string name() const override
Algorithm's name for identification.
std::map< std::string, std::string > validateInputs() override
Perform validation of ALL the input properties of the algorithm.
std::shared_ptr< IPeaksWorkspace > IPeaksWorkspace_sptr
shared pointer to Mantid::API::IPeaksWorkspace
std::shared_ptr< IAlgorithm > IAlgorithm_sptr
shared pointer to Mantid::API::IAlgorithm
std::complex< double > MANTID_API_DLL E1(std::complex< double > z)
Integral for Gamma.
std::shared_ptr< Workspace > Workspace_sptr
shared pointer to Mantid::API::Workspace
std::shared_ptr< IFunction > IFunction_sptr
shared pointer to the function base class
Definition IFunction.h:743
std::shared_ptr< PeaksWorkspace > PeaksWorkspace_sptr
Typedef for a shared pointer to a peaks workspace.
std::shared_ptr< Workspace2D > Workspace2D_sptr
shared pointer to Mantid::DataObjects::Workspace2D
std::array< double, PEAK_ELLIPSOID_DIMS > PeakEllipsoidExtent
long round(double x)
Custom rounding method for a double->long because none is portable in C++ (!)
Definition Utils.h:37
SpecialCoordinateSystem
Special coordinate systems for Q3D.
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.
double f_eval2(double x, void *params)
float coord_t
Typedef for the data type to use for coordinate axes in MD objects such as MDBox, MDEventWorkspace,...
Definition MDTypes.h:27
double signal_t
Typedef for the signal recorded in a MDBox, etc.
Definition MDTypes.h:36
@ Input
An input workspace.
Definition Property.h:53
@ Output
An output workspace.
Definition Property.h:54