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

Public Types

typedef npstat::PreciseType< double >::type precise_type
 
typedef std::pair< double, double > MeasurePoint
 
typedef GaussLegendreQuadrature PreciseQuadrature
 
typedef ContOrthoPoly1DDeg degree_functor
 

Public Member Functions

template<typename Numeric1 , typename Numeric2 >
 ContOrthoPoly1D (unsigned maxDegree, const std::vector< std::pair< Numeric1, Numeric2 > > &measure, OrthoPolyMethod m=OPOLY_STIELTJES, bool shiftMeasureCoords=true)
 
template<typename Numeric >
 ContOrthoPoly1D (unsigned maxDegree, const std::vector< Numeric > &coords, OrthoPolyMethod m=OPOLY_STIELTJES, bool shiftMeasureCoords=true)
 
template<typename Numeric1 , typename Numeric2 , class OrthoPoly >
 ContOrthoPoly1D (unsigned maxDegree, const std::vector< std::pair< Numeric1, Numeric2 > > &measure, const OrthoPoly &ops, bool shiftMeasureCoords=false)
 
template<typename Numeric , class OrthoPoly >
 ContOrthoPoly1D (unsigned maxDegree, const std::vector< Numeric > &coords, const OrthoPoly &ops, bool shiftMeasureCoords=false)
 
ContOrthoPoly1Dclone () const
 
unsigned maxDegree () const
 
unsigned long measureLength () const
 
double weight (const unsigned long index) const
 
precise_type sumOfWeights () const
 
precise_type sumOfWeightSquares () const
 
bool areAllWeightsEqual () const
 
double minCoordinate () const
 
double maxCoordinate () const
 
precise_type measureShift () const
 
precise_type meanCoordinate () const
 
precise_type location () const
 
precise_type scale () const
 
MeasurePoint measure (const unsigned long index) const
 
double coordinate (const unsigned long index) const
 
double effectiveSampleSize () const
 
double poly (unsigned deg, double x) const
 
std::pair< double, double > polyPair (unsigned deg1, unsigned deg2, double x) const
 
template<typename Numeric >
void allPolys (const double x, const unsigned maxdeg, Numeric *polyValues) const
 
void allpoly (const long double x, long double *values, const unsigned maxdeg) const
 
template<typename Numeric >
void sampleAverages (const Numeric *coords, unsigned long lenCoords, double *averages, unsigned maxdeg) const
 
template<class Pair >
void weightedPointsAverages (const Pair *points, unsigned long numPoints, double *averages, unsigned maxdeg) const
 
template<typename Numeric >
Matrix< double > sampleProductAverages (const Numeric *coords, unsigned long lenCoords, unsigned maxdeg) const
 
template<typename Numeric >
double series (const Numeric *coeffs, unsigned deg, double x) const
 
template<class Functor , typename Numeric >
void calculateCoeffs (const Functor &fcn, Numeric *coeffs, unsigned maxdeg) const
 
template<class Functor >
double constrainedCoeffs (const Functor &fcn, double *coeffs, unsigned maxdeg, double xmin, double xmax, double integralValue=1.0) const
 
void weightCoeffs (double *coeffs, unsigned maxdeg) const
 
template<class Functor >
double integratedSquaredError (const Functor &fcn, const double *coeffs, unsigned maxdeg, double xmin, double xmax, unsigned integrationRulePoints=0U) const
 
template<class Functor >
double weightedSquaredError (const Functor &fcn, const double *coeffs, unsigned maxdeg) const
 
precise_type empiricalKroneckerDelta (unsigned deg1, unsigned deg2) const
 
Matrix< precise_type > empiricalKroneckerMatrix (unsigned maxdeg) const
 
double empiricalKroneckerCovariance (unsigned deg1, unsigned deg2, unsigned deg3, unsigned deg4) const
 
Matrix< precise_type > jacobiMatrix (unsigned n) const
 
precise_type integratePoly (unsigned degree, unsigned power, double xmin, double xmax) const
 
precise_type integrateTripleProduct (unsigned deg1, unsigned deg2, unsigned deg3, double xmin, double xmax) const
 
Matrix< precise_type > tripleProductMatrix (const std::vector< std::pair< unsigned, unsigned > > &pairs, unsigned maxdeg, double xmin, double xmax) const
 
template<class Functor >
double integrateEWPolyProduct (const Functor &weight, unsigned deg1, unsigned deg2, double xmin, double xmax, unsigned integrationRulePoints) const
 
