Project
Loading...
Searching...
No Matches
Spline2DHelper.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
16
17#include "Spline2DHelper.h"
18#include "Spline1DHelper.h"
19
20#include "SymMatrixSolver.h"
21#include "GPUCommonDef.h"
22#include "GPUCommonLogger.h"
23
24#include "TMath.h"
25#include "TMatrixD.h"
26#include "TVectorD.h"
27#include "TDecompBK.h"
28
29#include <vector>
30#include "TRandom.h"
31#include "TMath.h"
32#include "TCanvas.h"
33#include "TNtuple.h"
34#include "TFile.h"
35#include "GPUCommonMath.h"
36#include <iostream>
37#include <chrono>
38
39using namespace o2::gpu;
40
41template <typename DataT>
42Spline2DHelper<DataT>::Spline2DHelper() : mError(), mFdimensions(0), mHelperU1(), mHelperU2()
43{
44}
45
46template <typename DataT>
47int32_t Spline2DHelper<DataT>::storeError(int32_t code, const char* msg)
48{
49 mError = msg;
50 return code;
51}
52
53template <typename DataT>
55 DataT* Fparameters, double x1Min, double x1Max, double x2Min, double x2Max,
56 std::function<void(double x1, double x2, double f[/*spline.getYdimensions()*/])> F) const
57{
60
61 std::vector<double> dataPointF(getNumberOfDataPoints() * mFdimensions);
62
63 double scaleX1 = (x1Max - x1Min) / ((double)mHelperU1.getSpline().getUmax());
64 double scaleX2 = (x2Max - x2Min) / ((double)mHelperU2.getSpline().getUmax());
65
66 for (int32_t iv = 0; iv < getNumberOfDataPointsU2(); iv++) {
67 double x2 = x2Min + mHelperU2.getDataPoint(iv).u * scaleX2;
68 for (int32_t iu = 0; iu < getNumberOfDataPointsU1(); iu++) {
69 double x1 = x1Min + mHelperU1.getDataPoint(iu).u * scaleX1;
70 F(x1, x2, &dataPointF[(iv * getNumberOfDataPointsU1() + iu) * mFdimensions]);
71 }
72 }
73 approximateFunction(Fparameters, dataPointF.data());
74}
75
76template <typename DataT>
78 DataT* Fparameters, double x1Min, double x1Max, double x2Min, double x2Max,
79 std::function<void(const std::vector<double>& x1, const std::vector<double>& x2, std::vector<double> f[/*mFdimensions*/])> F,
80 uint32_t batchsize) const
81{
85
86 std::vector<double> dataPointF(getNumberOfDataPoints() * mFdimensions);
87
88 double scaleX1 = (x1Max - x1Min) / ((double)mHelperU1.getSpline().getUmax());
89 double scaleX2 = (x2Max - x2Min) / ((double)mHelperU2.getSpline().getUmax());
90
91 std::vector<double> x1;
92 x1.reserve(batchsize);
93
94 std::vector<double> x2;
95 x2.reserve(batchsize);
96
97 std::vector<int32_t> index;
98 index.reserve(batchsize);
99
100 std::vector<double> dataPointFTmp[mFdimensions];
101 for (int32_t iDim = 0; iDim < mFdimensions; ++iDim) {
102 dataPointFTmp[iDim].reserve(batchsize);
103 }
104
105 uint32_t counter = 0;
106 for (int32_t iv = 0; iv < getNumberOfDataPointsU2(); iv++) {
107 double x2Tmp = x2Min + mHelperU2.getDataPoint(iv).u * scaleX2;
108 for (int32_t iu = 0; iu < getNumberOfDataPointsU1(); iu++) {
109 double x1Tmp = x1Min + mHelperU1.getDataPoint(iu).u * scaleX1;
110 x1.emplace_back(x1Tmp);
111 x2.emplace_back(x2Tmp);
112 index.emplace_back((iv * getNumberOfDataPointsU1() + iu) * mFdimensions);
113 ++counter;
114
115 if (counter == batchsize || (iu == (getNumberOfDataPointsU1() - 1) && (iv == (getNumberOfDataPointsU2() - 1)))) {
116 counter = 0;
117 F(x1, x2, dataPointFTmp);
118 uint32_t entries = index.size();
119
120 for (uint32_t i = 0; i < entries; ++i) {
121 const uint32_t indexTmp = index[i];
122 for (int32_t iDim = 0; iDim < mFdimensions; ++iDim) {
123 dataPointF[indexTmp + iDim] = dataPointFTmp[iDim][i];
124 }
125 }
126
127 x1.clear();
128 x2.clear();
129 index.clear();
130 for (int32_t iDim = 0; iDim < mFdimensions; ++iDim) {
131 dataPointFTmp[iDim].clear();
132 }
133 }
134 }
135 }
136 approximateFunction(Fparameters, dataPointF.data());
137}
138
139template <typename DataT>
141 DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim*/]) const
142{
144
145 const int32_t Ndim = mFdimensions;
146 const int32_t Ndim2 = 2 * Ndim;
147 const int32_t Ndim3 = 3 * Ndim;
148 const int32_t Ndim4 = 4 * Ndim;
149
150 int32_t nDataPointsU = getNumberOfDataPointsU1();
151 int32_t nDataPointsV = getNumberOfDataPointsU2();
152
153 int32_t nKnotsU = mHelperU1.getSpline().getNumberOfKnots();
154 int32_t nKnotsV = mHelperU2.getSpline().getNumberOfKnots();
155
156 std::unique_ptr<double[]> rotDataPointF(new double[nDataPointsU * nDataPointsV * Ndim]); // U DataPoints x V DataPoints : rotated DataPointF for one output dimension
157 std::unique_ptr<double[]> Dv(new double[nKnotsV * nDataPointsU * Ndim]); // V knots x U DataPoints
158
159 std::unique_ptr<DataT[]> parU(new DataT[mHelperU1.getSpline().calcNumberOfParameters(Ndim)]);
160 std::unique_ptr<DataT[]> parV(new DataT[mHelperU2.getSpline().calcNumberOfParameters(Ndim)]);
161
162 std::unique_ptr<double[]> parUdbl(new double[mHelperU1.getSpline().calcNumberOfParameters(Ndim)]);
163 std::unique_ptr<double[]> parVdbl(new double[mHelperU2.getSpline().calcNumberOfParameters(Ndim)]);
164
165 // rotated data points (u,v)->(v,u)
166
167 for (int32_t ipu = 0; ipu < nDataPointsU; ipu++) {
168 for (int32_t ipv = 0; ipv < nDataPointsV; ipv++) {
169 for (int32_t dim = 0; dim < Ndim; dim++) {
170 rotDataPointF[Ndim * (ipu * nDataPointsV + ipv) + dim] = DataPointF[Ndim * (ipv * nDataPointsU + ipu) + dim];
171 }
172 }
173 }
174
175 // get S and S'u at all the knots by interpolating along the U axis
176
177 for (int32_t iKnotV = 0; iKnotV < nKnotsV; ++iKnotV) {
178 int32_t ipv = mHelperU2.getKnotDataPoint(iKnotV);
179 const double* DataPointFrow = &(DataPointF[Ndim * ipv * nDataPointsU]);
180 mHelperU1.approximateFunctionGradually(parU.get(), DataPointFrow);
181
182 for (int32_t i = 0; i < mHelperU1.getSpline().calcNumberOfParameters(Ndim); i++) {
183 parUdbl[i] = parU[i];
184 }
185 for (int32_t iKnotU = 0; iKnotU < nKnotsU; ++iKnotU) {
186 DataT* knotPar = &Fparameters[Ndim4 * (iKnotV * nKnotsU + iKnotU)];
187 for (int32_t dim = 0; dim < Ndim; ++dim) {
188 knotPar[dim] = parU[Ndim * (2 * iKnotU) + dim]; // store S for all the knots
189 knotPar[Ndim2 + dim] = parU[Ndim * (2 * iKnotU) + Ndim + dim]; // store S'u for all the knots //SG!!!
190 }
191 }
192
193 // recalculate F values for all ipu DataPoints at V = ipv
194 for (int32_t ipu = 0; ipu < nDataPointsU; ipu++) {
195 double splineF[Ndim];
196 double u = mHelperU1.getDataPoint(ipu).u;
197 mHelperU1.getSpline().interpolateAtU(Ndim, parUdbl.get(), u, splineF);
198 for (int32_t dim = 0; dim < Ndim; dim++) {
199 rotDataPointF[(ipu * nDataPointsV + ipv) * Ndim + dim] = splineF[dim];
200 }
201 }
202 }
203
204 // calculate S'v at all data points with V == V of a knot
205
206 for (int32_t ipu = 0; ipu < nDataPointsU; ipu++) {
207 const double* DataPointFcol = &(rotDataPointF[ipu * nDataPointsV * Ndim]);
208 mHelperU2.approximateFunctionGradually(parV.get(), DataPointFcol);
209 for (int32_t iKnotV = 0; iKnotV < nKnotsV; iKnotV++) {
210 for (int32_t dim = 0; dim < Ndim; dim++) {
211 double dv = parV[(iKnotV * 2 + 1) * Ndim + dim];
212 Dv[(iKnotV * nDataPointsU + ipu) * Ndim + dim] = dv;
213 }
214 }
215 }
216
217 // fit S'v and S''_vu at all the knots
218
219 for (int32_t iKnotV = 0; iKnotV < nKnotsV; ++iKnotV) {
220 const double* Dvrow = &(Dv[iKnotV * nDataPointsU * Ndim]);
221 mHelperU1.approximateFunction(parU.get(), Dvrow);
222 for (int32_t iKnotU = 0; iKnotU < nKnotsU; ++iKnotU) {
223 for (int32_t dim = 0; dim < Ndim; ++dim) {
224 Fparameters[Ndim4 * (iKnotV * nKnotsU + iKnotU) + Ndim + dim] = parU[Ndim * 2 * iKnotU + dim]; // store S'v for all the knots
225 Fparameters[Ndim4 * (iKnotV * nKnotsU + iKnotU) + Ndim3 + dim] = parU[Ndim * 2 * iKnotU + Ndim + dim]; // store S''vu for all the knots
226 }
227 }
228 }
229}
230
231template <typename DataT>
234 double x1Min, double x1Max, double x2Min, double x2Max,
235 std::function<void(double x1, double x2, double f[/*spline.getYdimensions()*/])> F,
236 int32_t nAuxiliaryDataPointsU1, int32_t nAuxiliaryDataPointsU2)
237{
239
240 setSpline(spline, nAuxiliaryDataPointsU1, nAuxiliaryDataPointsU2);
241 mFdimensions = spline.getYdimensions();
242 std::vector<double> dataPointX1(getNumberOfDataPoints());
243 std::vector<double> dataPointX2(getNumberOfDataPoints());
244 std::vector<double> dataPointWeight(getNumberOfDataPoints(), 1.);
245 std::vector<double> dataPointF(getNumberOfDataPoints() * mFdimensions);
246
247 double scaleX1 = (x1Max - x1Min) / ((double)mHelperU1.getSpline().getUmax());
248 double scaleX2 = (x2Max - x2Min) / ((double)mHelperU2.getSpline().getUmax());
249
250 for (int32_t iv = 0; iv < getNumberOfDataPointsU2(); iv++) {
251 double x2 = x2Min + mHelperU2.getDataPoint(iv).u * scaleX2;
252 for (int32_t iu = 0; iu < getNumberOfDataPointsU1(); iu++) {
253 double x1 = x1Min + mHelperU1.getDataPoint(iu).u * scaleX1;
254 int32_t ind = iv * getNumberOfDataPointsU1() + iu;
255 dataPointX1[ind] = x1;
256 dataPointX2[ind] = x2;
257 F(x1, x2, &dataPointF[ind * mFdimensions]);
258 }
259 }
260 approximateDataPoints(spline, spline.getParameters(), x1Min, x1Max, x2Min, x2Max, dataPointX1.data(), dataPointX2.data(), dataPointF.data(),
261 dataPointWeight.data(), getNumberOfDataPoints());
262}
263
264template <typename DataT>
265void Spline2DHelper<DataT>::setGrid(Spline2DContainerBase<DataT, FlatObject>& spline, double x1Min, double x1Max, double x2Min, double x2Max)
266{
267 mFdimensions = spline.getYdimensions();
268 spline.setXrange(x1Min, x1Max, x2Min, x2Max);
269 {
270 std::vector<int32_t> knots;
271 for (int32_t i = 0; i < spline.getGridX1().getNumberOfKnots(); i++) {
272 knots.push_back(spline.getGridX1().getKnot(i).getU());
273 }
274 fGridU.recreate(0, knots.size(), knots.data());
275 fGridU.setXrange(x1Min, x1Max);
276 }
277 {
278 std::vector<int32_t> knots;
279 for (int32_t i = 0; i < spline.getGridX2().getNumberOfKnots(); i++) {
280 knots.push_back(spline.getGridX2().getKnot(i).getU());
281 }
282 fGridV.recreate(0, knots.size(), knots.data());
283 fGridV.setXrange(x2Min, x2Max);
284 }
285}
286
287template <typename DataT>
288void Spline2DHelper<DataT>::getScoefficients(int32_t iu, int32_t iv, double u, double v,
289 double coeff[16], int32_t indices[16])
290{
291 const Knot<double>& knotU = fGridU.getKnot(iu);
292 const Knot<double>& knotV = fGridV.getKnot(iv);
293 int32_t nu = fGridU.getNumberOfKnots();
294
295 // indices of parameters that are involved in spline calculation, 1D case
296 int32_t i00 = (nu * iv + iu) * 4; // values { S, S'v, S'u, S''vu } at {u0, v0}
297 int32_t i01 = i00 + 4 * nu; // values { ... } at {u0, v1}
298 int32_t i10 = i00 + 4;
299 int32_t i11 = i01 + 4;
300
301 double dSl, dDl, dSr, dDr;
302 Spline1DHelper<double>::getScoefficients(knotU, u, dSl, dDl, dSr, dDr);
303 double dSd, dDd, dSu, dDu;
304 Spline1DHelper<double>::getScoefficients(knotV, v, dSd, dDd, dSu, dDu);
305
306 // A = Parameters + i00, B = Parameters + i01
307 // S = dSl * (dSd * A[0] + dDd * A[1]) + dDl * (dSd * A[2] + dDd * A[3]) +
308 // dSr * (dSd * A[4] + dDd * A[5]) + dDr * (dSd * A[6] + dDd * A[7]) +
309 // dSl * (dSu * B[0] + dDu * B[1]) + dDl * (dSu * B[2] + dDu * B[3]) +
310 // dSr * (dSu * B[4] + dDu * B[5]) + dDr * (dSu * B[6] + dDu * B[7]);
311
312 double c[16] = {dSl * dSd, dSl * dDd, dDl * dSd, dDl * dDd,
313 dSr * dSd, dSr * dDd, dDr * dSd, dDr * dDd,
314 dSl * dSu, dSl * dDu, dDl * dSu, dDl * dDu,
315 dSr * dSu, dSr * dDu, dDr * dSu, dDr * dDu};
316 for (int32_t i = 0; i < 16; i++) {
317 coeff[i] = c[i];
318 }
319 for (int32_t i = 0; i < 4; i++) {
320 indices[0 + i] = i00 + i;
321 indices[4 + i] = i10 + i;
322 indices[8 + i] = i01 + i;
323 indices[12 + i] = i11 + i;
324 }
325}
326
327template <typename DataT>
329 Spline2DContainerBase<DataT, FlatObject>& spline, DataT* splineParameters, double x1Min, double x1Max, double x2Min, double x2Max,
330 const double dataPointX1[], const double dataPointX2[], const double dataPointF[/*getNumberOfDataPoints() x nFdim*/],
331 const double dataPointWeight[], int32_t nDataPoints)
332{
334
335 setGrid(spline, x1Min, x1Max, x2Min, x2Max);
336
337 int32_t nFdim = spline.getYdimensions();
338 int32_t nu = fGridU.getNumberOfKnots();
339 int32_t nv = fGridV.getNumberOfKnots();
340
341 const int32_t nPar = 4 * spline.getNumberOfKnots(); // n parameters for 1-dimensional F
342
343 SymMatrixSolver solver(nPar, nFdim);
344
345 for (int32_t iPoint = 0; iPoint < nDataPoints; ++iPoint) {
346 double u = fGridU.convXtoU(dataPointX1[iPoint]);
347 double v = fGridV.convXtoU(dataPointX2[iPoint]);
348 double weight = dataPointWeight[iPoint];
349 if (!(weight > 0.)) {
350 continue;
351 }
352 int32_t iu = fGridU.getLeftKnotIndexForU(u);
353 int32_t iv = fGridV.getLeftKnotIndexForU(v);
354 double c[16];
355 int32_t ind[16];
356 getScoefficients(iu, iv, u, v, c, ind);
357
358 // S(u,v) = sum c[i]*Parameters[ind[i]]
359
360 for (int32_t i = 0; i < 16; i++) {
361 for (int32_t j = i; j < 16; j++) {
362 solver.A(ind[i], ind[j]) += weight * c[i] * c[j];
363 }
364 }
365
366 for (int32_t iDim = 0; iDim < nFdim; iDim++) {
367 double f = (double)dataPointF[iPoint * nFdim + iDim];
368 for (int32_t i = 0; i < 16; i++) {
369 solver.B(ind[i], iDim) += weight * f * c[i];
370 }
371 }
372 } // data points
373
374 // add extra smoothness in case some data is missing
375
376 for (int32_t iu = 0; iu < nu - 1; iu++) {
377 for (int32_t iv = 0; iv < nv - 1; iv++) {
378 int32_t smoothPoint[4][2] = {
379 {-1, -1},
380 {-1, +1},
381 {+1, -1},
382 {+1, +1}};
383 for (int32_t iSet = 0; iSet < 4; iSet++) {
384 int32_t pu = iu + smoothPoint[iSet][0];
385 int32_t pv = iv + smoothPoint[iSet][1];
386 int32_t ip = (nu * pv + pu) * 4;
387 if (pu < 0 || pv < 0 || pu >= nu || pv >= nv) {
388 continue;
389 }
390 double c[17];
391 int32_t ind[17];
392 getScoefficients(iu, iv, fGridU.getKnot(pu).u, fGridV.getKnot(pv).u, c, ind);
393 c[16] = -1.;
394 ind[16] = ip;
395 // S = sum c[i]*Par[ind[i]]
396 double w = 1.e-8;
397 for (int32_t i = 0; i < 17; i++) {
398 for (int32_t j = i; j < 17; j++) {
399 solver.A(ind[i], ind[j]) += w * c[i] * c[j];
400 }
401 }
402 }
403 }
404 }
405
406 solver.solve();
407 for (int32_t i = 0; i < nPar; i++) {
408 for (int32_t iDim = 0; iDim < nFdim; iDim++) {
409 splineParameters[i * nFdim + iDim] = solver.B(i, iDim);
410 }
411 }
412}
413
414template <typename DataT>
415int32_t Spline2DHelper<DataT>::test(const bool draw, const bool drawDataPoints)
416{
417 using namespace std;
418
419 const int32_t Ndim = 3;
420
421 const int32_t Fdegree = 4;
422
423 double Fcoeff[Ndim][4 * (Fdegree + 1) * (Fdegree + 1)];
424
425 constexpr int32_t nKnots = 10;
426 constexpr int32_t nAuxiliaryPoints = 2;
427 constexpr int32_t uMax = nKnots; //* 3;
428
429 auto F = [&](double u, double v, double Fuv[]) {
430 const double scale = TMath::Pi() / uMax;
431 double uu = u * scale;
432 double vv = v * scale;
433 double cosu[Fdegree + 1], sinu[Fdegree + 1], cosv[Fdegree + 1], sinv[Fdegree + 1];
434 double ui = 0, vi = 0;
435 for (int32_t i = 0; i <= Fdegree; i++, ui += uu, vi += vv) {
436 GPUCommonMath::SinCosd(ui, sinu[i], cosu[i]);
437 GPUCommonMath::SinCosd(vi, sinv[i], cosv[i]);
438 }
439 for (int32_t dim = 0; dim < Ndim; dim++) {
440 double f = 0; // Fcoeff[dim][0]/2;
441 for (int32_t i = 1; i <= Fdegree; i++) {
442 for (int32_t j = 1; j <= Fdegree; j++) {
443 double* c = &(Fcoeff[dim][4 * (i * Fdegree + j)]);
444 f += c[0] * cosu[i] * cosv[j];
445 f += c[1] * cosu[i] * sinv[j];
446 f += c[2] * sinu[i] * cosv[j];
447 f += c[3] * sinu[i] * sinv[j];
448 }
449 }
450 Fuv[dim] = f;
451 }
452 };
453
454 TCanvas* canv = nullptr;
455 TNtuple* nt = nullptr;
456 TNtuple* knots = nullptr;
457
458 auto ask = [&]() -> bool {
459 if (!canv) {
460 return 0;
461 }
462 canv->Update();
463 LOG(info) << "type 'q ' to exit";
464 std::string str;
465 std::getline(std::cin, str);
466 return (str != "q" && str != ".q");
467 };
468
469 LOG(info) << "Test 2D interpolation with the compact spline";
470
471 int32_t nTries = 100;
472
473 if (draw) {
474 canv = new TCanvas("cQA", "Spline2D QA", 1500, 800);
475 nTries = 10000;
476 }
477
478 double statDf = 0;
479 double statDf1D = 0;
480 double statN = 0;
481
482 auto statTime = std::chrono::nanoseconds::zero();
483
484 for (int32_t seed = 1; seed < nTries + 1; seed++) {
485 // LOG(info) << "next try.." ;
486
487 gRandom->SetSeed(seed);
488
489 for (int32_t dim = 0; dim < Ndim; dim++) {
490 for (int32_t i = 0; i < 4 * (Fdegree + 1) * (Fdegree + 1); i++) {
491 Fcoeff[dim][i] = gRandom->Uniform(-1, 1);
492 }
493 }
494
496
497 int32_t knotsU[nKnots], knotsV[nKnots];
498 do {
499 knotsU[0] = 0;
500 knotsV[0] = 0;
501 double du = 1. * uMax / (nKnots - 1);
502 for (int32_t i = 1; i < nKnots; i++) {
503 knotsU[i] = (int32_t)(i * du); // + gRandom->Uniform(-du / 3, du / 3);
504 knotsV[i] = (int32_t)(i * du); // + gRandom->Uniform(-du / 3, du / 3);
505 }
506 knotsU[nKnots - 1] = uMax;
507 knotsV[nKnots - 1] = uMax;
508 spline.recreate(nKnots, knotsU, nKnots, knotsV);
509
510 if (nKnots != spline.getGridX1().getNumberOfKnots() ||
511 nKnots != spline.getGridX2().getNumberOfKnots()) {
512 LOG(info) << "warning: n knots changed during the initialisation " << nKnots
513 << " -> " << spline.getNumberOfKnots();
514 continue;
515 }
516 } while (0);
517
518 std::string err = FlatObject::stressTest(spline);
519 if (!err.empty()) {
520 LOG(info) << "error at FlatObject functionality: " << err;
521 return -1;
522 } else {
523 // LOG(info) << "flat object functionality is ok" ;
524 }
525
526 // Ndim-D spline
527
528 auto startTime = std::chrono::high_resolution_clock::now();
529 spline.approximateFunctionViaDataPoints(0., uMax, 0., uMax, F, nAuxiliaryPoints, nAuxiliaryPoints);
530 auto stopTime = std::chrono::high_resolution_clock::now();
531 statTime += std::chrono::duration_cast<std::chrono::nanoseconds>(stopTime - startTime);
532
533 // if (itry == 0)
534 if (0) {
535 TFile outf("testSpline2D.root", "recreate");
536 if (outf.IsZombie()) {
537 LOG(info) << "Failed to open output file testSpline2D.root ";
538 } else {
539 const char* name = "spline2Dtest";
540 spline.writeToFile(outf, name);
542 if (p == nullptr) {
543 LOG(info) << "Failed to read Spline1DOld from file testSpline1DOld.root ";
544 } else {
545 spline = *p;
546 }
547 outf.Close();
548 }
549 }
550
551 // 1-D splines for each of Ndim dimensions
552
553 Spline2D<DataT, 1> splines1D[Ndim];
554
555 for (int32_t dim = 0; dim < Ndim; dim++) {
556 auto F1 = [&](double x1, double x2, double f[]) {
557 double ff[Ndim];
558 F(x1, x2, ff);
559 f[0] = ff[dim];
560 };
561 splines1D[dim].recreate(nKnots, knotsU, nKnots, knotsV);
562 splines1D[dim].approximateFunctionViaDataPoints(0., uMax, 0., uMax, F1, nAuxiliaryPoints, nAuxiliaryPoints);
563 }
564
565 double stepU = .1;
566 for (double u = 0; u < uMax; u += stepU) {
567 for (double v = 0; v < uMax; v += stepU) {
568 double f[Ndim];
569 F(u, v, f);
570 DataT s[Ndim];
571 spline.interpolate(u, v, s);
572 for (int32_t dim = 0; dim < Ndim; dim++) {
573 statDf += (s[dim] - f[dim]) * (s[dim] - f[dim]);
574 DataT s1 = splines1D[dim].interpolate(u, v);
575 statDf1D += (s[dim] - s1) * (s[dim] - s1);
576 }
577 statN += Ndim;
578 // LOG(info) << u << " " << v << ": f " << f << " s " << s << " df "
579 // << s - f << " " << sqrt(statDf / statN) ;
580 }
581 }
582 // LOG(info) << "Spline2D standard deviation : " << sqrt(statDf / statN)
583 // ;
584
585 if (draw) {
586 delete nt;
587 delete knots;
588 nt = new TNtuple("nt", "nt", "u:v:f:s");
589 knots = new TNtuple("knots", "knots", "type:u:v:s");
590 double stepU = .3;
591 for (double u = 0; u < uMax; u += stepU) {
592 for (double v = 0; v < uMax; v += stepU) {
593 double f[Ndim];
594 F(u, v, f);
595 DataT s[Ndim];
596 spline.interpolate(u, v, s);
597 nt->Fill(u, v, f[0], s[0]);
598 }
599 }
600 nt->SetMarkerStyle(8);
601
602 nt->SetMarkerSize(.5);
603 nt->SetMarkerColor(kBlue);
604 nt->Draw("s:u:v", "", "");
605
606 nt->SetMarkerColor(kGray);
607 nt->SetMarkerSize(2.);
608 nt->Draw("f:u:v", "", "same");
609
610 nt->SetMarkerSize(.5);
611 nt->SetMarkerColor(kBlue);
612 nt->Draw("s:u:v", "", "same");
613
614 for (int32_t i = 0; i < nKnots; i++) {
615 for (int32_t j = 0; j < nKnots; j++) {
616 double u = spline.getGridX1().getKnot(i).u;
617 double v = spline.getGridX2().getKnot(j).u;
618 DataT s[Ndim];
619 spline.interpolate(u, v, s);
620 knots->Fill(1, u, v, s[0]);
621 }
622 }
623
624 knots->SetMarkerStyle(8);
625 knots->SetMarkerSize(1.5);
626 knots->SetMarkerColor(kRed);
627 knots->SetMarkerSize(1.5);
628 knots->Draw("s:u:v", "type==1", "same"); // knots
629
630 if (drawDataPoints) {
632 helper.setSpline(spline, 4, 4);
633 for (int32_t ipu = 0; ipu < helper.getHelperU1().getNumberOfDataPoints(); ipu++) {
634 const typename Spline1DHelperOld<DataT>::DataPoint& pu = helper.getHelperU1().getDataPoint(ipu);
635 for (int32_t ipv = 0; ipv < helper.getHelperU2().getNumberOfDataPoints(); ipv++) {
636 const typename Spline1DHelperOld<DataT>::DataPoint& pv = helper.getHelperU2().getDataPoint(ipv);
637 if (pu.isKnot && pv.isKnot) {
638 continue;
639 }
640 DataT s[Ndim];
641 spline.interpolate(pu.u, pv.u, s);
642 knots->Fill(2, pu.u, pv.u, s[0]);
643 }
644 }
645 knots->SetMarkerColor(kBlack);
646 knots->SetMarkerSize(1.);
647 knots->Draw("s:u:v", "type==2", "same"); // data points
648 }
649
650 if (!ask()) {
651 break;
652 }
653 }
654 }
655 // delete canv;
656 // delete nt;
657 // delete knots;
658
659 statDf = sqrt(statDf / statN);
660 statDf1D = sqrt(statDf1D / statN);
661
662 LOG(info) << "\n std dev for Spline2D : " << statDf;
663 LOG(info) << " mean difference between 1-D and " << Ndim
664 << "-D splines : " << statDf1D;
665 LOG(info) << " approximation time " << statTime.count() / 1000. / nTries << " ms";
666
667 if (statDf < 0.15 && statDf1D < 1.e-20) {
668 LOG(info) << "Everything is fine";
669 } else {
670 LOG(info) << "Something is wrong!!";
671 return -2;
672 }
673
674 return 0;
675}
676
677template class o2::gpu::Spline2DHelper<float>;
int32_t i
Definition of SymMatrixSolver class.
uint32_t j
Definition RawData.h:0
uint32_t c
Definition RawData.h:2
Definition of Spline1DHelper class.
Definition of Spline2DHelper class.
static std::string stressTest(T &obj)
Test the flat object functionality for a child class T.
Definition FlatObject.h:423
static void getScoefficients(const typename Spline1D< double >::KnotType &knotL, double u, double &cSl, double &cDl, double &cSr, double &cDr)
void recreate(int32_t nYdim, int32_t nKnotsX1, int32_t nKnotsX2)
Constructor for a regular spline.
const Spline1DHelperOld< DataT > & getHelperU2() const
void approximateFunctionViaDataPoints(Spline2DContainerBase< DataT, FlatObject > &spline, double x1Min, double x1Max, double x2Min, double x2Max, std::function< void(double x1, double x2, double f[])> F, int32_t nAuxiliaryDataPointsU1=4, int32_t nAuxiliaryDataPointsU2=4)
void approximateFunction(Spline2DContainerBase< DataT, FlatObject > &spline, double x1Min, double x1Max, double x2Min, double x2Max, std::function< void(double x1, double x2, double f[])> F, int32_t nAuxiliaryDataPointsU1=4, int32_t nAuxiliaryDataPointsU2=4)
_______________ Main functionality ________________________
int32_t setSpline(const Spline2DContainerBase< DataT, FlatObject > &spline, int32_t nAuxiliaryPointsU1, int32_t nAuxiliaryPointsU2)
_______________ Interface for a step-wise construction of the best-fit spline _______________________...
const Spline1DHelperOld< DataT > & getHelperU1() const
void approximateFunctionBatch(DataT *Fparameters, double x1Min, double x1Max, double x2Min, double x2Max, std::function< void(const std::vector< double > &x1, const std::vector< double > &x2, std::vector< double > f[])> F, uint32_t batchsize) const
approximate std::function, output in Fparameters. F calculates values for a batch of points.
void approximateDataPoints(Spline2DContainerBase< DataT, FlatObject > &spline, DataT *splineParameters, double x1Min, double x1Max, double x2Min, double x2Max, const double dataPointX1[], const double dataPointX2[], const double dataPointF[], const double dataPointWeight[], int32_t nDataPoints)
Create best-fit spline parameters for a given set of data points.
Spline2DHelper()
_____________ Constructors / destructors __________________________
static int32_t test(const bool draw=0, const bool drawDataPoints=1)
Test the Spline2D class functionality.
Forward declaration — specializations below select ClassDefNV based on FlatBase.
Definition Spline2D.h:104
double & B(int32_t i, int32_t j)
access to B elements
double & A(int32_t i, int32_t j)
access to A elements
const GLdouble * v
Definition glcorearb.h:832
GLuint GLfloat GLfloat GLfloat x1
Definition glcorearb.h:5034
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat s1
Definition glcorearb.h:5034
GLuint index
Definition glcorearb.h:781
GLuint const GLchar * name
Definition glcorearb.h:781
GLdouble f
Definition glcorearb.h:310
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLsizei GLenum const void * indices
Definition glcorearb.h:400
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
GLuint counter
Definition glcorearb.h:3987
TCanvas * drawDataPoints(TMultiGraph *mg, double min, double max)
Helper structure for 1D spline construction.
bool isKnot
is the point placed at a knot
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
const std::string str
uint64_t const void const *restrict const msg
Definition x9.h:153