npstat is hosted by Hepforge, IPPP Durham
NPStat  5.10.0
npstat::LegendreOrthoPoly1D Class Reference

#include <ClassicalOrthoPolys1D.hh>

Inheritance diagram for npstat::LegendreOrthoPoly1D:
npstat::AbsClassicalOrthoPoly1D

Public Member Functions

virtual LegendreOrthoPoly1Dclone () const
 
long double weight (const long double x) const
 
double xmin () const
 
double xmax () const
 
void integralCoeffs (const double *coeffs, unsigned degree, double *integCoeffs) const
 
void derivativeCoeffs (const double *coeffs, unsigned degree, double *derivCoeffs) const
 
void monomialCoeffs (const double *coeffs, unsigned degree, long double *monoCoeffs) const
 
void squaredPolyCoeffs (unsigned degree, long double *expansionCoeffs) const
 
virtual std::pair< long double, long double > monicRecurrenceCoeffs (const unsigned k) const
 
virtual std::pair< long double, long double > recurrenceCoeffs (const unsigned k) const
 
- Public Member Functions inherited from npstat::AbsClassicalOrthoPoly1D
virtual unsigned maxDegree () const
 
Matrix< long double > jacobiMatrix (unsigned n) const
 
long double poly (unsigned deg, long double x) const
 
long double monic (unsigned deg, long double x) const
 
std::pair< long double, long double > twopoly (unsigned deg1, unsigned deg2, long double x) const
 
void allpoly (long double x, long double *values, unsigned maxdeg) const
 
void allmonic (long double x, long double *values, unsigned maxdeg) const
 
double series (const double *coeffs, unsigned maxdeg, double x) const
 
double integrateSeries (const double *coeffs, unsigned maxdeg) const
 
template<class Functor , class Quadrature >
void calculateCoeffs (const Functor &fcn, const Quadrature &quad, double *coeffs, unsigned maxdeg) const
 
template<class Numeric >
void sampleCoeffs (const Numeric *coords, unsigned long lenCoords, double *coeffs, unsigned maxdeg) const
 
template<class Numeric >
void sampleCoeffVars (const Numeric *coords, unsigned long lenCoords, const double *coeffs, unsigned maxdeg, double *variances) const
 
template<class Numeric >
Matrix< double > sampleCoeffCovariance (const Numeric *coords, unsigned long lenCoords, const double *coeffs, unsigned maxdeg) const
 
template<class Numeric >
void sampleCoeffs (const Numeric *coords, unsigned long lenCoords, Numeric location, Numeric scale, double *coeffs, unsigned maxdeg) const
 
template<class Numeric >
void sampleCoeffVars (const Numeric *coords, unsigned long lenCoords, Numeric location, Numeric scale, const double *coeffs, unsigned maxdeg, double *variances) const
 
template<class Numeric >
Matrix< double > sampleCoeffCovariance (const Numeric *coords, unsigned long lenCoords, Numeric location, Numeric scale, const double *coeffs, unsigned maxdeg) const
 
template<class Pair >
void weightedSampleCoeffs (const Pair *points, unsigned long numPoints, double *coeffs, unsigned maxdeg) const
 
template<class Pair >
void weightedSampleCoeffVars (const Pair *points, unsigned long nPoints, const double *coeffs, unsigned maxdeg, double *variances) const
 
template<class Pair >
Matrix< double > weightedSampleCoeffCovariance (const Pair *points, unsigned long nPoints, const double *coeffs, unsigned maxdeg) const
 
template<class Pair >
void weightedSampleCoeffs (const Pair *points, unsigned long numPoints, typename Pair::first_type location, typename Pair::first_type scale, double *coeffs, unsigned maxdeg) const
 
template<class Pair >
void weightedSampleCoeffVars (const Pair *points, unsigned long nPoints, typename Pair::first_type location, typename Pair::first_type scale, const double *coeffs, unsigned maxdeg, double *variances) const
 
template<class Pair >
Matrix< double > weightedSampleCoeffCovariance (const Pair *points, unsigned long nPoints, typename Pair::first_type location, typename Pair::first_type scale, const double *coeffs, unsigned maxdeg) const
 
template<class Quadrature >
double empiricalKroneckerDelta (const Quadrature &quad, unsigned deg1, unsigned deg2) const
 
Matrix< double > empiricalKroneckerMatrix (const AbsIntervalQuadrature1D &quad, unsigned maxdeg) const
 
template<class Quadrature >
double jointAverage (const Quadrature &q, unsigned deg1, unsigned deg2, unsigned deg3, unsigned deg4) const
 
