npstat is hosted by Hepforge, IPPP Durham
NPStat  5.10.0

Public Member Functions

 OSDE1D (const AbsClassicalOrthoPoly1D &poly1d, double xmin, double xmax)
 
 OSDE1D (const OSDE1D &)
 
OSDE1Doperator= (const OSDE1D &)
 
const AbsClassicalOrthoPoly1Dclpoly () const
 
double xmin () const
 
double xmax () const
 
double weight (double x) const
 
double poly (const unsigned deg, const double x) const
 
std::pair< double, double > twopoly (unsigned deg1, unsigned deg2, double x) const
 
double series (const double *coeffs, unsigned maxdeg, double x) const
 
double integrateSeries (const double *coeffs, const unsigned maxdeg) const
 
double integratePoly (const unsigned degree, const unsigned power) const
 
double jointIntegral (const unsigned *degrees, unsigned nDegrees) const
 
Matrix< double > tMatrix (unsigned maxdeg) const
 
template<class Functor , class Quadrature >
void calculateCoeffs (const Functor &fcn, const Quadrature &quad, double *coeffs, unsigned maxdeg) const
 
void densityCoeffs (const AbsDistribution1D &distro, unsigned nIntegrationPoints, double *coeffs, unsigned maxdeg) const
 
void densityCoeffsVariance (const AbsDistribution1D &distro, unsigned nIntegrationPoints, const double *coeffs, unsigned maxdeg, double expectedSampleSize, double *variances) const
 
npstat::Matrix< double > densityCoeffsCovariance (const AbsDistribution1D &distro, unsigned nIntegrationPoints, const double *coeffs, unsigned maxdeg, double expectedSampleSize) const
 
template<class Numeric >
void sampleCoeffs (const Numeric *coords, unsigned long lenCoords, double *coeffs, unsigned maxdeg) const
 
template<class Numeric >
void sampleCoeffsVariance (const Numeric *coords, unsigned long lenCoords, const double *coeffs, unsigned maxdeg, double *variances) const
 
