Newer
Older
#include "MantidKernel/BoundedValidator.h"
#include <boost/shared_array.hpp>
#include <gsl/gsl_fft_complex.h>
namespace Mantid {
namespace Algorithms {
using Mantid::Kernel::Direction;
using Mantid::API::WorkspaceProperty;
using namespace API;
using namespace Kernel;
// Register the algorithm into the AlgorithmFactory
DECLARE_ALGORITHM(MaxEnt)
//----------------------------------------------------------------------------------------------
/** Constructor
*/
MaxEnt::MaxEnt() {}
//----------------------------------------------------------------------------------------------
/** Destructor
*/
MaxEnt::~MaxEnt() {}
//----------------------------------------------------------------------------------------------
/// Algorithm's name for identification. @see Algorithm::name
const std::string MaxEnt::name() const { return "MaxEnt"; }
/// Algorithm's version for identification. @see Algorithm::version
int MaxEnt::version() const { return 1; }
/// Algorithm's category for identification. @see Algorithm::category
const std::string MaxEnt::category() const { return "FFT"; }
/// Algorithm's summary for use in the GUI and help. @see Algorithm::summary
const std::string MaxEnt::summary() const {
return "Runs Maximum Entropy method on an input workspace";
}
//----------------------------------------------------------------------------------------------
/** Initialize the algorithm's properties.
*/
void MaxEnt::init() {
declareProperty(
new WorkspaceProperty<>("InputWorkspace", "", Direction::Input),
"An input workspace.");
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
auto mustBeNonNegative = boost::make_shared<BoundedValidator<double>>();
mustBeNonNegative->setLower(1E-12);
declareProperty(new PropertyWithValue<double>(
"Background", 0.4, mustBeNonNegative, Direction::Input),
"Default level above which the image is significant");
declareProperty(new PropertyWithValue<double>(
"ChiTarget", 100.0, mustBeNonNegative, Direction::Input),
"Target value of Chi-square");
declareProperty(new PropertyWithValue<double>(
"ChiEps", 0.001, mustBeNonNegative, Direction::Input),
"Required precision for Chi-square");
declareProperty(new PropertyWithValue<double>("DistancePenalty", 0.1,
mustBeNonNegative,
Direction::Input),
"Distance penalty applied to the current image");
declareProperty(new PropertyWithValue<double>(
"MaxAngle", 0.05, mustBeNonNegative, Direction::Input),
"Maximum degree of non-parallelism between S and C");
auto mustBePositive = boost::make_shared<BoundedValidator<size_t>>();
mustBePositive->setLower(0);
declareProperty(new PropertyWithValue<size_t>(
"MaxIterations", 20000, mustBePositive, Direction::Input),
"Maximum number of iterations");
declareProperty(new PropertyWithValue<size_t>("AlphaChopIterations", 500,
mustBePositive,
Direction::Input),
"Maximum number of iterations in alpha chop");
declareProperty(new WorkspaceProperty<>("EvolChi", "", Direction::Output),
"Output workspace containing the evolution of Chi-sq");
declareProperty(new WorkspaceProperty<>("EvolAngle", "", Direction::Output),
"Output workspace containing the evolution of "
"non-paralellism between S and C");
declareProperty(
new WorkspaceProperty<>("ReconstructedImage", "", Direction::Output),
"The output workspace containing the reconstructed image.");
new WorkspaceProperty<>("ReconstructedData", "", Direction::Output),
"The output workspace containing the reconstructed data.");
}
//----------------------------------------------------------------------------------------------
/** Validate the input properties.
*/
std::map<std::string, std::string> MaxEnt::validateInputs() {
std::map<std::string, std::string> result;
// X values in input workspace must be equally spaced
MatrixWorkspace_sptr inWS = getProperty("InputWorkspace");
const MantidVec &X = inWS->readX(0);
const double dx = X[1] - X[0];
for (size_t i = 1; i < X.size() - 2; i++) {
if (std::abs(dx - X[i + 1] + X[i]) / dx > 1e-7) {
result["InputWorkspace"] =
"X axis must be linear (all bins must have the same width)";
}
}
return result;
}
//----------------------------------------------------------------------------------------------
/** Execute the algorithm.
*/
void MaxEnt::exec() {
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
179
180
181
182
183
184
185
186
187
188
189
// Read input workspace
MatrixWorkspace_sptr inWS = getProperty("InputWorkspace");
// Background (default level, sky background, etc)
double background = getProperty("DefaultLevel");
// Chi target
double chiTarget = getProperty("ChiTarget");
// Required precision for Chi arget
double chiEps = getProperty("ChiEps");
// Maximum degree of non-parallelism between S and C
double angle = getProperty("MaxAngle");
// Distance penalty for current image
double distEps = getProperty("DistancePenalty");
// Maximum number of iterations
size_t niter = getProperty("MaxIterations");
// Maximum number of iterations in alpha chop
size_t alphaIter = getProperty("AlphaChopIterations");
// Number of spectra
size_t nspec = inWS->getNumberHistograms();
// Number of data points
size_t npoints = inWS->blocksize();
// Create output workspaces
MatrixWorkspace_sptr outImageWS = WorkspaceFactory::Instance().create(inWS);
MatrixWorkspace_sptr outDataWS = WorkspaceFactory::Instance().create(inWS);
// Create evol workspaces
MatrixWorkspace_sptr outEvolChi =
WorkspaceFactory::Instance().create(inWS, nspec, niter, niter);
MatrixWorkspace_sptr outEvolTest =
WorkspaceFactory::Instance().create(inWS, nspec, niter, niter);
// Start distribution (flat background)
std::vector<double> image(npoints, background * 1.01);
for (size_t s = 0; s < nspec; s++) {
// Progress
Progress progress(this, 0, 1, niter);
// Run maxent algorithm
for (size_t it = 0; it < niter; it++) {
// Check for canceling the algorithm
if (!(it % 1000)) {
interruption_point();
}
progress.report();
} // iterations
// Populate the output workspaces
// Populate workspaces recording the evolution of Chi and Test
} // Next spectrum
setProperty("EvolChi", outEvolChi);
setProperty("EvolAngle", outEvolTest);
setProperty("ReconstructedImage", outImageWS);
setProperty("ReconstructedData", outDataWS);
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
//----------------------------------------------------------------------------------------------
/**
* Transforms from solution-space to data-space
* @param input :: [input] An input vector in image space
* @return :: The input vector transformed to data space
*/
std::vector<double> MaxEnt::opus(const std::vector<double> &input) {
/* Performs backward Fourier transform */
size_t n = input.size();
/* Prepare the data */
boost::shared_array<double> result(new double[2 * n]);
for (size_t i = 0; i < n; i++) {
result[2 * i] = input[i]; // Real part
result[2 * i + 1] = 0.; // Imaginary part
}
/* Backward FT */
gsl_fft_complex_wavetable *wavetable = gsl_fft_complex_wavetable_alloc(n);
gsl_fft_complex_workspace *workspace = gsl_fft_complex_workspace_alloc(n);
gsl_fft_complex_inverse(result.get(), 1, n, wavetable, workspace);
gsl_fft_complex_wavetable_free(wavetable);
gsl_fft_complex_workspace_free(workspace);
std::vector<double> output(n);
for (size_t i = 0; i < n; i++) {
output[i] = result[2 * i];
}
return output;
}
//----------------------------------------------------------------------------------------------
/**
* Transforms from data-space to solution-space
* @param input :: [input] An input vector in data space
* @return :: The input vector converted to image space
*/
std::vector<double> MaxEnt::tropus(const std::vector<double> &input) {
/* Performs forward Fourier transform */
size_t n = input.size();
/* Prepare the data */
boost::shared_array<double> result(new double[n * 2]);
for (size_t i = 0; i < n; i++) {
result[2 * i] = input[i]; // even indexes filled with the real part
result[2 * i + 1] = 0.; // odd indexes filled with the imaginary part
}
/* Fourier transofrm */
gsl_fft_complex_wavetable *wavetable = gsl_fft_complex_wavetable_alloc(n);
gsl_fft_complex_workspace *workspace = gsl_fft_complex_workspace_alloc(n);
gsl_fft_complex_forward(result.get(), 1, n, wavetable, workspace);
gsl_fft_complex_wavetable_free(wavetable);
gsl_fft_complex_workspace_free(workspace);
/* Get the data */
std::vector<double> output(n);
for (size_t i = 0; i < n; i++) {
output[i] = result[2 * i];
}
return output;
}
} // namespace Algorithms
} // namespace Mantid