template<class Functor >
void extWeightAverages (const Functor &weight, const AbsIntervalQuadrature1D &quad, unsigned maxdeg, double xmin, double xmax, double *averages) const
 
template<class Functor >
Matrix< double > extWeightProductAverages (const Functor &weight, const AbsIntervalQuadrature1D &quad, unsigned maxdeg, double xmin, double xmax) const
 
double powerAverage (unsigned deg, unsigned power) const
 
double jointPowerAverage (unsigned deg1, unsigned power1, unsigned deg2, unsigned power2) const
 
double jointAverage (const unsigned *degrees, unsigned nDegrees, bool degreesSortedIncreasingOrder=false) const
 
double cachedJointAverage (unsigned deg1, unsigned deg2, unsigned deg3, unsigned deg4) const
 
double cachedJointAverage (unsigned deg1, unsigned deg2, unsigned deg3, unsigned deg4, unsigned deg5, unsigned deg6) const
 
double cachedJointAverage (unsigned deg1, unsigned deg2, unsigned deg3, unsigned deg4, unsigned deg5, unsigned deg6, unsigned deg7, unsigned deg8) const
 
double cov4 (unsigned a, unsigned b, unsigned c, unsigned d) const
 
double cov6 (unsigned a, unsigned b, unsigned c, unsigned d, unsigned e, unsigned f) const
 
double cov8 (unsigned a, unsigned b, unsigned c, unsigned d, unsigned e, unsigned f, unsigned g, unsigned h) const
 
double covCov4 (unsigned a, unsigned b, unsigned c, unsigned d, unsigned e, unsigned f, unsigned g, unsigned h) const
 
double slowCov8 (unsigned a, unsigned b, unsigned c, unsigned d, unsigned e, unsigned f, unsigned g, unsigned h) const
 
double epsExpectation (unsigned m, unsigned n, bool highOrder) const
 
double epsCovariance (unsigned i, unsigned j, unsigned m, unsigned n, bool highOrder) const
 
Matrix< double > epsCovarianceMatrix (const std::vector< std::pair< unsigned, unsigned > > &pairs, bool highOrder) const
 
CPP11_auto_ptr< StorablePolySeries1DmakeStorablePolySeries (double xmin, double xmax, const double *coeffs=0, unsigned maxdeg=0) const
 
std::pair< precise_type, precise_type > recurrenceCoeffs (const unsigned k) const
 
std::pair< precise_type, precise_type > monicRecurrenceCoeffs (const unsigned k) const
 

Friends

class PolyFcn
 
class TripleProdFcn
 
template<class Funct >
class EWPolyProductFcn
 

Constructor & Destructor Documentation

◆ ContOrthoPoly1D() [1/4]

template<typename Numeric1 , typename Numeric2 >
npstat::ContOrthoPoly1D::ContOrthoPoly1D ( unsigned  maxDegree,
const std::vector< std::pair< Numeric1, Numeric2 > > &  measure,
OrthoPolyMethod  m = OPOLY_STIELTJES,
bool  shiftMeasureCoords = true 
)

Main constructor, with obvious arguments. The first element of the measure pair is the coordinate and the second element is the weight (all weights must be non-negative). Internally, the measure will be sorted in the order of increasing weight and the measure coordinates will normally be shifted so that their weighted mean is at 0.

◆ ContOrthoPoly1D() [2/4]

template<typename Numeric >
npstat::ContOrthoPoly1D::ContOrthoPoly1D ( unsigned  maxDegree,
const std::vector< Numeric > &  coords,
OrthoPolyMethod  m = OPOLY_STIELTJES,
bool  shiftMeasureCoords = true 
)

Constructor which assumes that all weights are 1.0

◆ ContOrthoPoly1D() [3/4]

template<typename Numeric1 , typename Numeric2 , class OrthoPoly >
npstat::ContOrthoPoly1D::ContOrthoPoly1D ( unsigned  maxDegree,
const std::vector< std::pair< Numeric1, Numeric2 > > &  measure,
const OrthoPoly &  ops,
bool  shiftMeasureCoords = false 
)

Constructor that uses the modified Chebyshev algorithm

◆ ContOrthoPoly1D() [4/4]

template<typename Numeric , class OrthoPoly >
npstat::ContOrthoPoly1D::ContOrthoPoly1D ( unsigned  maxDegree,
const std::vector< Numeric > &  coords,
const OrthoPoly &  ops,
bool  shiftMeasureCoords = false 
)

