1 #ifndef NPSTAT_STORABLEPOLYSERIES1D_HH_
2 #define NPSTAT_STORABLEPOLYSERIES1D_HH_
17 #include "geners/ClassId.hh"
34 template <
typename Real>
36 const std::vector<std::pair<Real,Real> >& rCoeffs,
37 double xmin,
double xmax,
double shift = 0.0,
38 const double *coeffs = 0,
unsigned maxdeg = 0);
42 inline virtual double operator()(
const double& x)
const
45 &coeffs_[0], coeffs_.size()-1U, x-shift_);
50 inline double xmin()
const {
return xmin_;}
51 inline double xmax()
const {
return xmax_;}
52 inline unsigned maxDegCoeffs()
const {
return coeffs_.size()-1U;}
53 inline double getCoefficient(
const unsigned i) {
return coeffs_.at(i);}
54 inline unsigned maxDegPoly()
const {
return this->
maxDegree();}
55 inline double getShift()
const {
return shift_;}
56 inline virtual unsigned maxDegree()
const {
return recur_.size() - 1U;}
59 inline double poly(
const unsigned deg,
const double x)
const
62 inline double series(
const double* coeffs,
const unsigned maxdeg,
66 inline double integratePoly(
const unsigned degree,
67 const unsigned power)
const
68 {
return integratePoly(degree, power, xmin_, xmax_);}
70 double integratePoly(
unsigned degree,
unsigned power,
71 double xmin,
double xmax)
const;
74 void setCoeffs(
const double *coeffs,
unsigned maxdeg);
77 {
return (
typeid(*
this) ==
typeid(r)) && this->isEqual(r);}
79 {
return !(*
this == r);}
82 virtual gs::ClassId classId()
const {
return gs::ClassId(*
this);}
83 virtual bool write(std::ostream& os)
const;
85 static inline const char* classname()
86 {
return "npstat::StorablePolySeries1D";}
87 static inline unsigned version() {
return 1;}
90 inline virtual std::pair<long double,long double>
91 recurrenceCoeffs(
const unsigned deg)
const {
return recur_.at(deg);}
93 inline virtual std::pair<long double,long double>
94 monicRecurrenceCoeffs(
const unsigned deg)
const
96 const std::pair<long double,long double>& c = recur_.at(deg);
97 return std::pair<long double,long double>(c.first, c.second*c.second);
101 virtual bool isEqual(
const StorablePolySeries1D& other)
const;
109 virtual long double weight(
long double x)
const;
111 std::vector<double> coeffs_;
112 std::vector<std::pair<long double,long double> > recur_;
119 inline std::pair<double,double> monicRecurrenceCoeffs2(
const unsigned deg)
const
121 const std::pair<long double,long double>& p = monicRecurrenceCoeffs(deg);
122 return std::pair<double,double>(p.first, p.second);
125 inline std::pair<double,double> recurrenceCoeffs2(
const unsigned deg)
const
127 const std::pair<long double,long double>& p = recurrenceCoeffs(deg);
128 return std::pair<double,double>(p.first, p.second);
134 #include "npstat/nm/StorablePolySeries1D.icc"
Base class for classical continuous orthonormal polynomials.
Interface definitions and concrete simple functors for a variety of functor-based calculations.
Definition: AbsClassicalOrthoPoly1D.hh:32
double series(const double *coeffs, unsigned maxdeg, double x) const
long double poly(unsigned deg, long double x) const
Definition: StorablePolySeries1D.hh:27
double xmin() const
Definition: StorablePolySeries1D.hh:50
virtual unsigned maxDegree() const
Definition: StorablePolySeries1D.hh:56
void setCoeffs(const double *coeffs, unsigned maxdeg)
StorablePolySeries1D(const std::vector< std::pair< Real, Real > > &rCoeffs, double xmin, double xmax, double shift=0.0, const double *coeffs=0, unsigned maxdeg=0)
Definition: AbsArrayProjector.hh:14
Definition: SimpleFunctors.hh:58