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
30{
31namespace gpu
32{
36template <typename DataT>
38{
39 public:
43 struct DataPoint {
44 double u;
45 double cS0;
46 double cZ0;
47 double cS1;
48 double cZ1;
49 int32_t iKnot;
50 bool isKnot;
51 };
52
54
57
60
63
65 ~Spline1DHelperOld() = default;
66
68
69 void bandGauss(double A[], double b[], int32_t n);
70
73 double xMin, double xMax,
74 double x[], double f[], int32_t nDataPoints);
75
78 double xMin, double xMax, std::function<void(double x, double f[/*spline.getFdimensions()*/])> F,
79 int32_t nAuxiliaryDataPoints = 4);
80
83 double xMin, double xMax, std::function<void(double x, double f[/*spline.getFdimensions()*/])> F,
84 int32_t nAuxiliaryDataPoints = 4);
85
88 double xMin, double xMax, std::function<void(double x, double f[/*spline.getFdimensions()*/])> F);
89
91
93 int32_t setSpline(const Spline1DContainer<DataT>& spline, int32_t nFdimensions, int32_t nAuxiliaryDataPoints);
94
96 int32_t setSpline(const Spline1DContainer<DataT>& spline, int32_t nFdimensions, double xMin, double xMax, double vx[], int32_t nDataPoints);
97
99 void approximateFunction(DataT* Fparameters, double xMin, double xMax, std::function<void(double x, double f[])> F) const;
100
102 void approximateFunctionGradually(DataT* Fparameters, double xMin, double xMax, std::function<void(double x, double f[])> F) const;
103
105 int32_t getNumberOfDataPoints() const { return mDataPoints.size(); }
106
108 void approximateFunction(DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim*/]) const;
109
111 void approximateFunctionGradually(DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim */]) const;
112
114 void copySfromDataPoints(DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim*/]) const;
115
118 void approximateDerivatives(DataT* Fparameters, const double DataPointF[/*getNumberOfDataPoints() x nFdim*/]) const;
119
121
122 const Spline1D<double>& getSpline() const { return mSpline; }
123
124 int32_t getKnotDataPoint(int32_t iknot) const { return mKnotDataPoints[iknot]; }
125
126 const DataPoint& getDataPoint(int32_t ip) const { return mDataPoints[ip]; }
127
130 static void getScoefficients(const typename Spline1D<double>::Knot& knotL, double u,
131 double& cSl, double& cDl, double& cSr, double& cDr);
132
133 static void getDScoefficients(const typename Spline1D<double>::Knot& knotL, double u,
134 double& cSl, double& cDl, double& cSr, double& cDr);
135
136 static void getDDScoefficients(const typename Spline1D<double>::Knot& knotL, double u,
137 double& cSl, double& cDl, double& cSr, double& cDr);
138
139 static void getDDScoefficientsLeft(const typename Spline1D<double>::Knot& knotL,
140 double& cSl, double& cDl, double& cSr, double& cDr);
141
142 static void getDDScoefficientsRight(const typename Spline1D<double>::Knot& knotL,
143 double& cSl, double& cDl, double& cSr, double& cDr);
144 static void getDDDScoefficients(const typename Spline1D<double>::Knot& knotL,
145 double& cSl, double& cDl, double& cSr, double& cDr);
146
148 const char* getLastError() const { return mError.c_str(); }
149
150#if !defined(GPUCA_GPUCODE) && !defined(GPUCA_STANDALONE) // code invisible on GPU and in the standalone compilation
152 static int32_t test(const bool draw = 0, const bool drawDataPoints = 1);
153#endif
154
155 private:
157 int32_t storeError(int32_t code, const char* msg);
158
159 std::string mError = "";
160
162
163 Spline1D<double> mSpline;
164 int32_t mFdimensions;
165 std::vector<DataPoint> mDataPoints;
166 std::vector<int32_t> mKnotDataPoints;
167 std::vector<double> mLSMmatrixFull;
168 std::vector<double> mLSMmatrixSderivatives;
169 std::vector<double> mLSMmatrixSvalues;
170
171 ClassDefNV(Spline1DHelperOld, 0);
172};
173
174} // namespace gpu
175} // namespace o2
176
177#endif
Definition of Spline1D class.
Definition A.h:16
static void getDDScoefficientsLeft(const typename Spline1D< double >::Knot &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.
~Spline1DHelperOld()=default
Destructor.
Spline1DHelperOld()
_____________ Constructors / destructors __________________________
void approximateFunction(Spline1DContainer< DataT > &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.
static void getScoefficients(const typename Spline1D< double >::Knot &knotL, double u, double &cSl, double &cDl, double &cSr, double &cDr)
static void getDScoefficients(const typename Spline1D< double >::Knot &knotL, double u, double &cSl, double &cDl, double &cSr, double &cDr)
int32_t getKnotDataPoint(int32_t iknot) const
static void getDDScoefficientsRight(const typename Spline1D< double >::Knot &knotL, double &cSl, double &cDl, double &cSr, double &cDr)
void approximateDataPoints(Spline1DContainer< DataT > &spline, double xMin, double xMax, double x[], double f[], int32_t nDataPoints)
Create best-fit spline parameters for a given input function F.
void approximateFunctionGradually(Spline1DContainer< DataT > &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.
void copySfromDataPoints(DataT *Fparameters, const double DataPointF[]) const
a tool for the gradual approximation: set spline values S_i at knots == function values
void approximateFunctionClassic(Spline1DContainer< DataT > &spline, double xMin, double xMax, std::function< void(double x, double f[])> F)
Create classic spline parameters for a given input function F.
static void getDDScoefficients(const typename Spline1D< double >::Knot &knotL, double u, double &cSl, double &cDl, double &cSr, double &cDr)
const Spline1D< double > & getSpline() const
_______________ Utilities ________________________
static void getDDDScoefficients(const typename Spline1D< double >::Knot &knotL, double &cSl, double &cDl, double &cSr, double &cDr)
int32_t setSpline(const Spline1DContainer< DataT > &spline, int32_t nFdimensions, int32_t nAuxiliaryDataPoints)
_______________ Interface for a step-wise construction of the best-fit spline _______________________...
Spline1DHelperOld & operator=(const Spline1DHelperOld &)=default
Assignment operator: disabled.
void approximateDerivatives(DataT *Fparameters, const double DataPointF[]) const
int32_t getNumberOfDataPoints() const
number of data points
const DataPoint & getDataPoint(int32_t ip) const
Spline1DHelperOld(const Spline1DHelperOld &)=default
Copy constructor: disabled.
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
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
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