votca 2024.2-dev
Loading...
Searching...
No Matches
votca::tools::CubicSpline Class Reference

A cubic spline class. More...

#include <cubicspline.h>

Inheritance diagram for votca::tools::CubicSpline:
Collaboration diagram for votca::tools::CubicSpline:

Public Member Functions

 CubicSpline ()=default
 
 ~CubicSpline () override=default
 
void Interpolate (const Eigen::VectorXd &x, const Eigen::VectorXd &y) override
 Calculate interpolating spline for given (x,y) values. Points on resulting spline can be obtained via Calculate().
 
void Fit (const Eigen::VectorXd &x, const Eigen::VectorXd &y) override
 Fit spline through noisy (x,y) values. Points on resulting fitted spline can be obtained via Calculate().
 
double Calculate (double r) override
 Calculate spline function value for a given x value on the spline created by Interpolate() or Fit()
 
double CalculateDerivative (double r) override
 Calculate y value for a given x value on the derivative of the spline created by function Interpolate or Fit.
 
void setSplineData (const Eigen::VectorXd &f, const Eigen::VectorXd &f2)
 
template<typename matrix_type >
void AddToFitMatrix (matrix_type &M, double x, Index offset1, Index offset2=0, double scale=1)
 Add a point (one entry) to fitting matrix.
 
template<typename matrix_type >
void AddToFitMatrix (matrix_type &M, double x, Index offset1, Index offset2, double scale1, double scale2)
 Add a point (one entry) to fitting matrix.
 
template<typename matrix_type , typename vector_type >
void AddToFitMatrix (matrix_type &M, vector_type &x, Index offset1, Index offset2=0)
 Add a vector of points to fitting matrix.
 
template<typename matrix_type >
void AddBCSumZeroToFitMatrix (matrix_type &M, Index offset1, Index offset2=0)
 Add boundary condition of sum_i f_i =0 to fitting matrix.
 
template<typename matrix_type >
void AddBCToFitMatrix (matrix_type &M, Index offset1, Index offset2=0)
 Add boundary conditions to fitting matrix.
 
Eigen::VectorXd Calculate (const Eigen::VectorXd &x)
 Calculate spline function values for given x values on the spline created by Interpolate() or Fit()
 
Eigen::VectorXd CalculateDerivative (const Eigen::VectorXd &x)
 Calculate y values for given x values on the derivative of the spline created by function Interpolate or Fit.
 
- Public Member Functions inherited from votca::tools::Spline
 Spline ()=default
 
virtual ~Spline ()=default
 
void setBC (eBoundary bc)
 Set the boundary type of the spline.
 
void setBCInt (Index bc)
 Set the boundary type of the spline.
 
double getGridPoint (int i)
 Get the grid point of certain index.
 
Eigen::VectorXd Calculate (const Eigen::VectorXd &x)
 Calculate spline function values for given x values on the spline created by Interpolate() or Fit()
 
Eigen::VectorXd CalculateDerivative (const Eigen::VectorXd &x)
 Calculate y values for given x values on the derivative of the spline created by function Interpolate or Fit.
 
void Print (std::ostream &out, double interval)
 Print spline values (using Calculate()) on output "out" on the entire grid in steps of "interval".
 
Index getInterval (double r)
 Determine the index of the interval containing value r.
 
Index GenerateGrid (double min, double max, double h)
 Generate the grid for fitting from "min" to "max" in steps of "h".
 
Eigen::VectorXd & getX ()
 Get the grid array x.
 
const Eigen::VectorXd & getX () const
 

Private Member Functions

double A (double r)
 
double B (double r)
 
double C (double r)
 
double D (double r)
 
double Aprime (double r)
 
double Bprime (double r)
 
double Cprime (double r)
 
double Dprime (double r)
 
double A_prime_l (Index i)
 
double A_prime_r (Index i)
 
double B_prime_l (Index i)
 
double B_prime_r (Index i)
 
double C_prime_l (Index i)
 
double C_prime_r (Index i)
 
double D_prime_l (Index i)
 
double D_prime_r (Index i)
 

Private Attributes

Eigen::VectorXd f_
 
Eigen::VectorXd f2_
 

Additional Inherited Members

- Public Types inherited from votca::tools::Spline
enum  eBoundary { splineNormal = 0 , splinePeriodic = 1 , splineDerivativeZero = 2 }
 enum for type of boundary condition More...
 
- Protected Attributes inherited from votca::tools::Spline
eBoundary boundaries_ = eBoundary::splineNormal
 Get the spline data f_.
 