Constructor that uses the modified Chebyshev algorithm and assumes that all weights are 1.0

Member Function Documentation

◆ allpoly()

void npstat::ContOrthoPoly1D::allpoly ( const long double  x,
long double *  values,
const unsigned  maxdeg 
) const
inline

Function added for interface compatibility with AbsClassicalOrthoPoly1D

◆ allPolys()

template<typename Numeric >
void npstat::ContOrthoPoly1D::allPolys ( const double  x,
const unsigned  maxdeg,
Numeric *  polyValues 
) const
inline

Return the values of all orthonormal polynomials up to some maximum degree. Size of the "polyValues" array should be at least maxdeg + 1.

◆ cachedJointAverage()

double npstat::ContOrthoPoly1D::cachedJointAverage ( unsigned  deg1,
unsigned  deg2,
unsigned  deg3,
unsigned  deg4 
) const

A measure-weighted average that will be remembered internally so that another call to this function with the same arguments will return the answer quickly

◆ calculateCoeffs()

template<class Functor , typename Numeric >
void npstat::ContOrthoPoly1D::calculateCoeffs ( const Functor &  fcn,
Numeric *  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).

◆ clone()

ContOrthoPoly1D* npstat::ContOrthoPoly1D::clone ( ) const
inline

Function added for interface compatibility with AbsClassicalOrthoPoly1D

◆ constrainedCoeffs()

template<class Functor >
double npstat::ContOrthoPoly1D::constrainedCoeffs ( const Functor &  fcn,
double *  coeffs,
unsigned  maxdeg,
double  xmin,
double  xmax,
double  integralValue = 1.0 
) const

Build the coefficients of the orthogonal polynomial series for the given function in such a way that the integral of the truncated series on the [xmin, xmax] interval is constrained to "integralValue". 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). The construction minimizes the ISE weighted by the empirical density.

The function returns the chi-square from the corresponding constrained least squares problem. This is the sum of (coeffs[i] - c[i])^2, 0 <= i <= maxdeg, where c[i] are determined by the "calculateCoeffs" method.

◆ cov4()

double npstat::ContOrthoPoly1D::cov4 ( unsigned  a,
unsigned  b,
unsigned  c,
unsigned  d 
) const

Covariance between v_{ab} and v_{cd}

◆ cov6()

double npstat::ContOrthoPoly1D::cov6 ( unsigned  a,
unsigned  b,
unsigned  c,
unsigned  d,
unsigned  e,
unsigned  f 
) const

Covariance between v_{ab}*v_{cd} and v_{ef}

◆ cov8()

double npstat::ContOrthoPoly1D::cov8 ( unsigned  a,
unsigned  b,
unsigned  c,
unsigned  d,
unsigned  e,
unsigned  f,
unsigned  g,
unsigned  h 
) const

Covariance between v_{ab}*v_{cd} and v_{ef}*v_{gh}

◆ covCov4()

double npstat::ContOrthoPoly1D::covCov4 ( unsigned  a,
unsigned  b,
unsigned  c,
unsigned  d,
unsigned  e,
unsigned  f,
unsigned  g,
unsigned  h 
) const

Covariance between two cov4 estimates (i.e., error on the error)

◆ effectiveSampleSize()

double npstat::ContOrthoPoly1D::effectiveSampleSize ( ) const

Kish's effective sample size for the measure

◆ empiricalKroneckerCovariance()

double npstat::ContOrthoPoly1D::empiricalKroneckerCovariance ( unsigned  deg1,
unsigned  deg2,
unsigned  deg3,
unsigned  deg4 
) const

If the Kronecker deltas were calculated from a sample, they would be random and would have a covariance matrix. This is an estimate of the terms in that covariance matrix. The covariance is between delta_{deg1,deg2} and delta_{deg3,deg4}.

◆ empiricalKroneckerDelta()

precise_type npstat::ContOrthoPoly1D::empiricalKroneckerDelta ( unsigned  deg1,
unsigned  deg2 
) const

This method is useful for testing the numerical precision of the orthonormalization procedure. It returns the scalar products between various polynomials.

◆ empiricalKroneckerMatrix()

Matrix<precise_type> npstat::ContOrthoPoly1D::empiricalKroneckerMatrix ( unsigned  maxdeg) const

A faster way to generate Kronecker deltas if the whole matrix of them is needed. The argument must not exceed the return value of the "maxDegree" function. The resulting matrix will be dimensioned (maxdeg + 1) x (maxdeg + 1).