template<class Numeric >
npstat::Matrix< double > sampleCoeffsCovariance (const Numeric *coords, unsigned long lenCoords, 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 weightedSampleCoeffsVariance (const Pair *points, unsigned long nPoints, const double *coeffs, unsigned maxdeg, double *variances) const
 
template<class Pair >
npstat::Matrix< double > weightedSampleCoeffsCovariance (const Pair *points, unsigned long nPoints, const double *coeffs, unsigned maxdeg) const
 
template<class Functor , class Quadrature >
double integratedSquaredError (const Functor &fcn, const Quadrature &quad, const double *coeffs, unsigned maxdeg) const
 
template<class Functor , class Quadrature >
double integratedSquaredError (const Functor &fcn, const Quadrature &quad, const double *coeffs, unsigned maxdeg, double xmin, double xmax) const
 

Static Public Member Functions

static std::pair< double, double > supportRegion (const AbsClassicalOrthoPoly1D &poly1d, double leftLimit, bool leftIsFromSample, double rightLimit, bool rightIsFromSample, unsigned long nPoints)
 
static unsigned optimalDegreeHart (const double *coeffs, const double *variances, unsigned maxdeg, double k=2.0)
 

Constructor & Destructor Documentation

◆ OSDE1D()

npstat::OSDE1D::OSDE1D ( const AbsClassicalOrthoPoly1D poly1d,
double  xmin,
double  xmax 
)

Main constructor. The arguments are as follows:

poly1d – classical orthogonal polynomial system to use

xmin – lower bound of the density support region

xmax – upper bound of the density support region

Note that the density estimate will be automatically normalized for Legendre polynomials only.

Member Function Documentation

◆ calculateCoeffs()

template<class Functor , class Quadrature >
void npstat::OSDE1D::calculateCoeffs ( const Functor &  fcn,
const Quadrature &  quad,
double *  coeffs,
unsigned  maxdeg 
) const

Build the coefficients of the orthogonal polynomial series for the given function. The length of the array "coeffs" should be at least maxdeg + 1. Note that the coefficients are returned in the order of increasing degree (same order is used by the "series" function).

◆ clpoly()

const AbsClassicalOrthoPoly1D& npstat::OSDE1D::clpoly ( ) const
inline

Basic inspection of object properties

◆ densityCoeffs()

void npstat::OSDE1D::densityCoeffs ( const AbsDistribution1D distro,
unsigned  nIntegrationPoints,
double *  coeffs,
unsigned  maxdeg 
) const

Build the coefficients of the orthogonal polynomial series for the given density using Gauss-Legendre quadratures. The length of the array "coeffs" should be at least maxdeg + 1.

◆ densityCoeffsCovariance()

npstat::Matrix<double> npstat::OSDE1D::densityCoeffsCovariance ( const AbsDistribution1D distro,
unsigned  nIntegrationPoints,
const double *  coeffs,
unsigned  maxdeg,
double  expectedSampleSize 
) const

Expected covariances of the coefficients obtained previously with the "densityCoeffs" method.

◆ densityCoeffsVariance()

void npstat::OSDE1D::densityCoeffsVariance ( const AbsDistribution1D distro,
unsigned  nIntegrationPoints,
const double *  coeffs,
unsigned  maxdeg,
double  expectedSampleSize,
double *  variances 
) const

Expected variances of the coefficients obtained previously with the "densityCoeffs" method. The length of the array "variances" should be at least maxdeg + 1.

◆ integratedSquaredError() [1/2]

template<class Functor , class Quadrature >
double npstat::OSDE1D::integratedSquaredError ( const Functor &  fcn,
const Quadrature &  quad,
const double *  coeffs,
unsigned  maxdeg 
) const

Integrated squared error between the given function and the polynomial series

◆ integratedSquaredError() [2/2]

template<class Functor , class Quadrature >
double npstat::OSDE1D::integratedSquaredError ( const Functor &  fcn,
const Quadrature &  quad,
const double *  coeffs,
unsigned  maxdeg,
double  xmin,
double  xmax 
) const

Integrated squared error between the given function and the polynomial series on an arbitrary interval. If the given interval has regions beyond xmin() and xmax(), the integration will assume that on those regions the value of the series is 0, and the quadrature will be carried out there separately.

◆ integratePoly()

double npstat::OSDE1D::integratePoly ( const unsigned  degree,
const unsigned  power 
) const
inline

Unweighted polynomial integral

◆ integrateSeries()

double npstat::OSDE1D::integrateSeries ( const double *  coeffs,
const unsigned  maxdeg 
) const
inline

Integral of the series

◆ jointIntegral()

double npstat::OSDE1D::jointIntegral ( const unsigned *  degrees,
unsigned  nDegrees 
) const
inline

Unweighted integral of a product of multiple polynomials

◆ optimalDegreeHart()

static unsigned npstat::OSDE1D::optimalDegreeHart ( const double *  coeffs,
const double *  variances,
unsigned  maxdeg,
double  k = 2.0 
)
static

Optimal OSDE truncation degree according to the Hart scheme. Minimizes Sum_{i=0}^M (k*variances[i] - coeffs[i]^2) over M. The default value of k = 2 corresponds to the original scheme. Argument arrays "coeffs" and "variances" should have at least maxdeg + 1 elements and should be calculated in advance with "sampleCoeffs"/"sampleCoeffsVariance" or with "weightedSampleCoeffs"/"weightedSampleCoeffsVariance" methods.

◆ poly()

double npstat::OSDE1D::poly ( const unsigned  deg,
const double  x 
) const
inline

Polynomial values

◆ sampleCoeffs()

template<class Numeric >
void npstat::OSDE1D::sampleCoeffs ( const Numeric *  coords,
unsigned long  lenCoords,
double *  coeffs,
unsigned  maxdeg 
) const

Build the coefficients of the orthogonal polynomial series for the given sample of points (empirical density function). The length of the array "coeffs" should be at least maxdeg + 1. Note that the coefficients are returned in the order of increasing degree (same order is used by the "series" function).

◆ sampleCoeffsCovariance()

template<class Numeric >
npstat::Matrix<double> npstat::OSDE1D::sampleCoeffsCovariance ( const Numeric *  coords,
unsigned long  lenCoords,
const double *  coeffs,
unsigned  maxdeg 
) const

Estimate covariances of the coefficients obtained previously with the "sampleCoeffs" method.

◆ sampleCoeffsVariance()

template<class Numeric >
void npstat::OSDE1D::sampleCoeffsVariance ( const Numeric *  coords,
unsigned long  lenCoords,
const double *  coeffs,
unsigned  maxdeg,
double *  variances 
) const

Estimate variances of the coefficients obtained previously with the "sampleCoeffs" method.

◆ series()

double npstat::OSDE1D::series ( const double *  coeffs,
unsigned  maxdeg,
double  x 
) const

Polynomial series representing the fitted density

◆ supportRegion()

static std::pair<double,double> npstat::OSDE1D::supportRegion ( const AbsClassicalOrthoPoly1D poly1d,
double  leftLimit,
bool  leftIsFromSample,
double  rightLimit,
bool  rightIsFromSample,
unsigned long  nPoints 
)
static

A helper function for choosing the density support interval if it is unknown and has to be estimated from the sample. If the limit is not from the sample (as indicated by the corresponding boolean argument), it is left unchanged. The first element of the returned pair will correspond to the lower limit of the interval and the second element to the upper.

◆ tMatrix()

Matrix<double> npstat::OSDE1D::tMatrix ( unsigned  maxdeg) const

Matrix of unweighted pairwise integrals. The returned matrix will be dimensioned (maxdeg + 1) x (maxdeg + 1).

◆ twopoly()

std::pair<double,double> npstat::OSDE1D::twopoly ( unsigned  deg1,
unsigned  deg2,
double  x 
) const

A pair of polynomial values. Faster than calling "poly" two times.

◆ weightedSampleCoeffs()

template<class Pair >
void npstat::OSDE1D::weightedSampleCoeffs ( const Pair *  points,
unsigned long  numPoints,
double *  coeffs,
unsigned  maxdeg 
) const

Build the coefficients of the orthogonal polynomial series for the given sample of weighted points (empirical density function). The first element of the pair will be treated as the coordinate and the second element will be treated as weight. Weights must not be negative. The length of the array "coeffs" should be at least maxdeg + 1. Note that the coefficients are returned in the order of increasing degree (same order is used by the "series" function).

◆ weightedSampleCoeffsCovariance()

template<class Pair >
npstat::Matrix<double> npstat::OSDE1D::weightedSampleCoeffsCovariance ( const Pair *  points,
unsigned long  nPoints,
const double *  coeffs,
unsigned  maxdeg 
) const

Estimate covariances of the coefficients obtained previously with the "weightedSampleCoeffs" method. This code assumes that weights and coordinates are statistically independent from each other.

◆ weightedSampleCoeffsVariance()

template<class Pair >
void npstat::OSDE1D::weightedSampleCoeffsVariance ( const Pair *  points,
unsigned long  nPoints,
const double *  coeffs,
unsigned  maxdeg,
double *  variances 
) const

Estimate variances of the coefficients obtained previously with the "weightedSampleCoeffs" method. This code assumes that weights and coordinates are statistically independent from each other.


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