Eigen::VectorXd r_
 

Detailed Description

A cubic spline class.

This class does cubic piecewise spline interpolation and spline fitting. As representation of a single spline, the general form

\[ S_i(x) = A(x,h_i) f_i + B(x,h_i) f_{i+1} + C(x,h_i) f''_i + d(x,h_i) f''_{i+1} \]

with

\[ x_i \le x < x_{i+1}\,,\\ h_i = x_{i+1} - x_{i} \]

The \(f_i\,,\,,f'' i_\) are the function values and second derivates at point \(x_i\).

The parameters \(f''_i\) are no free parameters, they are determined by the smoothing condition that the first derivatives are continuous. So the only free paramers are the grid points x_i as well as the functon values f_i at these points. A spline can be generated in several ways:

  • Interpolation spline
  • Fitting spline (fit to noisy data)
  • calculate the parameters elsewere and fill the spline class

Definition at line 54 of file cubicspline.h.

Constructor & Destructor Documentation

◆ CubicSpline()

votca::tools::CubicSpline::CubicSpline ( )
default

◆ ~CubicSpline()

votca::tools::CubicSpline::~CubicSpline ( )
overridedefault

Member Function Documentation

◆ A()

double votca::tools::CubicSpline::A ( double r)
private

Definition at line 141 of file cubicspline.cc.

◆ A_prime_l()

double votca::tools::CubicSpline::A_prime_l ( Index i)
private

Definition at line 190 of file cubicspline.cc.

◆ A_prime_r()

double votca::tools::CubicSpline::A_prime_r ( Index i)
private

Definition at line 202 of file cubicspline.cc.

◆ AddBCSumZeroToFitMatrix()

template<typename matrix_type >
void votca::tools::CubicSpline::AddBCSumZeroToFitMatrix ( matrix_type & M,
Index offset1,
Index offset2 = 0 )
inline

Add boundary condition of sum_i f_i =0 to fitting matrix.

Parameters
Mpointer to matrix
offset1relative to getInterval(x)
offset2relative to getInterval(x)

Definition at line 212 of file cubicspline.h.

◆ AddBCToFitMatrix()

template<typename matrix_type >
void votca::tools::CubicSpline::AddBCToFitMatrix ( matrix_type & M,
Index offset1,
Index offset2 = 0 )
inline

Add boundary conditions to fitting matrix.

Parameters
Mpointer to matrix
offset1relative to getInterval(x)
offset2relative to getInterval(x)

Definition at line 221 of file cubicspline.h.

◆ AddToFitMatrix() [1/3]

template<typename matrix_type >
void votca::tools::CubicSpline::AddToFitMatrix ( matrix_type & M,
double x,
Index offset1,
Index offset2,
double scale1,
double scale2 )
inline

Add a point (one entry) to fitting matrix.

Parameters
Mpointer to matrix [in] [out]
xvalue [in]
offset1relative to getInterval(x) [in]
offset2relative to getInterval(x) [in]
scale1parameters for terms "A,B,C,D" [in]
scale2parameters for terms "AA,BB,CC,DD" [in] When creating a matrix to fit data with a spline, this function creates one entry in that fitting matrix.

Definition at line 187 of file cubicspline.h.

◆ AddToFitMatrix() [2/3]

template<typename matrix_type >
void votca::tools::CubicSpline::AddToFitMatrix ( matrix_type & M,
double x,
Index offset1,
Index offset2 = 0,
double scale = 1 )
inline

Add a point (one entry) to fitting matrix.

Parameters
Mpointer to matrix
xvalue
offset1relative to getInterval(x)
offset2relative to getInterval(x)
scaleparameters for terms "A,B,C,D" When creating a matrix to fit data with a spline, this function creates one entry in that fitting matrix.

Definition at line 176 of file cubicspline.h.

◆ AddToFitMatrix() [3/3]

template<typename matrix_type , typename vector_type >
void votca::tools::CubicSpline::AddToFitMatrix ( matrix_type & M,
vector_type & x,
Index offset1,
Index offset2 = 0 )
inline

Add a vector of points to fitting matrix.

Parameters
Mpointer to matrix
xvector of x values
offset1relative to getInterval(x)
offset2relative to getInterval(x) Same as previous function, but vector-valued and with scale=1.0

Definition at line 200 of file cubicspline.h.

◆ Aprime()

double votca::tools::CubicSpline::Aprime ( double r)
private

Definition at line 146 of file cubicspline.cc.

