-
Notifications
You must be signed in to change notification settings - Fork 122
/
SplineBackground.cpp
179 lines (148 loc) · 5.05 KB
/
SplineBackground.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
//----------------------------------------------------------------------
// Includes
//----------------------------------------------------------------------
#include "MantidCurveFitting/Algorithms/SplineBackground.h"
#include "MantidKernel/cow_ptr.h"
#include "MantidKernel/BoundedValidator.h"
#include "MantidAPI/MatrixWorkspace.h"
#include "MantidAPI/WorkspaceFactory.h"
#include <gsl/gsl_bspline.h>
#include <gsl/gsl_multifit.h>
#include <gsl/gsl_statistics.h>
namespace Mantid {
namespace CurveFitting {
namespace Algorithms {
DECLARE_ALGORITHM(SplineBackground)
using namespace Kernel;
using namespace API;
/// Initialisation method
void SplineBackground::init() {
declareProperty(make_unique<WorkspaceProperty<API::MatrixWorkspace>>(
"InputWorkspace", "", Direction::Input),
"The name of the input workspace.");
declareProperty(make_unique<WorkspaceProperty<API::MatrixWorkspace>>(
"OutputWorkspace", "", Direction::Output),
"The name to use for the output workspace.");
auto mustBePositive = boost::make_shared<BoundedValidator<int>>();
mustBePositive->setLower(0);
declareProperty("WorkspaceIndex", 0, mustBePositive,
"The index of the spectrum for fitting.");
declareProperty("NCoeff", 10, "The number of b-spline coefficients.");
}
/** Executes the algorithm
*
*/
void SplineBackground::exec() {
API::MatrixWorkspace_sptr inWS = getProperty("InputWorkspace");
int spec = getProperty("WorkspaceIndex");
if (spec > static_cast<int>(inWS->getNumberHistograms()))
throw std::out_of_range("WorkspaceIndex is out of range.");
const MantidVec &X = inWS->readX(spec);
const MantidVec &Y = inWS->readY(spec);
const MantidVec &E = inWS->readE(spec);
const bool isHistogram = inWS->isHistogramData();
const int ncoeffs = getProperty("NCoeff");
const int k = 4; // order of the spline + 1 (cubic)
const int nbreak = ncoeffs - (k - 2);
if (nbreak <= 0)
throw std::out_of_range("Too low NCoeff");
gsl_bspline_workspace *bw;
gsl_vector *B;
gsl_vector *c, *w, *x, *y;
gsl_matrix *Z, *cov;
gsl_multifit_linear_workspace *mw;
double chisq;
int n = static_cast<int>(Y.size());
bool isMasked = inWS->hasMaskedBins(spec);
std::vector<int> masked(Y.size());
if (isMasked) {
auto maskedBins = inWS->maskedBins(spec);
for (auto &maskedBin : maskedBins)
masked[maskedBin.first] = 1;
n -= static_cast<int>(inWS->maskedBins(spec).size());
}
if (n < ncoeffs) {
g_log.error("Too many basis functions (NCoeff)");
throw std::out_of_range("Too many basis functions (NCoeff)");
}
/* allocate a cubic bspline workspace (k = 4) */
bw = gsl_bspline_alloc(k, nbreak);
B = gsl_vector_alloc(ncoeffs);
x = gsl_vector_alloc(n);
y = gsl_vector_alloc(n);
Z = gsl_matrix_alloc(n, ncoeffs);
c = gsl_vector_alloc(ncoeffs);
w = gsl_vector_alloc(n);
cov = gsl_matrix_alloc(ncoeffs, ncoeffs);
mw = gsl_multifit_linear_alloc(n, ncoeffs);
/* this is the data to be fitted */
int j = 0;
for (MantidVec::size_type i = 0; i < Y.size(); ++i) {
if (isMasked && masked[i])
continue;
gsl_vector_set(x, j,
(isHistogram ? (0.5 * (X[i] + X[i + 1]))
: X[i])); // Middle of the bins, if a histogram
gsl_vector_set(y, j, Y[i]);
gsl_vector_set(w, j, E[i] > 0. ? 1. / (E[i] * E[i]) : 0.);
++j;
}
if (n != j) {
gsl_bspline_free(bw);
gsl_vector_free(B);
gsl_vector_free(x);
gsl_vector_free(y);
gsl_matrix_free(Z);
gsl_vector_free(c);
gsl_vector_free(w);
gsl_matrix_free(cov);
gsl_multifit_linear_free(mw);
throw std::runtime_error("Assertion failed: n != j");
}
double xStart = X.front();
double xEnd = X.back();
/* use uniform breakpoints */
gsl_bspline_knots_uniform(xStart, xEnd, bw);
/* construct the fit matrix X */
for (int i = 0; i < n; ++i) {
double xi = gsl_vector_get(x, i);
/* compute B_j(xi) for all j */
gsl_bspline_eval(xi, B, bw);
/* fill in row i of X */
for (j = 0; j < ncoeffs; ++j) {
double Bj = gsl_vector_get(B, j);
gsl_matrix_set(Z, i, j, Bj);
}
}
/* do the fit */
gsl_multifit_wlinear(Z, w, y, c, cov, &chisq, mw);
/* output the smoothed curve */
API::MatrixWorkspace_sptr outWS =
WorkspaceFactory::Instance().create(inWS, 1, X.size(), Y.size());
{
outWS->getSpectrum(0)
.setSpectrumNo(inWS->getSpectrum(spec).getSpectrumNo());
double yi, yerr;
for (MantidVec::size_type i = 0; i < Y.size(); i++) {
double xi = X[i];
gsl_bspline_eval(xi, B, bw);
gsl_multifit_linear_est(B, c, cov, &yi, &yerr);
outWS->dataY(0)[i] = yi;
outWS->dataE(0)[i] = yerr;
}
outWS->dataX(0) = X;
}
gsl_bspline_free(bw);
gsl_vector_free(B);
gsl_vector_free(x);
gsl_vector_free(y);
gsl_matrix_free(Z);
gsl_vector_free(c);
gsl_vector_free(w);
gsl_matrix_free(cov);
gsl_multifit_linear_free(mw);
setProperty("OutputWorkspace", outWS);
}
} // namespace Algorithms
} // namespace CurveFitting
} // namespace Mantid