Loading...
Searching...
No Matches
rbfspline.C
Go to the documentation of this file.
1/*
2 * This file is part of the Multivariate Splines library.
3 * Copyright (C) 2012 Bjarne Grimstad (bjarne.grimstad@gmail.com)
4 *
5 * This Source Code Form is subject to the terms of the Mozilla Public
6 * License, v. 2.0. If a copy of the MPL was not distributed with this
7 * file, You can obtain one at http://mozilla.org/MPL/2.0/.
8*/
9
10
11#include "rbfspline.h"
12#include "linearsolvers.h"
13#include <Eigen/Eigen>
14
15namespace SPLINTER
16{
17
18RBFSpline::RBFSpline(const DataTable& samples, RadialBasisFunctionType type,
19 double e)
20 : RBFSpline(samples, type, false)
21{
22}
23
24RBFSpline::RBFSpline(const DataTable& samples, RadialBasisFunctionType type,
25 DenseMatrix w, double e)
26 : samples(samples),
27 normalized(false),
28 precondition(false),
29 dim(samples.getNumVariables()),
30 numSamples(samples.getNumSamples())
31{
32 if (type == RadialBasisFunctionType::THIN_PLATE_SPLINE)
33 {
34 fn = std::shared_ptr<RadialBasisFunction>(new ThinPlateSpline());
35 fn->e = e;
36 }
37 else if (type == RadialBasisFunctionType::MULTIQUADRIC)
38 {
39 fn = std::shared_ptr<RadialBasisFunction>(new Multiquadric());
40 }
41 else if (type == RadialBasisFunctionType::INVERSE_QUADRIC)
42 {
43 fn = std::shared_ptr<RadialBasisFunction>(new InverseQuadric());
44 }
45 else if (type == RadialBasisFunctionType::INVERSE_MULTIQUADRIC)
46 {
47 fn = std::shared_ptr<RadialBasisFunction>(new InverseMultiquadric());
48 }
49 else if (type == RadialBasisFunctionType::GAUSSIAN)
50 {
51 fn = std::shared_ptr<RadialBasisFunction>(new Gaussian());
52 fn->e = e;
53 }
54 else
55 {
56 fn = std::shared_ptr<RadialBasisFunction>(new ThinPlateSpline());
57 }
58
59 weights = w;
60}
61
62RBFSpline::RBFSpline(const DataTable& samples, RadialBasisFunctionType type,
63 bool normalized, double e)
64 : samples(samples),
65 normalized(normalized),
66 precondition(false),
67 dim(samples.getNumVariables()),
68 numSamples(samples.getNumSamples())
69{
70 if (type == RadialBasisFunctionType::THIN_PLATE_SPLINE)
71 {
72 fn = std::shared_ptr<RadialBasisFunction>(new ThinPlateSpline());
73 }
74 else if (type == RadialBasisFunctionType::MULTIQUADRIC)
75 {
76 fn = std::shared_ptr<RadialBasisFunction>(new Multiquadric());
77 }
78 else if (type == RadialBasisFunctionType::INVERSE_QUADRIC)
79 {
80 fn = std::shared_ptr<RadialBasisFunction>(new InverseQuadric());
81 }
82 else if (type == RadialBasisFunctionType::INVERSE_MULTIQUADRIC)
83 {
84 fn = std::shared_ptr<RadialBasisFunction>(new InverseMultiquadric());
85 }
86 else if (type == RadialBasisFunctionType::GAUSSIAN)
87 {
88 fn = std::shared_ptr<RadialBasisFunction>(new Gaussian());
89 fn->e = e;
90 }
91 else
92 {
93 fn = std::shared_ptr<RadialBasisFunction>(new ThinPlateSpline());
94 }
95
96 /* Want to solve the linear system A*w = b,
97 * where w is the vector of weights.
98 * NOTE: the system is dense and by default badly conditioned.
99 * It should be solved by a specialized solver such as GMRES
100 * with preconditioning (e.g. ACBF) as in matlab.
101 * NOTE: Consider trying the Łukaszyk–Karmowski metric (for two variables)
102 */
103 //SparseMatrix A(numSamples,numSamples);
104 //A.reserve(numSamples*numSamples);
105 DenseMatrix A;
106 A.setZero(numSamples, numSamples);
107 DenseMatrix b;
108 b.setZero(numSamples, 1);
109 int i = 0;
110
111 for (auto it1 = samples.cbegin(); it1 != samples.cend(); ++it1, ++i)
112 {
113 double sum = 0;
114 int j = 0;
115
116 for (auto it2 = samples.cbegin(); it2 != samples.cend(); ++it2, ++j)
117 {
118 double val = fn->eval(dist(*it1, *it2));
119
120 if (val != 0)
121 {
122 //A.insert(i,j) = val;
123 A(i, j) = val;
124 sum += val;
125 }
126 }
127
128 double y = it1->getY();
129
130 if (normalized)
131 {
132 b(i) = sum * y;
133 }
134 else
135 {
136 b(i) = y;
137 }
138 }
139
140 //A.makeCompressed();
141
142 if (precondition)
143 {
144 // Calcualte precondition matrix P
145 DenseMatrix P = computePreconditionMatrix();
146 // Preconditioned A and b
147 DenseMatrix Ap = P * A;
148 DenseMatrix bp = P * b;
149 A = Ap;
150 b = bp;
151 }
152
153#ifndef NDEBUG
154 std::cout << "Computing RBF weights using dense solver." << std::endl;
155 std::cout << "The radius of the RBF is equal to " << e << std::endl;
156#endif // NDEBUG
157 // SVD analysis
158 // Eigen::JacobiSVD<DenseMatrix> svd(A, Eigen::ComputeThinU | Eigen::ComputeThinV);
159 // auto svals = svd.singularValues();
160 // double svalmax = svals(0);
161 // double svalmin = svals(svals.rows() - 1);
162 // double rcondnum = (svalmax <= 0.0 || svalmin <= 0.0) ? 0.0 : svalmin / svalmax;
163 // #ifndef NDEBUG
164 // std::cout << "The reciprocal of the condition number is: " << rcondnum <<
165 // std::endl;
166 // std::cout << "Largest/smallest singular value: " << svalmax << " / " << svalmin
167 // << std::endl;
168 // #endif // NDEBUG
169 // // Solve for weights
170 // weights = svd.solve(b);
171 weights = A.colPivHouseholderQr().solve(b);
172#ifndef NDEBUG
173 // Compute error. If it is used later on, move this statement above the NDEBUG
174 double err = (A * weights - b).norm() / b.norm();
175 std::cout << "Error: " << std::setprecision(10) << err << std::endl;
176#endif // NDEBUG
177 // // Alternative solver
178 // DenseQR s;
179 // bool success = s.solve(A,b,weights);
180 // assert(success);
181 // NOTE: Tried using experimental GMRES solver in Eigen, but it did not work very well.
182}
183
184double RBFSpline::eval(DenseVector x) const
185{
186 std::vector<double> y;
187
188 for (int i = 0; i < x.rows(); i++)
189 {
190 y.push_back(x(i));
191 }
192
193 return eval(y);
194}
195
196double RBFSpline::eval(std::vector<double> x) const
197{
198 assert(x.size() == dim);
199 double fval, sum = 0, sumw = 0;
200 int i = 0;
201
202 for (auto it = samples.cbegin(); it != samples.cend(); ++it, ++i)
203 {
204 fval = fn->eval(dist(x, it->getX()));
205 sumw += weights(i) * fval;
206 sum += fval;
207 }
208
209 return normalized ? sumw / sum : sumw;
210}
211
212/*
213 * TODO: test for errors
214 */
215//DenseMatrix RBFSpline::evalJacobian(DenseVector x) const
216//{
217// std::vector<double> x_vec;
218// for(unsigned int i = 0; i<x.size(); i++)
219// x_vec.push_back(x(i));
220
221// DenseMatrix jac;
222// jac.setZero(1,dim);
223
224// for(unsigned int i = 0; i < dim; i++)
225// {
226// double sumw = 0;
227// double sumw_d = 0;
228// double sum = 0;
229// double sum_d = 0;
230
231// int j = 0;
232// for(auto it = samples.cbegin(); it != samples.cend(); ++it, ++j)
233// {
234// // Sample
235// auto s_vec = it->getX();
236
237// // Distance from sample
238// double r = dist(x_vec, s_vec);
239// double ri = x_vec.at(i) - s_vec.at(i);
240
241// // Evaluate RBF and its derivative at r
242// double f = fn->eval(r);
243// double dfdr = fn->evalDerivative(r);
244
245// sum += f;
246// sumw += weights(j)*f;
247
248// // TODO: check if this assumption is correct
249// if(r != 0)
250// {
251// sum_d += dfdr*ri/r;
252// sumw_d += weights(j)*dfdr*ri/r;
253// }
254// }
255
256// if(normalized)
257// jac(i) = (sum*sumw_d - sum_d*sumw)/(sum*sum);
258// else
259// jac(i) = sumw_d;
260// }
261// return jac;
262//}
263
264/*
265 * Calculate precondition matrix
266 */
267DenseMatrix RBFSpline::computePreconditionMatrix() const
268{
269 DenseMatrix P;
270 P.setZero(numSamples, numSamples);
271 // Calculate precondition matrix P based on
272 // purely local approximate cardinal basis functions (ACBF)
273 int sigma = std::max(1.0,
274 std::floor(0.1 * numSamples)); // Local points to consider
275 int i = 0;
276
277 for (auto it1 = samples.cbegin(); it1 != samples.cend(); ++it1, ++i)
278 {
279 Point p1(it1->getX());
280 // Shift data using p1 as origin
281 std::vector<Point> shifted_points;
282 int j = 0;
283
284 for (auto it2 = samples.cbegin(); it2 != samples.cend(); ++it2, ++j)
285 {
286 Point p2(it2->getX());
287 Point p3(p2 - p1);
288 p3.setIndex(j); // store index with point
289 shifted_points.push_back(p3);
290 }
291
292 std::sort(shifted_points.begin(), shifted_points.end());
293 // Find sigma closest points to p1
294 std::vector<Point> points;
295 std::vector<int> indices;
296
297 for (int j = 0; j < sigma; j++)
298 {
299 Point p(shifted_points.at(j));
300 indices.push_back(p.getIndex());
301 Point p2 = p + p1;
302 p2.setIndex(p.getIndex());
303 points.push_back(
304 p2); // The resulting point has a different index than that of p
305 // cout << p.getIndex() << "/" << p1.getIndex() << "/" << p2.getIndex() << endl;
306 // assert(p.getIndex() == p2.getIndex());
307 }
308
309 // Add some points far away from p1 and preferably scattered around the domain boundary
310 for (int k = 0; k < 1; k++)
311 {
312 Point p(shifted_points.at(shifted_points.size() - 1 - k));
313 indices.push_back(p.getIndex());
314 Point p2 = p + p1;
315 p2.setIndex(p.getIndex());
316 points.push_back(p2);
317 }
318
319 // Build and solve linear system
320 int m = points.size();
321 DenseMatrix e;
322 e.setZero(m, 1);
323 e(0, 0) = 1;
324 DenseMatrix B;
325 B.setZero(m, m);
326 DenseMatrix w;
327 w.setZero(m, 1);
328 assert(points.front().getIndex() == i);
329
330 for (int k1 = 0; k1 < m; k1++)
331 {
332 for (int k2 = 0; k2 < m; k2++)
333 {
334 Point p = points.at(k1) - points.at(k2);
335 B(k1, k2) = fn->eval(p.dist());
336 }
337 }
338
339 // DenseQR s;
340 // bool success = s.solve(B,e,w);
341 // assert(success);
342 Eigen::JacobiSVD<DenseMatrix> svd(B, Eigen::ComputeThinU | Eigen::ComputeThinV);
343 w = svd.solve(e);
344 //assert(svd.info() == Eigen::Success);
345
346 for (unsigned int j = 0; j < numSamples; j++)
347 {
348 auto it = find(indices.begin(), indices.end(), j);
349
350 if (it != indices.end())
351 {
352 int k = it - indices.begin();
353 P(i, j) = w(k, 0);
354 //cout << "j/k/g " << j << "/" << k << "/" << points.at(k).getIndex() << endl;
355 assert(points.at(k).getIndex() == j);
356 }
357 }
358 }
359
360 return P;
361}
362
363/*
364 * Computes Euclidean distance ||x-y||
365 */
366double RBFSpline::dist(std::vector<double> x, std::vector<double> y) const
367{
368 assert(x.size() == y.size());
369 double sum = 0.0;
370
371 for (unsigned int i = 0; i < x.size(); i++)
372 {
373 sum += (x.at(i) - y.at(i)) * (x.at(i) - y.at(i));
374 }
375
376 return std::sqrt(sum);
377}
378
379/*
380 * Computes Euclidean distance ||x-y||
381 */
382double RBFSpline::dist(DataPoint x, DataPoint y) const
383{
384 return dist(x.getX(), y.getX());
385}
386
387bool RBFSpline::dist_sort(DataPoint x, DataPoint y) const
388{
389 std::vector<double> zeros(x.getDimX(), 0);
390 DataPoint origin(zeros, 0.0);
391 double x_dist = dist(x, origin);
392 double y_dist = dist(y, origin);
393 return (x_dist < y_dist);
394}
395
396} // namespace MultivariateSplines
volScalarField & A
double dist(const std::vector< double > x, const std::vector< double > y)
Definition datapoint.C:72
Info<< endl;Info<< "*********************************************************"<< endl;Info<< "Performing test for the parameterized BC inverse solver"<< endl;Info<< endl;word outputFolder="./ITHACAoutput/parameterizedBCtest_RBFparameter/";volScalarField gTrueField=example_paramBC.list2Field(example_paramBC.gTrue);ITHACAstream::exportSolution(gTrueField, "1", outputFolder, "gTrue");Eigen::VectorXd rbfWidth=EigenFunctions::ExpSpaced(0.001, 100, rbfWidthTest_size);ITHACAstream::exportMatrix(rbfWidth, "rbfShapeParameters", "eigen", outputFolder);Eigen::VectorXd residualNorms;residualNorms.resize(rbfWidthTest_size);Eigen::VectorXd heatFluxL2norm(rbfWidthTest_size);Eigen::VectorXd heatFluxLinfNorm=heatFluxL2norm;Eigen::VectorXd condNumber=heatFluxL2norm;Eigen::MatrixXd singVal;for(int i=0;i< rbfWidthTest_size;i++){ Info<< "*********************************************************"<< endl;Info<< "RBF parameter "<< rbfWidth(i)<< endl;Info<< "Test "<< i+1<< endl;Info<< endl;example_paramBC.set_gParametrized("rbf", rbfWidth(i));example_paramBC.parameterizedBCoffline(1);example_paramBC.parameterizedBC("fullPivLU");volScalarField gParametrizedField=example_paramBC.list2Field(example_paramBC.g);ITHACAstream::exportSolution(gParametrizedField, std::to_string(i+1), outputFolder, "gParametrized");volScalarField &T(example_paramBC._T());ITHACAstream::exportSolution(T, std::to_string(i+1), outputFolder, "T");residualNorms(i)=Foam::sqrt(example_paramBC.residual.squaredNorm());Eigen::MatrixXd A=example_paramBC.Theta.transpose() *example_paramBC.Theta;Eigen::JacobiSVD< Eigen::MatrixXd > svd(A, Eigen::ComputeThinU|Eigen::ComputeThinV)
volScalarField & p
label i
Definition pEqn.H:46