◆ B()

double votca::tools::CubicSpline::B ( double r)
private

Definition at line 150 of file cubicspline.cc.

◆ B_prime_l()

double votca::tools::CubicSpline::B_prime_l ( Index i)
private

Definition at line 192 of file cubicspline.cc.

◆ B_prime_r()

double votca::tools::CubicSpline::B_prime_r ( Index i)
private

Definition at line 206 of file cubicspline.cc.

◆ Bprime()

double votca::tools::CubicSpline::Bprime ( double r)
private

Definition at line 155 of file cubicspline.cc.

◆ C()

double votca::tools::CubicSpline::C ( double r)
private

Definition at line 159 of file cubicspline.cc.

◆ C_prime_l()

double votca::tools::CubicSpline::C_prime_l ( Index i)
private

Definition at line 194 of file cubicspline.cc.

◆ C_prime_r()

double votca::tools::CubicSpline::C_prime_r ( Index i)
private

Definition at line 208 of file cubicspline.cc.

◆ Calculate() [1/2]

Eigen::VectorXd votca::tools::Spline::Calculate ( const Eigen::VectorXd & x)

Calculate spline function values for given x values on the spline created by Interpolate() or Fit()

Parameters
xvector of data values
Returns
vector of y value

Definition at line 116 of file spline.cc.

◆ Calculate() [2/2]

double votca::tools::CubicSpline::Calculate ( double x)
overridevirtual

Calculate spline function value for a given x value on the spline created by Interpolate() or Fit()

Parameters
xdata value
Returns
y value

Implements votca::tools::Spline.

Definition at line 129 of file cubicspline.cc.

◆ CalculateDerivative() [1/2]

Eigen::VectorXd votca::tools::Spline::CalculateDerivative ( const Eigen::VectorXd & x)

Calculate y values for given x values on the derivative of the spline created by function Interpolate or Fit.

Parameters
xvector of data values
Returns
vector of y value

Definition at line 124 of file spline.cc.

◆ CalculateDerivative() [2/2]

double votca::tools::CubicSpline::CalculateDerivative ( double x)
overridevirtual

Calculate y value for a given x value on the derivative of the spline created by function Interpolate or Fit.

Parameters
xdata value
Returns
y value of derivative

Implements votca::tools::Spline.

Definition at line 135 of file cubicspline.cc.

◆ Cprime()

double votca::tools::CubicSpline::Cprime ( double r)
private

Definition at line 168 of file cubicspline.cc.

◆ D()

double votca::tools::CubicSpline::D ( double r)
private

Definition at line 175 of file cubicspline.cc.

◆ D_prime_l()

double votca::tools::CubicSpline::D_prime_l ( Index i)
private

Definition at line 198 of file cubicspline.cc.

◆ D_prime_r()

double votca::tools::CubicSpline::D_prime_r ( Index i)
private

Definition at line 212 of file cubicspline.cc.

◆ Dprime()

double votca::tools::CubicSpline::Dprime ( double r)
private

Definition at line 183 of file cubicspline.cc.

◆ Fit()

void votca::tools::CubicSpline::Fit ( const Eigen::VectorXd & x,
const Eigen::VectorXd & y )
overridevirtual

Fit spline through noisy (x,y) values. Points on resulting fitted spline can be obtained via Calculate().

Parameters
xvalues of data to be fitted
yvalues of data to be fitted both vectors must be of same size

Implements votca::tools::Spline.

Definition at line 86 of file cubicspline.cc.

◆ Interpolate()

void votca::tools::CubicSpline::Interpolate ( const Eigen::VectorXd & x,
const Eigen::VectorXd & y )
overridevirtual

Calculate interpolating spline for given (x,y) values. Points on resulting spline can be obtained via Calculate().

Parameters
xvalues of data to be interpolated
yvalues of data to be interpolated both vectors must be of same size

Implements votca::tools::Spline.

Definition at line 31 of file cubicspline.cc.

◆ setSplineData()

void votca::tools::CubicSpline::setSplineData ( const Eigen::VectorXd & f,
const Eigen::VectorXd & f2 )
inline

Definition at line 81 of file cubicspline.h.

Member Data Documentation

◆ f2_

Eigen::VectorXd votca::tools::CubicSpline::f2_
private

Definition at line 150 of file cubicspline.h.

◆ f_

Eigen::VectorXd votca::tools::CubicSpline::f_
private

Definition at line 148 of file cubicspline.h.


The documentation for this class was generated from the following files: