Project
Loading...
Searching...
No Matches
Spline1DHelperOld.h
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
14
16
17#ifndef ALICEO2_GPUCOMMON_TPCFASTTRANSFORMATION_Spline1DHelperOld_H
18#define ALICEO2_GPUCOMMON_TPCFASTTRANSFORMATION_Spline1DHelperOld_H
19
20#include <cmath>
21#include <vector>
22
23#include "GPUCommonDef.h"
24#include "GPUCommonRtypes.h"
25#include "Spline1D.h"
26#include <functional>
27#include <string>
28
29namespace o2::gpu
30{
34template <typename DataT>
36{
37 public:
41 struct DataPoint {
42 double u;
43 double cS0;
44 double cZ0;
45 double cS1;
46 double cZ1;
47 int32_t iKnot;
48 bool isKnot;
49 };
50
52
55
58
61
63 ~Spline1DHelperOld() = default;
64
66
67 void bandGauss(double A[], double b[], int32_t n);
68
71 double xMin, double xMax,
72 double x[], double f[], int32_t nDataPoints);
73
76 double xMin, double xMax, std::function<void(double x, double f[/*spline.getFdimensions()*/])> F,
77 int32_t nAuxiliaryDataPoints = 4);
78
81 double xMin, double xMax, std::function<void(double x, double f[/*spline.getFdimensions()*/])> F,
82 int32_t nAuxiliaryDataPoints = 4);
83
86 double xMin, double xMax, std::function<void(double x, double f[/*spline.getFdimensions()*/])> F);
87
89
91 int32_t setSpline(const Spline1DContainer<DataT, FlatObject>& spline, int32_t nFdimensions, int32_t nAuxiliaryDataPoints);
92
94 int32_t setSpline(const Spline1DContainer<DataT, FlatObject>& spline, int32_t nFdimensions, double xMin, double xMax, double vx[], int32_t nDataPoints);
95
97 void approximateFunction(DataT* Fparameters, double xMin, double xMax, std::function<void(double x, double f[])> F) const;
98
100 void approximateFunctionGradually(DataT* Fparameters, double xMin, double xMax, std::function<void(double x, double f[])> F) const;
101
103 int32_t getNumberOfDataPoints() const { return mDataPoints.size(); }
104
106 void approximateFunction(DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim*/]) const;
107
109 void approximateFunctionGradually(DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim */]) const;
110
112 void copySfromDataPoints(DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim*/]) const;
113
116 void approximateDerivatives(DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim*/]) const;
117
119
120 const Spline1D<double>& getSpline() const { return mSpline; }
121
122 int32_t getKnotDataPoint(int32_t iknot) const { return mKnotDataPoints[iknot]; }
123
124 const DataPoint& getDataPoint(int32_t ip) const { return mDataPoints[ip]; }
125
128 static void getScoefficients(const typename Spline1D<double>::KnotType& knotL, double u,
129 double& cSl, double& cDl, double& cSr, double& cDr);
130
131 static void getDScoefficients(const typename Spline1D<double>::KnotType& knotL, double u,
132 double& cSl, double& cDl, double& cSr, double& cDr);
133
134 static void getDDScoefficients(const typename Spline1D<double>::KnotType& knotL, double u,
135 double& cSl, double& cDl, double& cSr, double& cDr);
136
137 static void getDDScoefficientsLeft(const typename Spline1D<double>::KnotType& knotL,
138 double& cSl, double& cDl, double& cSr, double& cDr);
139
140 static void getDDScoefficientsRight(const typename Spline1D<double>::KnotType& knotL,
141 double& cSl, double& cDl, double& cSr, double& cDr);
142 static void getDDDScoefficients(const typename Spline1D<double>::KnotType& knotL,
143 double& cSl, double& cDl, double& cSr, double& cDr);
144
146 const char* getLastError() const { return mError.c_str(); }
147
148#if !defined(GPUCA_GPUCODE) && !defined(GPUCA_STANDALONE) // code invisible on GPU and in the standalone compilation
150 static int32_t test(const bool draw = 0, const bool drawDataPoints = 1);
151#endif
152
153 private:
155 int32_t storeError(int32_t code, const char* msg);
156
157 std::string mError = "";
158
160
161 Spline1D<double> mSpline;
162 int32_t mFdimensions;
163 std::vector<DataPoint> mDataPoints;
164 std::vector<int32_t> mKnotDataPoints;
165 std::vector<double> mLSMmatrixFull;
166 std::vector<double> mLSMmatrixSderivatives;
167 std::vector<double> mLSMmatrixSvalues;
168
169 ClassDefNV(Spline1DHelperOld, 0);
170};
171
172} // namespace o2::gpu
173
174#endif
Definition of Spline1D class.
Definition A.h:16
static void getDDDScoefficients(const typename Spline1D< double >::KnotType &knotL, double &cSl, double &cDl, double &cSr, double &cDr)
void bandGauss(double A[], double b[], int32_t n)
_______________ Main functionality ________________________
const char * getLastError() const
Gives error string.
void approximateFunctionClassic(Spline1DContainer< DataT, FlatObject > &spline, double xMin, double xMax, std::function< void(double x, double f[])> F)
Create classic spline parameters for a given input function F.
~Spline1DHelperOld()=default
Destructor.
Spline1DHelperOld()
_____________ Constructors / destructors __________________________
static void getDDScoefficients(const typename Spline1D< double >::KnotType &knotL, double u, double &cSl, double &cDl, double &cSr, double &cDr)
void approximateFunctionGradually(Spline1DContainer< DataT, FlatObject > &spline, double xMin, double xMax, std::function< void(double x, double f[])> F, int32_t nAuxiliaryDataPoints=4)
Create best-fit spline parameters gradually for a given input function F.
static void getDScoefficients(const typename Spline1D< double >::KnotType &knotL, double u, double &cSl, double &cDl, double &cSr, double &cDr)
int32_t getKnotDataPoint(int32_t iknot) const
void copySfromDataPoints(DataT *Fparameters, const double DataPointF[]) const
a tool for the gradual approximation: set spline values S_i at knots == function values
static void getScoefficients(const typename Spline1D< double >::KnotType &knotL, double u, double &cSl, double &cDl, double &cSr, double &cDr)
void approximateDataPoints(Spline1DContainer< DataT, FlatObject > &spline, double xMin, double xMax, double x[], double f[], int32_t nDataPoints)
Create best-fit spline parameters for a given input function F.
const Spline1D< double > & getSpline() const
_______________ Utilities ________________________
static void getDDScoefficientsRight(const typename Spline1D< double >::KnotType &knotL, double &cSl, double &cDl, double &cSr, double &cDr)
Spline1DHelperOld & operator=(const Spline1DHelperOld &)=default
Assignment operator: disabled.
static void getDDScoefficientsLeft(const typename Spline1D< double >::KnotType &knotL, double &cSl, double &cDl, double &cSr, double &cDr)
void approximateDerivatives(DataT *Fparameters, const double DataPointF[]) const
void approximateFunction(Spline1DContainer< DataT, FlatObject > &spline, double xMin, double xMax, std::function< void(double x, double f[])> F, int32_t nAuxiliaryDataPoints=4)
Create best-fit spline parameters for a given input function F.
int32_t getNumberOfDataPoints() const
number of data points
const DataPoint & getDataPoint(int32_t ip) const
int32_t setSpline(const Spline1DContainer< DataT, FlatObject > &spline, int32_t nFdimensions, int32_t nAuxiliaryDataPoints)
_______________ Interface for a step-wise construction of the best-fit spline _______________________...
Spline1DHelperOld(const Spline1DHelperOld &)=default
Copy constructor: disabled.
Forward declaration — specializations below select ClassDefNV based on FlatBase.
Definition Spline1D.h:171
GLdouble n
Definition glcorearb.h:1982
GLint GLenum GLint x
Definition glcorearb.h:403
GLdouble f
Definition glcorearb.h:310
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
Helper structure for 1D spline construction.
bool isKnot
is the point placed at a knot
int32_t iKnot
index of the left knot of the segment
uint64_t const void const *restrict const msg
Definition x9.h:153