◆ epsCovariance()

double npstat::ContOrthoPoly1D::epsCovariance ( unsigned  i,
unsigned  j,
unsigned  m,
unsigned  n,
bool  highOrder 
) const

Estimate of covariance between eps_{ij} and eps_{mn}.

The result should be identical with "empiricalKroneckerCovariance" in case "highOrder" argument is "false". However, this method caches results of various calculations of averages and should perform better inside cycles over indices.

◆ epsCovarianceMatrix()

Matrix<double> npstat::ContOrthoPoly1D::epsCovarianceMatrix ( const std::vector< std::pair< unsigned, unsigned > > &  pairs,
bool  highOrder 
) const

Covariance matrix created by multiple calls to "epsCovariance"

◆ epsExpectation()

double npstat::ContOrthoPoly1D::epsExpectation ( unsigned  m,
unsigned  n,
bool  highOrder 
) const

Estimate of eps_{mn} expectation

◆ extWeightAverages()

template<class Functor >
void npstat::ContOrthoPoly1D::extWeightAverages ( const Functor &  weight,
const AbsIntervalQuadrature1D quad,
unsigned  maxdeg,
double  xmin,
double  xmax,
double *  averages 
) const

Integral of the polynomials with some external weight function (e.g., reconstructed or oracle density). The length of array "averages" should be at least maxdeg + 1.

◆ extWeightProductAverages()

template<class Functor >
Matrix<double> npstat::ContOrthoPoly1D::extWeightProductAverages ( const Functor &  weight,
const AbsIntervalQuadrature1D quad,
unsigned  maxdeg,
double  xmin,
double  xmax 
) const

Integral of all products of two orthonormal polynomials with some external weight function (e.g., oracle density).

◆ integratedSquaredError()

template<class Functor >
double npstat::ContOrthoPoly1D::integratedSquaredError ( const Functor &  fcn,
const double *  coeffs,
unsigned  maxdeg,
double  xmin,
double  xmax,
unsigned  integrationRulePoints = 0U 
) const

Integrated squared error between the given function and the polynomial series. Parameter "integrationRulePoints" specifies the parameter "npoints" for the GaussLegendreQuadrature class. If "integrationRulePoints" is 0, the rule will be picked automatically in such a way that it integrates a polynomial of degree maxdeg*2 exactly.

◆ integrateEWPolyProduct()

template<class Functor >
double npstat::ContOrthoPoly1D::integrateEWPolyProduct ( const Functor &  weight,
unsigned  deg1,
unsigned  deg2,
double  xmin,
double  xmax,
unsigned  integrationRulePoints 
) const

Integral of the product of two orthonormal polynomials with some external weight function (e.g., oracle density). "EW" in the method name stands for "externally weighted".

"integrationRulePoints" argument will be passed to the GaussLegendreQuadrature constructor and must be supported by that class.

◆ integratePoly()

precise_type npstat::ContOrthoPoly1D::integratePoly ( unsigned  degree,
unsigned  power,
double  xmin,
double  xmax 
) const

Integral of an orthonormal polynomials to some power without the weight (so this is not an inner product with the poly of degree 0).

◆ integrateTripleProduct()

precise_type npstat::ContOrthoPoly1D::integrateTripleProduct ( unsigned  deg1,
unsigned  deg2,
unsigned  deg3,
double  xmin,
double  xmax 
) const

Integral of the product of three orthonormal polynomials without the weight (so this is not an inner product)

◆ jacobiMatrix()

Matrix<precise_type> npstat::ContOrthoPoly1D::jacobiMatrix ( unsigned  n) const

Generate principal minor of order n of the Jacobi matrix. n must not exceed "maxDegree". Note that, if you calculate roots using this Jacobi matrix, they will be shifted by meanCoordinate() (that is, add meanCoordinate() to all such roots).

◆ jointAverage()

double npstat::ContOrthoPoly1D::jointAverage ( const unsigned *  degrees,
unsigned  nDegrees,
bool  degreesSortedIncreasingOrder = false 
) const

A measure-weighted average of the product of several orthonormal poly values

◆ jointPowerAverage()

double npstat::ContOrthoPoly1D::jointPowerAverage ( unsigned  deg1,
unsigned  power1,
unsigned  deg2,
unsigned  power2 
) const

A measure-weighted average of the product of two orthonormal poly values to some powers

◆ makeStorablePolySeries()

CPP11_auto_ptr<StorablePolySeries1D> npstat::ContOrthoPoly1D::makeStorablePolySeries ( double  xmin,
double  xmax,
const double *  coeffs = 0,
unsigned  maxdeg = 0 
) const

Construct a StorablePolySeries1D object out of this one

◆ maxDegree()

unsigned npstat::ContOrthoPoly1D::maxDegree ( ) const
inline

A simple inspector of object properties

◆ measure()

MeasurePoint npstat::ContOrthoPoly1D::measure ( const unsigned long  index) const
inline

Note that this returns the internal, shifted coordinate. Add the mean coordinate to restore the original coordinate.

◆ monicRecurrenceCoeffs()

std::pair<precise_type,precise_type> npstat::ContOrthoPoly1D::monicRecurrenceCoeffs ( const unsigned  k) const
inline

Examine recurrence coefficients for the corresponding monic polynomials

◆ poly()

double npstat::ContOrthoPoly1D::poly ( unsigned  deg,
double  x 
) const

Return the value of one of the normalized polynomials

◆ polyPair()

std::pair<double,double> npstat::ContOrthoPoly1D::polyPair ( unsigned  deg1,
unsigned  deg2,
double  x 
) const

Return the values of two of the normalized polynomials

◆ powerAverage()

double npstat::ContOrthoPoly1D::powerAverage ( unsigned  deg,
unsigned  power 
) const

A measure-weighted average of orthonormal poly values to some power

◆ recurrenceCoeffs()

std::pair<precise_type,precise_type> npstat::ContOrthoPoly1D::recurrenceCoeffs ( const unsigned  k) const
inline

Examine recurrence coefficients

◆ sampleAverages()

template<typename Numeric >
void npstat::ContOrthoPoly1D::sampleAverages ( const Numeric *  coords,
unsigned long  lenCoords,
double *  averages,
unsigned  maxdeg 
) const

Average values of the polynomials for a sample of points. Size of the "averages" array should be at least maxdeg + 1. Note that the resulting averages are not weighted by the measure.

◆ sampleProductAverages()

template<typename Numeric >
Matrix<double> npstat::ContOrthoPoly1D::sampleProductAverages ( const Numeric *  coords,
unsigned long  lenCoords,
unsigned  maxdeg 
) const

Average values of the pairwise polynomial products for a sample of points. The returned matrix will be symmetric and will have the dimensions (maxdeg + 1) x (maxdeg + 1). Note that the resulting averages are not weighted by the measure.

◆ series()

template<typename Numeric >
double npstat::ContOrthoPoly1D::series ( const Numeric *  coeffs,
unsigned  deg,
double  x 
) const

Return the value of the orthonormal polynomial series

◆ slowCov8()

double npstat::ContOrthoPoly1D::slowCov8 ( unsigned  a,
unsigned  b,
unsigned  c,
unsigned  d,
unsigned  e,
unsigned  f,
unsigned  g,
unsigned  h 
) const

Alternative implementation of "cov8". It is slower than "cov8" but easier to program and verify. Useful mainly for testing "cov8".

◆ tripleProductMatrix()

Matrix<precise_type> npstat::ContOrthoPoly1D::tripleProductMatrix ( const std::vector< std::pair< unsigned, unsigned > > &  pairs,
unsigned  maxdeg,
double  xmin,
double  xmax 
) const

Unweighted triple product integrals arranged in a matrix. The returned matrix will be dimensioned pairs.size() x (maxdeg+1).

◆ weightCoeffs()

void npstat::ContOrthoPoly1D::weightCoeffs ( double *  coeffs,
unsigned  maxdeg 
) const

Build the coefficients of the orthogonal polynomial series for the discrete weight values (that is, fcn(x_i) = w_i, using the terminology of the "calculateCoeffs" method). This can sometimes be useful for weighted measures. Of course, for unweighted measures this results in just delta_{deg,0}.

◆ weightedPointsAverages()

template<class Pair >
void npstat::ContOrthoPoly1D::weightedPointsAverages ( const Pair *  points,
unsigned long  numPoints,
double *  averages,
unsigned  maxdeg 
) const

Averages of the polynomial values for the given sample of weighted points (not weighting by the measure). 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 array "averages" should be at least maxdeg + 1.

◆ weightedSquaredError()

template<class Functor >
double npstat::ContOrthoPoly1D::weightedSquaredError ( const Functor &  fcn,
const double *  coeffs,
unsigned  maxdeg 
) const

Squared error between the given function and the polynomial series, weighted according to the measure


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