template<class Quadrature >
double jointAverage (const Quadrature &q, const unsigned *degrees, unsigned nDegrees, bool checkedForZeros=false) const
 
template<class Quadrature >
double directQuadrature (const Quadrature &q, unsigned deg1, unsigned deg2, unsigned deg3, unsigned deg4) const
 
template<class Quadrature >
double directQuadrature (const Quadrature &q, const unsigned *degrees, unsigned nDegrees, bool checkedForZeros=false) const
 
double integratePoly (unsigned degree, unsigned power) const
 
double jointIntegral (const unsigned *degrees, unsigned nDegrees) const
 
template<typename Numeric , typename Real >
void sampleAverages (const Numeric *coords, unsigned long lenCoords, Real *averages, unsigned maxdeg) const
 
template<class Pair , typename Real >
void weightedPointsAverages (const Pair *points, unsigned long numPoints, Real *averages, unsigned maxdeg) const
 
template<class Functor , typename Real >
void extWeightAverages (const Functor &extWeight, const AbsIntervalQuadrature1D &quad, Real *averages, unsigned maxdeg) const
 
template<class Numeric >
Matrix< double > sampleProductAverages (const Numeric *coords, unsigned long lenCoords, unsigned maxdeg) const
 
template<class Functor >
Matrix< double > extWeightProductAverages (const Functor &extWeight, const AbsIntervalQuadrature1D &quad, unsigned maxdeg) const
 
template<class Numeric >
double jointSampleAverage (const Numeric *coords, unsigned long lenCoords, const unsigned *degrees, unsigned lenDegrees) const
 
template<class Numeric >
std::pair< double, double > twoJointSampleAverages (const Numeric *coords, unsigned long lenCoords, const unsigned *degrees1, unsigned lenDegrees1, const unsigned *degrees2, unsigned lenDegrees2) const
 
CPP11_auto_ptr< StorablePolySeries1DmakeStorablePolySeries (unsigned maxPolyDeg, const double *coeffs=0, unsigned maxdeg=0, double shift=0.0) const
 
long double location () const
 
long double scale () const
 

Additional Inherited Members

- Public Types inherited from npstat::AbsClassicalOrthoPoly1D
typedef OrthoPoly1DDeg degree_functor
 
typedef OrthoPoly1DWeight weight_functor
 
- Static Protected Member Functions inherited from npstat::AbsClassicalOrthoPoly1D
static int kdelta (const unsigned i, const unsigned j)
 

Detailed Description

Orthonormal Legendre polynomials on [-1, 1]

Member Function Documentation

◆ clone()

virtual LegendreOrthoPoly1D* npstat::LegendreOrthoPoly1D::clone ( ) const
inlinevirtual

Virtual copy constructor

Implements npstat::AbsClassicalOrthoPoly1D.

◆ derivativeCoeffs()

void npstat::LegendreOrthoPoly1D::derivativeCoeffs ( const double *  coeffs,
unsigned  degree,
double *  derivCoeffs 
) const

Derive the series coefficients for the derivative of the argument series. The array of coefficients must be at least degree+1 long, and the buffer for the derivative coefficients must be at least degree long.

◆ integralCoeffs()

void npstat::LegendreOrthoPoly1D::integralCoeffs ( const double *  coeffs,
unsigned  degree,
double *  integCoeffs 
) const

Derive the series coefficients for the integral of the argument series from -1 to x. The array of coefficients must be at least degree+1 long, and the buffer for the integral coefficients must be at least degree+2 long.

◆ monomialCoeffs()

void npstat::LegendreOrthoPoly1D::monomialCoeffs ( const double *  coeffs,
unsigned  degree,
long double *  monoCoeffs 
) const

Derive the series coefficients for the monomial expansion. The arrays of coefficients must be at least degree+1 long.

◆ squaredPolyCoeffs()

void npstat::LegendreOrthoPoly1D::squaredPolyCoeffs ( unsigned  degree,
long double *  expansionCoeffs 
) const

Expand a squared polynomial into polynomial series. The number of coefficients will be 2*degree + 1.

◆ weight()

long double npstat::LegendreOrthoPoly1D::weight ( const long double  x) const
inlinevirtual

The weight function should be normalized

Implements npstat::AbsClassicalOrthoPoly1D.

◆ xmin()

double npstat::LegendreOrthoPoly1D::xmin ( ) const
inlinevirtual

Support of the polynomial

Implements npstat::AbsClassicalOrthoPoly1D.


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