SG++
sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT > Class Template Reference

B-spline basis on Clenshaw-Curtis grids. More...

#include <BsplineModifiedClenshawCurtisBasis.hpp>

Inheritance diagram for sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >:
sgpp::base::Basis< LT, IT >

Public Member Functions

 BsplineModifiedClenshawCurtisBasis ()
 Default constructor. More...
 
 BsplineModifiedClenshawCurtisBasis (size_t degree)
 Constructor. More...
 
double clenshawCurtisPoint (LT l, IT i) const
 
double clenshawCurtisPointNegativeIndex (LT l, IT ni) const
 
double eval (LT l, IT i, double x) override
 
double evalDx (LT l, IT i, double x)
 
double evalDxDx (LT l, IT i, double x)
 
size_t getDegree () const
 
 ~BsplineModifiedClenshawCurtisBasis () override
 Destructor. More...
 
- Public Member Functions inherited from sgpp::base::Basis< LT, IT >
virtual ~Basis ()
 Destructor. More...
 

Protected Member Functions

double clenshawCurtisPoint (LT l, IT i, IT hInv) const
 
double clenshawCurtisPointNegativeIndex (LT l, IT ni, IT hInv) const
 
void constructKnots (LT l, IT i, IT hInv)
 Construct the (p+2) Clenshaw-Curtis knots of a B-spline basis function and save them in xi. More...
 
void constructKnotsNegativeIndex (LT l, IT ni, IT hInv)
 Construct the (p+2) Clenshaw-Curtis knots of a B-spline basis function and save them in xi. More...
 
double modifiedBSpline (LT l, IT hInv, double x, size_t p)
 
double modifiedBSplineDx (LT l, IT hInv, double x, size_t p)
 
double modifiedBSplineDxDx (LT l, IT hInv, double x, size_t p)
 
double nonUniformBSpline (double x, size_t p, size_t k) const
 
double nonUniformBSplineDx (double x, size_t p, size_t k) const
 
double nonUniformBSplineDxDx (double x, size_t p, size_t k) const
 

Protected Attributes

ClenshawCurtisTableclenshawCurtisTable
 reference to the Clenshaw-Curtis cache table More...
 
size_t degree
 degree of the B-spline More...
 
std::vector< double > xi
 temporary helper vector of fixed size p+2 containing B-spline knots More...
 

Detailed Description

template<class LT, class IT>
class sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >

B-spline basis on Clenshaw-Curtis grids.

Constructor & Destructor Documentation

template<class LT, class IT>
sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::BsplineModifiedClenshawCurtisBasis ( )
inline

Default constructor.

template<class LT, class IT>
sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::BsplineModifiedClenshawCurtisBasis ( size_t  degree)
inlineexplicit

Constructor.

Parameters
degreeB-spline degree, must be odd (if it's even, degree - 1 is used)
template<class LT, class IT>
sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::~BsplineModifiedClenshawCurtisBasis ( )
inlineoverride

Destructor.

Member Function Documentation

template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::clenshawCurtisPoint ( LT  l,
IT  i 
) const
inline
template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::clenshawCurtisPoint ( LT  l,
IT  i,
IT  hInv 
) const
inlineprotected
Parameters
llevel of the grid point
iindex of the grid point
hInv2^l
Returns
i-th Clenshaw-Curtis grid point with level l
template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::clenshawCurtisPointNegativeIndex ( LT  l,
IT  ni 
) const
inline
template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::clenshawCurtisPointNegativeIndex ( LT  l,
IT  ni,
IT  hInv 
) const
inlineprotected
Parameters
llevel of the grid point
ninegative index -i of the grid point
hInv2^l
Returns
(-ni)-th Clenshaw-Curtis grid point with level l
template<class LT, class IT>
void sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::constructKnotsNegativeIndex ( LT  l,
IT  ni,
IT  hInv 
)
inlineprotected

Construct the (p+2) Clenshaw-Curtis knots of a B-spline basis function and save them in xi.

Parameters
llevel of basis function
ninegative index -i of basis function
hInv2^l

Referenced by sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::modifiedBSpline(), sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::modifiedBSplineDx(), and sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::modifiedBSplineDxDx().

template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::evalDx ( LT  l,
IT  i,
double  x 
)
inline
Parameters
llevel of basis function
iindex of basis function
xevaluation point
Returns
value of derivative of modified Clenshaw-Curtis B-spline basis function

Referenced by sgpp::base::OperationEvalGradientModBsplineClenshawCurtisNaive::evalGradient(), sgpp::base::OperationEvalHessianModBsplineClenshawCurtisNaive::evalHessian(), and sgpp::base::OperationEvalPartialDerivativeModBsplineClenshawCurtisNaive::evalPartialDerivative().

template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::evalDxDx ( LT  l,
IT  i,
double  x 
)
inline
Parameters
llevel of basis function
iindex of basis function
xevaluation point
Returns
value of 2nd derivative of modified Clenshaw-Curtis B-spline basis function

Referenced by sgpp::base::OperationEvalHessianModBsplineClenshawCurtisNaive::evalHessian().

template<class LT, class IT>
size_t sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::getDegree ( ) const
inline
Returns
B-spline degree
template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::modifiedBSpline ( LT  l,
IT  hInv,
double  x,
size_t  p 
)
inlineprotected
Parameters
llevel of basis function
hInv2^l
xevaluation point
pB-spline degree
Returns
value of modified Clenshaw-Curtis B-spline (e.g. index == 1)

Referenced by sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::eval().

template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::modifiedBSplineDx ( LT  l,
IT  hInv,
double  x,
size_t  p 
)
inlineprotected
Parameters
llevel of basis function
hInv2^l
xevaluation point
pB-spline degree
Returns
value of derivative of modified Clenshaw-Curtis B-spline (e.g. index == 1)

Referenced by sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::evalDx().

template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::modifiedBSplineDxDx ( LT  l,
IT  hInv,
double  x,
size_t  p 
)
inlineprotected
Parameters
llevel of basis function
hInv2^l
xevaluation point
pB-spline degree
Returns
value of 2nd derivative of modified Clenshaw-Curtis B-spline (e.g. index == 1)

Referenced by sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::evalDxDx().

template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::nonUniformBSplineDx ( double  x,
size_t  p,
size_t  k 
) const
inlineprotected
Parameters
xevaluation point
pB-spline degree
kindex of B-spline in the knot sequence
Returns
value of derivative of non-uniform B-spline with knots \(\{\xi_k, ... \xi_{k+p+1}\}\)

Referenced by sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::evalDx(), and sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::modifiedBSplineDx().

template<class LT, class IT>
double sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::nonUniformBSplineDxDx ( double  x,
size_t  p,
size_t  k 
) const
inlineprotected
Parameters
xevaluation point
pB-spline degree
kindex of B-spline in the knot sequence
Returns
value of 2nd derivative of non-uniform B-spline with knots \(\{\xi_k, ... \xi_{k+p+1}\}\)

Referenced by sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::evalDxDx(), and sgpp::base::BsplineModifiedClenshawCurtisBasis< unsigned int, unsigned int >::modifiedBSplineDxDx().

Member Data Documentation

template<class LT, class IT>
ClenshawCurtisTable& sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::clenshawCurtisTable
protected

reference to the Clenshaw-Curtis cache table

template<class LT, class IT>
std::vector<double> sgpp::base::BsplineModifiedClenshawCurtisBasis< LT, IT >::xi
protected

temporary helper vector of fixed size p+2 containing B-spline knots


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