-
Notifications
You must be signed in to change notification settings - Fork 122
/
MaxEnt.cpp
856 lines (718 loc) · 30.3 KB
/
MaxEnt.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
179
180
181
182
183
184
185
186
187
188
189
190
191
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
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
#include "MantidAlgorithms/MaxEnt.h"
#include "MantidAPI/EqualBinSizesValidator.h"
#include "MantidAPI/MatrixWorkspace.h"
#include "MantidAPI/TextAxis.h"
#include "MantidAPI/WorkspaceFactory.h"
#include "MantidAlgorithms/MaxEnt/MaxentEntropyNegativeValues.h"
#include "MantidAlgorithms/MaxEnt/MaxentEntropyPositiveValues.h"
#include "MantidAlgorithms/MaxEnt/MaxentSpaceComplex.h"
#include "MantidAlgorithms/MaxEnt/MaxentSpaceReal.h"
#include "MantidAlgorithms/MaxEnt/MaxentTransformFourier.h"
#include "MantidHistogramData/LinearGenerator.h"
#include "MantidKernel/BoundedValidator.h"
#include "MantidKernel/ListValidator.h"
#include "MantidKernel/UnitFactory.h"
#include "MantidKernel/VectorHelper.h"
#include <gsl/gsl_linalg.h>
#include <algorithm>
#include <numeric>
namespace Mantid {
namespace Algorithms {
using Mantid::HistogramData::LinearGenerator;
using Mantid::HistogramData::Points;
using namespace API;
using namespace Kernel;
// Register the algorithm into the AlgorithmFactory
DECLARE_ALGORITHM(MaxEnt)
namespace {
// Maps defining the inverse caption and label for the reconstructed image
// Example:
// The input workspaces (X axis) is in (Time, s)
// The output image should be in (Frequency, Hz)
// Defines the new caption
std::map<std::string, std::string> inverseCaption = {{"Time", "Frequency"},
{"Frequency", "Time"},
{"d-Spacing", "q"},
{"q", "d-Spacing"}};
// Defines the new label
std::map<std::string, std::string> inverseLabel = {{"s", "Hz"},
{"microsecond", "MHz"},
{"Hz", "s"},
{"MHz", "microsecond"},
{"Angstrom", "Angstrom^-1"},
{"Angstrom^-1", "Angstrom"}};
// A threshold for small singular values
const double THRESHOLD = 1E-6;
/** removes zeros from converged results
* @param ws :: [input] The input workspace with zeros
* @param maxIt :: [input] The max number of iteration this alg used
* @param yLabel :: [input] y-label to use for returned ws
* @return : ws cut down in lenght to maxIt
*/
MatrixWorkspace_sptr removeZeros(const MatrixWorkspace_sptr &ws,
const size_t maxIt,
const ::std::string yLabel) {
ws->setYUnitLabel(yLabel);
try {
ws->getAxis(0)->unit() =
UnitFactory::Instance().create("Number of Iterations");
} catch (Exception::NotFoundError &) {
ws->getAxis(0)->unit() = UnitFactory::Instance().create("Label");
Unit_sptr unit = ws->getAxis(0)->unit();
boost::shared_ptr<Units::Label> label =
boost::dynamic_pointer_cast<Units::Label>(unit);
label->setLabel("Number of Iterations", "");
}
if (maxIt == ws->readY(0).size()) {
return ws;
}
const size_t nspec = ws->getNumberHistograms();
MatrixWorkspace_sptr outWS =
WorkspaceFactory::Instance().create(ws, nspec, maxIt, maxIt);
for (size_t spec = 0; spec < nspec; spec++) {
outWS->setPoints(spec, Points(maxIt, LinearGenerator(0.0, 1.0)));
auto Data = ws->readY(spec);
outWS->setCounts(spec,
std::vector<double>(Data.begin(), Data.begin() + maxIt));
}
return outWS;
}
}
//----------------------------------------------------------------------------------------------
/// 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 "Arithmetic\\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 every spectrum of an input workspace. "
"It currently works for the case where data and image are related by a"
" 1D Fourier transform.";
}
//----------------------------------------------------------------------------------------------
/** Initialize the algorithm's properties.
*/
void MaxEnt::init() {
// X values in input workspace must be (almost) equally spaced
const double warningLevel = 0.01;
const double errorLevel = 0.5;
declareProperty(
make_unique<WorkspaceProperty<>>(
"InputWorkspace", "", Direction::Input,
boost::make_shared<EqualBinSizesValidator>(errorLevel, warningLevel)),
"An input workspace.");
declareProperty("ComplexData", false,
"If true, the input data is assumed to be complex and the "
"input workspace is expected to have an even number of "
"histograms (2N). Spectrum numbers S and S+N are assumed to "
"be the real and imaginary part of the complex signal "
"respectively.");
declareProperty("ComplexImage", true,
"If true, the algorithm will use complex images for the "
"calculations. This is the recommended option when there is "
"no prior knowledge about the image. If the image is known "
"to be real, this option can be set to false and the "
"algorithm will only consider the real part for "
"calculations.");
declareProperty("PositiveImage", false,
"If true, the reconstructed image is only allowed to take "
"positive values. It can take negative values otherwise. "
"This option defines the entropy formula that will be used "
"for the calculations (see next section for more details).");
declareProperty("AutoShift", false,
"Automatically calculate and apply phase shift. Zero on the "
"X axis is assumed to be in the first bin. If it is not, "
"setting this property will automatically correct for this.");
auto mustBePositive = boost::make_shared<BoundedValidator<size_t>>();
mustBePositive->setLower(0);
declareProperty(make_unique<PropertyWithValue<size_t>>(
"ResolutionFactor", 1, mustBePositive, Direction::Input),
"An integer number indicating the factor by which the number "
"of points will be increased in the image and reconstructed "
"data");
auto mustBeNonNegative = boost::make_shared<BoundedValidator<double>>();
mustBeNonNegative->setLower(1E-12);
declareProperty(
make_unique<PropertyWithValue<double>>("A", 0.4, mustBeNonNegative,
Direction::Input),
"A maximum entropy constant. This algorithm was first developed for the "
"ISIS muon group where the default 0.4 was found to give good "
"reconstructions. "
"In general the user will need to experiment with this value. Choosing a "
"small value may lead to unphysical spiky reconstructions and choosing "
"an increasingly large "
"value the reconstruction will start to resamble that of a direct "
"fourier "
"transform reconstruction. However, where the data contain a "
"zero Fourier data point with a small error the "
"reconstruction will be insensitive to the choice "
"of this property (and increasing so the more well determined "
"this data point is).");
declareProperty(
make_unique<PropertyWithValue<double>>(
"ChiTargetOverN", 1.0, mustBeNonNegative, Direction::Input),
"Target value of Chi-square divided by the number of data points (N)");
declareProperty(make_unique<PropertyWithValue<double>>(
"ChiEps", 0.001, mustBeNonNegative, Direction::Input),
"Required precision for Chi-square");
declareProperty(
make_unique<PropertyWithValue<double>>(
"DistancePenalty", 0.1, mustBeNonNegative, Direction::Input),
"Distance penalty applied to the current image at each iteration.");
declareProperty(
make_unique<PropertyWithValue<double>>(
"MaxAngle", 0.001, mustBeNonNegative, Direction::Input),
"Maximum degree of non-parallelism between S (the entropy) and C "
"(chi-squared). These needs to be parallel. Chosing a smaller "
"shouldn't change the output. However, if you find this is the "
"case please let the Mantid team know since this indicates that "
"the default value of this proporty may need changing or "
"other changes to this implementation are required.");
mustBePositive = boost::make_shared<BoundedValidator<size_t>>();
mustBePositive->setLower(1);
declareProperty(make_unique<PropertyWithValue<size_t>>(
"MaxIterations", 20000, mustBePositive, Direction::Input),
"Maximum number of iterations.");
declareProperty(make_unique<PropertyWithValue<size_t>>("AlphaChopIterations",
500, mustBePositive,
Direction::Input),
"Maximum number of iterations in alpha chop.");
declareProperty(
make_unique<WorkspaceProperty<>>("EvolChi", "", Direction::Output),
"Output workspace containing the evolution of Chi-sq.");
declareProperty(
make_unique<WorkspaceProperty<>>("EvolAngle", "", Direction::Output),
"Output workspace containing the evolution of "
"non-paralellism between S and C.");
declareProperty(make_unique<WorkspaceProperty<>>("ReconstructedImage", "",
Direction::Output),
"The output workspace containing the reconstructed image.");
declareProperty(make_unique<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;
MatrixWorkspace_sptr inWS = getProperty("InputWorkspace");
if (inWS) {
// If the input signal is complex, we expect an even number of histograms
// in the input workspace
size_t nhistograms = inWS->getNumberHistograms();
bool complex = getProperty("ComplexData");
if (complex && (nhistograms % 2))
result["InputWorkspace"] = "The number of histograms in the input "
"workspace must be even for complex data";
}
return result;
}
//----------------------------------------------------------------------------------------------
/** Execute the algorithm.
*/
void MaxEnt::exec() {
// MaxEnt parameters
// Complex data?
const bool complexData = getProperty("ComplexData");
// Complex image?
const bool complexImage = getProperty("ComplexImage");
// Image must be positive?
const bool positiveImage = getProperty("PositiveImage");
// Autoshift
const bool autoShift = getProperty("AutoShift");
// Increase the number of points in the image by this factor
const size_t resolutionFactor = getProperty("ResolutionFactor");
// Background (default level, sky background, etc)
const double background = getProperty("A");
// Chi target
const double ChiTargetOverN = getProperty("ChiTargetOverN");
// Required precision for Chi arget
const double chiEps = getProperty("ChiEps");
// Maximum degree of non-parallelism between S and C
const double angle = getProperty("MaxAngle");
// Distance penalty for current image
const double distEps = getProperty("DistancePenalty");
// Maximum number of iterations
const size_t niter = getProperty("MaxIterations");
// Maximum number of iterations in alpha chop
const size_t alphaIter = getProperty("AlphaChopIterations");
// Number of spectra and datapoints
// Read input workspace
MatrixWorkspace_const_sptr inWS = getProperty("InputWorkspace");
// Number of spectra
size_t nspec = inWS->getNumberHistograms();
// Number of data points - assumed to be constant between spectra or
// this will throw an exception
size_t npoints = inWS->blocksize() * resolutionFactor;
// Number of X bins
const size_t npointsX = inWS->isHistogramData() ? npoints + 1 : npoints;
// For now have the requirement that data must have non-zero
// (and positive!) errors
for (size_t s = 0; s < nspec; s++) {
auto errors = inWS->e(s).rawData();
size_t npoints = errors.size();
for (size_t i = 0; i < npoints; i++) {
if (errors[i] <= 0.0) {
throw std::invalid_argument("Input data must have non-zero errors.");
}
}
}
// Is our data space real or complex?
MaxentSpace_sptr dataSpace;
if (complexData) {
dataSpace = std::make_shared<MaxentSpaceComplex>();
} else {
dataSpace = std::make_shared<MaxentSpaceReal>();
}
// Is our image space real or complex?
MaxentSpace_sptr imageSpace;
if (complexImage) {
imageSpace = std::make_shared<MaxentSpaceComplex>();
} else {
imageSpace = std::make_shared<MaxentSpaceReal>();
}
// The type of transform. Currently a 1D Fourier Transform
MaxentTransform_sptr transform =
std::make_shared<MaxentTransformFourier>(dataSpace, imageSpace);
// The type of entropy we are going to use (depends on the type of image,
// positive only, or positive and/or negative)
MaxentEntropy_sptr entropy;
if (positiveImage) {
entropy = std::make_shared<MaxentEntropyPositiveValues>();
} else {
entropy = std::make_shared<MaxentEntropyNegativeValues>();
}
// Entropy and transform is all we need to set up a calculator
MaxentCalculator maxentCalculator = MaxentCalculator(entropy, transform);
// Output workspaces
MatrixWorkspace_sptr outImageWS;
MatrixWorkspace_sptr outDataWS;
MatrixWorkspace_sptr outEvolChi;
MatrixWorkspace_sptr outEvolTest;
nspec = complexData ? nspec / 2 : nspec;
outImageWS =
WorkspaceFactory::Instance().create(inWS, 2 * nspec, npoints, npoints);
for (size_t i = 0; i < outImageWS->getNumberHistograms(); ++i)
outImageWS->getSpectrum(i).setDetectorID(static_cast<detid_t>(i + 1));
outDataWS =
WorkspaceFactory::Instance().create(inWS, 2 * nspec, npointsX, npoints);
for (size_t i = 0; i < outDataWS->getNumberHistograms(); ++i)
outDataWS->getSpectrum(i).setDetectorID(static_cast<detid_t>(i + 1));
outEvolChi = WorkspaceFactory::Instance().create(inWS, nspec, niter, niter);
outEvolTest = WorkspaceFactory::Instance().create(inWS, nspec, niter, niter);
npoints = complexImage ? npoints * 2 : npoints;
size_t maxIt = 0; // used to determine max iterations used by alg
outEvolChi->setPoints(0, Points(niter, LinearGenerator(0.0, 1.0)));
for (size_t s = 0; s < nspec; s++) {
// Start distribution (flat background)
std::vector<double> image(npoints, background);
std::vector<double> data;
std::vector<double> errors;
if (complexData) {
data = toComplex(inWS, s, false); // false -> data
errors = toComplex(inWS, s, true); // true -> errors
} else {
data = inWS->y(s).rawData();
errors = inWS->e(s).rawData();
}
// To record the algorithm's progress
std::vector<double> evolChi(niter, 0.);
std::vector<double> evolTest(niter, 0.);
// Progress
Progress progress(this, 0.0, 1.0, niter);
// Run maxent algorithm
for (size_t it = 0; it < niter; it++) {
// Iterates one step towards the solution. This means calculating
// quadratic coefficients, search directions, angle and chi-sq
maxentCalculator.iterate(data, errors, image, background);
// Calculate delta to construct new image (SB eq. 25)
double currChisq = maxentCalculator.getChisq();
auto coeffs = maxentCalculator.getQuadraticCoefficients();
auto delta = move(coeffs, ChiTargetOverN / currChisq, chiEps, alphaIter);
// Apply distance penalty (SB eq. 33)
delta = applyDistancePenalty(delta, coeffs, image, background, distEps);
// Update image
auto dirs = maxentCalculator.getSearchDirections();
image = updateImage(image, delta, dirs);
// Record the evolution of Chi-square and angle(S,C)
double currAngle = maxentCalculator.getAngle();
evolChi[it] = currChisq;
evolTest[it] = currAngle;
if (it > maxIt) {
maxIt = it;
}
// Stop condition, solution found
if ((std::abs(currChisq / ChiTargetOverN - 1.) < chiEps) &&
(currAngle < angle)) {
g_log.information() << "Stopped after " << it << " iterations"
<< std::endl;
break;
}
// Check for canceling the algorithm
if (!(it % 1000)) {
interruption_point();
}
progress.report();
} // iterations
// Get calculated data
auto solData = maxentCalculator.getReconstructedData();
auto solImage = maxentCalculator.getImage();
// Populate the output workspaces
populateDataWS(inWS, s, nspec, solData, complexData, outDataWS);
populateImageWS(inWS, s, nspec, solImage, complexImage, outImageWS,
autoShift);
// Populate workspaces recording the evolution of Chi and Test
// X values
outEvolChi->setSharedX(s, outEvolChi->sharedX(0));
outEvolTest->setSharedX(s, outEvolChi->sharedX(0));
// Y values
outEvolChi->setCounts(s, std::move(evolChi));
outEvolTest->setCounts(s, std::move(evolTest));
// No errors
} // Next spectrum
// add 1 to maxIt to account for starting at 0
maxIt++;
setProperty("EvolChi", removeZeros(outEvolChi, maxIt, "Chi squared"));
setProperty("EvolAngle", removeZeros(outEvolTest, maxIt, "Maximum Angle"));
setProperty("ReconstructedImage", outImageWS);
setProperty("ReconstructedData", outDataWS);
}
//----------------------------------------------------------------------------------------------
/** Returns a given spectrum as a complex number
* @param inWS :: [input] The input workspace containing all the spectra
* @param spec :: [input] The spectrum of interest
* @param errors :: [input] If true, returns the errors, otherwise returns the
* counts
* @return : Spectrum 'spec' as a complex vector
*/
std::vector<double> MaxEnt::toComplex(API::MatrixWorkspace_const_sptr &inWS,
size_t spec, bool errors) {
const size_t numBins = inWS->y(0).size();
std::vector<double> result(numBins * 2);
if (inWS->getNumberHistograms() % 2)
throw std::invalid_argument(
"Cannot convert input workspace to complex data");
size_t nspec = inWS->getNumberHistograms() / 2;
if (!errors) {
for (size_t i = 0; i < numBins; i++) {
result[2 * i] = inWS->y(spec)[i];
result[2 * i + 1] = inWS->y(spec + nspec)[i];
}
} else {
for (size_t i = 0; i < numBins; i++) {
result[2 * i] = inWS->e(spec)[i];
result[2 * i + 1] = inWS->e(spec + nspec)[i];
}
}
return result;
}
/** Bisection method to move delta one step closer towards the solution
* @param coeffs :: [input] The current quadratic coefficients
* @param ChiTargetOverN :: [input] The requested Chi target over N
* (data points)
* @param chiEps :: [input] Precision required for Chi target
* @param alphaIter :: [input] Maximum number of iterations in the bisection
* method (alpha chop)
* @return : The increment length to be added to the current image
*/
std::vector<double> MaxEnt::move(const QuadraticCoefficients &coeffs,
double ChiTargetOverN, double chiEps,
size_t alphaIter) {
double aMin = 0.; // Minimum alpha
double aMax = 1.; // Maximum alpha
// Dimension, number of search directions
size_t dim = coeffs.c2.size().first;
std::vector<double> deltaMin(dim, 0); // delta at alpha min
std::vector<double> deltaMax(dim, 0); // delta at alpha max
double chiMin = calculateChi(coeffs, aMin, deltaMin); // Chi at alpha min
double chiMax = calculateChi(coeffs, aMax, deltaMax); // Chi at alpha max
double dchiMin = chiMin - ChiTargetOverN; // max - target
double dchiMax = chiMax - ChiTargetOverN; // min - target
if (dchiMin * dchiMax > 0) {
// ChiTargetOverN could be outside the range [chiMin, chiMax]
if (fabs(dchiMin) < fabs(dchiMax)) {
return deltaMin;
} else {
return deltaMax;
}
// throw std::runtime_error("Error in alpha chop\n");
}
// Initial values of eps and iter to start while loop
double eps = 2. * chiEps;
size_t iter = 0;
// Bisection method
std::vector<double> delta(dim, 0); // delta at current alpha
while ((fabs(eps / ChiTargetOverN) > chiEps) && (iter < alphaIter)) {
double aMid = 0.5 * (aMin + aMax);
double chiMid = calculateChi(coeffs, aMid, delta);
eps = chiMid - ChiTargetOverN;
if (dchiMin * eps > 0) {
aMin = aMid;
dchiMin = eps;
}
if (dchiMax * eps > 0) {
aMax = aMid;
dchiMax = eps;
}
iter++;
}
// Check if move was successful
if ((fabs(eps / ChiTargetOverN) > chiEps) || (iter > alphaIter)) {
throw std::runtime_error("Error encountered when calculating solution "
"image. No convergence in alpha chop.\n");
}
return delta;
}
/** Calculates Chi given the quadratic coefficients and an alpha value by
* solving the matrix equation A*b = B
* @param coeffs :: [input] The quadratic coefficients
* @param a :: [input] The alpha value
* @param b :: [output] The solution
* @return :: The calculated chi-square
*/
double MaxEnt::calculateChi(const QuadraticCoefficients &coeffs, double a,
std::vector<double> &b) {
size_t dim = coeffs.c2.size().first;
double ax = a;
double bx = 1 - ax;
Kernel::DblMatrix A(dim, dim);
Kernel::DblMatrix B(dim, 1);
// Construct the matrix A and vector B such that Ax=B
for (size_t k = 0; k < dim; k++) {
for (size_t l = 0; l < dim; l++) {
A[k][l] = bx * coeffs.c2[k][l] - ax * coeffs.s2[k][l];
}
B[k][0] = -bx * coeffs.c1[k][0] + ax * coeffs.s1[k][0];
}
// Alternatives I have tried:
// Gauss-Jordan
// LU
// SVD seems to work better
// Solve using SVD
b = solveSVD(A, B);
// Now compute Chi
double ww = 0.;
for (size_t k = 0; k < dim; k++) {
double z = 0.;
for (size_t l = 0; l < dim; l++) {
z += coeffs.c2[k][l] * b[l];
}
ww += b[k] * (coeffs.c1[k][0] + 0.5 * z);
}
// Return chi
return ww + 1.;
}
/** Solves A*x = B using SVD
* @param A :: [input] The matrix A
* @param B :: [input] The vector B
* @return :: The solution x
*/
std::vector<double> MaxEnt::solveSVD(DblMatrix &A, const DblMatrix &B) {
size_t dim = A.size().first;
auto a = gsl_matrix_view_array(A[0], dim, dim);
auto b = gsl_vector_const_view_array(B[0], dim);
std::vector<double> vVec(dim * dim), sVec(dim), wVec(dim), delta(dim);
auto v = gsl_matrix_view_array(vVec.data(), dim, dim);
auto s = gsl_vector_view_array(sVec.data(), dim);
auto w = gsl_vector_view_array(wVec.data(), dim);
auto x = gsl_vector_view_array(delta.data(), dim);
// Singular value decomposition
gsl_linalg_SV_decomp(&a.matrix, &v.matrix, &s.vector, &w.vector);
// A could be singular or ill-conditioned. We can use SVD to obtain a least
// squares solution by setting the small (compared to the maximum) singular
// values to zero
// Find largest sing value
double max = *std::max_element(sVec.begin(), sVec.end());
// Apply a threshold to small singular values
double threshold = THRESHOLD * max;
std::transform(sVec.begin(), sVec.end(), sVec.begin(),
[&threshold](double el) { return el > threshold ? el : 0.0; });
// Solve A*x = B
gsl_linalg_SV_solve(&a.matrix, &v.matrix, &s.vector, &b.vector, &x.vector);
return delta;
}
/** Applies a distance penalty
* @param delta :: [input] The current increment
* @param coeffs :: [input] The quadratic coefficients
* @param image :: [input] The current image
* @param background :: [input] The background
* @param distEps :: [input] The distance constraint
* @return :: The new increment
*/
std::vector<double> MaxEnt::applyDistancePenalty(
const std::vector<double> &delta, const QuadraticCoefficients &coeffs,
const std::vector<double> &image, double background, double distEps) {
double sum = 0.;
for (double point : image)
sum += fabs(point);
size_t dim = coeffs.s2.size().first;
double dist = 0.;
for (size_t k = 0; k < dim; k++) {
double sum = 0.0;
for (size_t l = 0; l < dim; l++)
sum -= coeffs.s2[k][l] * delta[l];
dist += delta[k] * sum;
}
if (dist > distEps * sum / background) {
auto newDelta = delta;
for (size_t k = 0; k < delta.size(); k++) {
newDelta[k] *= sqrt(distEps * sum / dist / background);
}
return newDelta;
}
return delta;
}
/**
* Updates the image according to an increment delta
* @param image : [input] The current image as a vector (can be real or complex)
* @param delta : [input] The increment delta as a vector (can be real or
* complex)
* @param dirs : [input] The search directions
* @return : The new image
*/
std::vector<double>
MaxEnt::updateImage(const std::vector<double> &image,
const std::vector<double> &delta,
const std::vector<std::vector<double>> dirs) {
if (image.empty() || dirs.empty() || (delta.size() != dirs.size())) {
throw std::runtime_error("Cannot calculate new image");
}
std::vector<double> newImage = image;
// Calculate the new image
for (size_t i = 0; i < image.size(); i++) {
for (size_t k = 0; k < delta.size(); k++) {
newImage[i] += delta[k] * dirs[k][i];
}
}
return newImage;
}
/** Populates the output workspaces
* @param inWS :: [input] The input workspace
* @param spec :: [input] The current spectrum being analyzed
* @param nspec :: [input] The total number of histograms in the input workspace
* @param result :: [input] The image to be written in the output workspace (can
* be real or complex vector)
* @param complex :: [input] True if the result is a complex vector, false
* otherwise
* @param outWS :: [input] The output workspace to populate
* @param autoShift :: [input] Whether or not to correct the phase shift
*/
void MaxEnt::populateImageWS(MatrixWorkspace_const_sptr &inWS, size_t spec,
size_t nspec, const std::vector<double> &result,
bool complex, MatrixWorkspace_sptr &outWS,
bool autoShift) {
if (complex && result.size() % 2)
throw std::invalid_argument("Cannot write results to output workspaces");
int npoints = complex ? static_cast<int>(result.size() / 2)
: static_cast<int>(result.size());
MantidVec X(npoints);
MantidVec YR(npoints);
MantidVec YI(npoints);
MantidVec E(npoints, 0.);
auto dataPoints = inWS->points(spec);
double x0 = dataPoints[0];
double dx = dataPoints[1] - x0;
double delta = 1. / dx / npoints;
const int isOdd = (inWS->y(0).size() % 2) ? 1 : 0;
double shift = x0 * 2. * M_PI;
if (!autoShift)
shift = 0.;
// X values
for (int i = 0; i < npoints; i++) {
X[i] = delta * (-npoints / 2 + i);
}
// Y values
if (complex) {
for (int i = 0; i < npoints; i++) {
int j = (npoints / 2 + i + isOdd) % npoints;
double xShift = X[i] * shift;
double c = cos(xShift);
double s = sin(xShift);
YR[i] = result[2 * j] * c - result[2 * j + 1] * s;
YI[i] = result[2 * j] * s + result[2 * j + 1] * c;
YR[i] *= dx;
YI[i] *= dx;
}
} else {
for (int i = 0; i < npoints; i++) {
int j = (npoints / 2 + i + isOdd) % npoints;
double xShift = X[i] * shift;
double c = cos(xShift);
double s = sin(xShift);
YR[i] = result[j] * c;
YI[i] = result[j] * s;
YR[i] *= dx;
YI[i] *= dx;
}
}
// X caption & label
auto inputUnit = inWS->getAxis(0)->unit();
if (inputUnit) {
boost::shared_ptr<Kernel::Units::Label> lblUnit =
boost::dynamic_pointer_cast<Kernel::Units::Label>(
UnitFactory::Instance().create("Label"));
if (lblUnit) {
lblUnit->setLabel(
inverseCaption[inWS->getAxis(0)->unit()->caption()],
inverseLabel[inWS->getAxis(0)->unit()->label().ascii()]);
outWS->getAxis(0)->unit() = lblUnit;
}
}
outWS->mutableX(spec) = std::move(X);
outWS->mutableY(spec) = std::move(YR);
outWS->mutableE(spec) = std::move(E);
outWS->setSharedX(nspec + spec, outWS->sharedX(spec));
outWS->mutableY(nspec + spec) = std::move(YI);
outWS->setSharedE(nspec + spec, outWS->sharedE(spec));
}
/** Populates the output workspaces
* @param inWS :: [input] The input workspace
* @param spec :: [input] The current spectrum being analyzed
* @param nspec :: [input] The total number of histograms in the input workspace
* @param result :: [input] The reconstructed data to be written in the output
* workspace (can be a real or complex vector)
* @param complex :: [input] True if result is a complex vector, false otherwise
* @param outWS :: [input] The output workspace to populate
*/
void MaxEnt::populateDataWS(MatrixWorkspace_const_sptr &inWS, size_t spec,
size_t nspec, const std::vector<double> &result,
bool complex, MatrixWorkspace_sptr &outWS) {
if (complex && result.size() % 2)
throw std::invalid_argument("Cannot write results to output workspaces");
int npoints = complex ? static_cast<int>(result.size() / 2)
: static_cast<int>(result.size());
int npointsX = inWS->isHistogramData() ? npoints + 1 : npoints;
MantidVec X(npointsX);
MantidVec YR(npoints);
MantidVec YI(npoints);
MantidVec E(npoints, 0.);
double x0 = inWS->x(spec)[0];
double dx = inWS->x(spec)[1] - x0;
// X values
for (int i = 0; i < npointsX; i++) {
X[i] = x0 + i * dx;
}
// Y values
if (complex) {
for (int i = 0; i < npoints; i++) {
YR[i] = result[2 * i];
YI[i] = result[2 * i + 1];
}
} else {
for (int i = 0; i < npoints; i++) {
YR[i] = result[i];
YI[i] = 0.;
}
}
outWS->mutableX(spec) = std::move(X);
outWS->mutableY(spec) = std::move(YR);
outWS->mutableE(spec) = std::move(E);
outWS->setSharedX(nspec + spec, outWS->sharedX(spec));
outWS->mutableY(nspec + spec) = std::move(YI);
outWS->setSharedE(nspec + spec, outWS->sharedE(spec));
}
} // namespace Algorithms
} // namespace Mantid