-
Notifications
You must be signed in to change notification settings - Fork 122
/
ParallaxCorrection.cpp
178 lines (164 loc) · 7.52 KB
/
ParallaxCorrection.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
// Mantid Repository : https://github.com/mantidproject/mantid
//
// Copyright © 2018 ISIS Rutherford Appleton Laboratory UKRI,
// NScD Oak Ridge National Laboratory, European Spallation Source,
// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
// SPDX - License - Identifier: GPL - 3.0 +
#include "MantidAlgorithms/ParallaxCorrection.h"
#include "MantidAPI/InstrumentValidator.h"
#include "MantidAPI/MatrixWorkspace.h"
#include "MantidAPI/Progress.h"
#include "MantidAPI/WorkspaceUnitValidator.h"
#include "MantidGeometry/Instrument.h"
#include "MantidGeometry/Instrument/ComponentInfo.h"
#include "MantidGeometry/Instrument/DetectorInfo.h"
#include "MantidHistogramData/HistogramMath.h"
#include "MantidKernel/ArrayLengthValidator.h"
#include "MantidKernel/ArrayProperty.h"
#include "MantidKernel/CompositeValidator.h"
#include <math.h>
#include <muParser.h>
namespace {
static const std::string PARALLAX_PARAMETER = "parallax";
static const std::string DIRECTION_PARAMETER = "direction";
/**
* @brief validateFormula Checks if the formula is evaluable
* @param parallax : formula
* @param direction : x or y
* @return empty string if valid, error message otherwise
*/
std::string validateFormula(const std::string ¶llax, const std::string &direction) {
if (direction != "x" && direction != "y") {
return "Direction must be x or y";
}
mu::Parser muParser;
double t = 0.;
muParser.DefineVar("t", &t);
muParser.SetExpr(parallax);
try {
muParser.Eval();
} catch (mu::Parser::exception_type &e) {
return e.GetMsg();
}
return "";
}
} // namespace
namespace Mantid {
namespace Algorithms {
// Register the algorithm into the AlgorithmFactory
DECLARE_ALGORITHM(ParallaxCorrection)
//----------------------------------------------------------------------------------------------
/// Algorithms name for identification. @see Algorithm::name
const std::string ParallaxCorrection::name() const { return "ParallaxCorrection"; }
/// Algorithm's version for identification. @see Algorithm::version
int ParallaxCorrection::version() const { return 1; }
/// Algorithm's category for identification. @see Algorithm::category
const std::string ParallaxCorrection::category() const { return "SANS"; }
/// Algorithm's summary for use in the GUI and help. @see Algorithm::summary
const std::string ParallaxCorrection::summary() const {
return "Performs parallax correction for tube based SANS instruments.";
}
//----------------------------------------------------------------------------------------------
/** Initialize the algorithm's properties.
*/
void ParallaxCorrection::init() {
auto validator = std::make_shared<Kernel::CompositeValidator>();
validator->add(std::make_unique<API::InstrumentValidator>());
validator->add(std::make_unique<API::WorkspaceUnitValidator>("Wavelength"));
auto lengthValidator = std::make_shared<Kernel::ArrayLengthValidator<std::string>>();
lengthValidator->setLengthMin(1);
declareProperty(std::make_unique<API::WorkspaceProperty<API::MatrixWorkspace>>("InputWorkspace", "",
Kernel::Direction::Input, validator),
"An input workspace.");
declareProperty(std::make_unique<Kernel::ArrayProperty<std::string>>("ComponentNames", lengthValidator),
"List of instrument components to perform the corrections for.");
declareProperty(
std::make_unique<API::WorkspaceProperty<API::MatrixWorkspace>>("OutputWorkspace", "", Kernel::Direction::Output),
"An output workspace.");
}
//----------------------------------------------------------------------------------------------
/**
* @brief ParallaxCorrection::performCorrection for the given bank
* @param outWS : the corrected workspace
* @param indices : the workspaces indices corresponding to the bank
* @param parallax : the correction formula for the bank
* @param direction : the tube direction in the bank
*/
void ParallaxCorrection::performCorrection(const API::MatrixWorkspace_sptr &outWS, const std::vector<size_t> &indices,
const std::string ¶llax, const std::string &direction) {
double t;
mu::Parser muParser;
muParser.DefineVar("t", &t);
muParser.SetExpr(parallax);
const auto &detectorInfo = outWS->detectorInfo();
// note that this is intenionally serial
for (const auto wsIndex : indices) {
const Kernel::V3D pos = detectorInfo.position(wsIndex);
if (direction == "y") {
t = std::fabs(std::atan2(pos.X(), pos.Z()));
} else {
t = std::fabs(std::atan2(pos.Y(), pos.Z()));
}
const double correction = muParser.Eval();
if (correction > 0.) {
auto &spectrum = outWS->mutableY(wsIndex);
auto &errors = outWS->mutableE(wsIndex);
spectrum /= correction;
errors /= correction;
} else {
g_log.warning() << "Correction is <=0 for workspace index " << wsIndex << ". Skipping the correction.\n";
}
}
}
//----------------------------------------------------------------------------------------------
/** Execute the algorithm.
*/
void ParallaxCorrection::exec() {
API::MatrixWorkspace_const_sptr inputWorkspace = getProperty("InputWorkspace");
API::MatrixWorkspace_sptr outputWorkspace = getProperty("OutputWorkspace");
if (inputWorkspace != outputWorkspace) {
outputWorkspace = inputWorkspace->clone();
}
const std::vector<std::string> componentNames = getProperty("ComponentNames");
const auto &instrument = inputWorkspace->getInstrument();
const auto &detectorInfo = outputWorkspace->detectorInfo();
const auto &allDetIDs = detectorInfo.detectorIDs();
const auto &componentInfo = outputWorkspace->componentInfo();
auto progress = std::make_unique<API::Progress>(this, 0., 1., componentNames.size());
for (const auto &componentName : componentNames) {
progress->report("Performing parallax correction for component " + componentName);
const auto component = instrument->getComponentByName(componentName);
if (!component) {
g_log.error() << "No component defined with name " << componentName << "\n";
continue;
}
if (!component->hasParameter(PARALLAX_PARAMETER) || !component->hasParameter(DIRECTION_PARAMETER)) {
g_log.error() << "No parallax correction defined in IPF for component " << componentName << "\n";
continue;
}
const std::string parallax = component->getStringParameter(PARALLAX_PARAMETER)[0];
const std::string direction = component->getStringParameter(DIRECTION_PARAMETER)[0];
const auto valid = validateFormula(parallax, direction);
if (!valid.empty()) {
g_log.error() << "Unable to parse the parallax formula and direction "
"for component "
<< componentName << ". Reason: " << valid << "\n";
continue;
}
const auto componentIndex = componentInfo.indexOfAny(componentName);
const auto &detectorIndices = componentInfo.detectorsInSubtree(componentIndex);
if (detectorIndices.empty()) {
g_log.error() << "No detectors found in component " << componentName << "\n";
continue;
}
std::vector<detid_t> detIDs;
detIDs.reserve(detectorIndices.size());
std::transform(detectorIndices.cbegin(), detectorIndices.cend(), std::back_inserter(detIDs),
[&allDetIDs](size_t i) { return allDetIDs[i]; });
const auto indices = outputWorkspace->getIndicesFromDetectorIDs(detIDs);
performCorrection(outputWorkspace, indices, parallax, direction);
}
setProperty("OutputWorkspace", outputWorkspace);
}
} // namespace Algorithms
} // namespace Mantid