npstat is hosted by Hepforge, IPPP Durham
NPStat  5.10.0
StorablePolySeries1D.hh
Go to the documentation of this file.
1 #ifndef NPSTAT_STORABLEPOLYSERIES1D_HH_
2 #define NPSTAT_STORABLEPOLYSERIES1D_HH_
3 
4 /*!
5 // \file StorablePolySeries1D.hh
6 //
7 // \brief Storable functor for orthogonal polynomial series
8 //
9 // Author: I. Volobouev
10 //
11 // November 2018
12 */
13 
14 #include <vector>
15 #include <typeinfo>
16 
17 #include "geners/ClassId.hh"
18 
21 
22 namespace npstat {
23  class StorablePolySeries1D : public Functor1<double,double>
24 #ifndef SWIGBUG
25  , private AbsClassicalOrthoPoly1D
26 #endif
27  {
28  public:
29  /**
30  // The best way to create an object of this class is to call
31  // "makeStorablePolySeries" method of either ContOrthoPoly1D
32  // or AbsClassicalOrthoPoly1D class
33  */
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);
39 
40  inline virtual ~StorablePolySeries1D() {}
41 
42  inline virtual double operator()(const double& x) const
43  {
45  &coeffs_[0], coeffs_.size()-1U, x-shift_);
46  }
47 
48  //@{
49  /** A simple inspector of object properties */
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;}
57  //@}
58 
59  inline double poly(const unsigned deg, const double x) const
60  {return AbsClassicalOrthoPoly1D::poly(deg, x-shift_);}
61 
62  inline double series(const double* coeffs, const unsigned maxdeg,
63  const double x) const
64  {return AbsClassicalOrthoPoly1D::series(coeffs, maxdeg, x-shift_);}
65 
66  inline double integratePoly(const unsigned degree,
67  const unsigned power) const
68  {return integratePoly(degree, power, xmin_, xmax_);}
69 
70  double integratePoly(unsigned degree, unsigned power,
71  double xmin, double xmax) const;
72 
73  /** Set the series coefficients */
74  void setCoeffs(const double *coeffs, unsigned maxdeg);
75 
76  inline bool operator==(const StorablePolySeries1D& r) const
77  {return (typeid(*this) == typeid(r)) && this->isEqual(r);}
78  inline bool operator!=(const StorablePolySeries1D& r) const
79  {return !(*this == r);}
80 
81  // Methods needed for I/O
82  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
83  virtual bool write(std::ostream& os) const;
84 
85  static inline const char* classname()
86  {return "npstat::StorablePolySeries1D";}
87  static inline unsigned version() {return 1;}
88  static StorablePolySeries1D* read(const gs::ClassId& id, std::istream& in);
89 
90  inline virtual std::pair<long double,long double>
91  recurrenceCoeffs(const unsigned deg) const {return recur_.at(deg);}
92 
93  inline virtual std::pair<long double,long double>
94  monicRecurrenceCoeffs(const unsigned deg) const
95  {
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);
98  }
99 
100  protected:
101  virtual bool isEqual(const StorablePolySeries1D& other) const;
102 
103  private:
104  inline StorablePolySeries1D() {}
105 
106  inline virtual StorablePolySeries1D* clone() const
107  {return new StorablePolySeries1D(*this);}
108 
109  virtual long double weight(long double x) const;
110 
111  std::vector<double> coeffs_;
112  std::vector<std::pair<long double,long double> > recur_;
113  double xmin_;
114  double xmax_;
115  double shift_;
116 
117 #ifdef SWIG
118  public:
119  inline std::pair<double,double> monicRecurrenceCoeffs2(const unsigned deg) const
120  {
121  const std::pair<long double,long double>& p = monicRecurrenceCoeffs(deg);
122  return std::pair<double,double>(p.first, p.second);
123  }
124 
125  inline std::pair<double,double> recurrenceCoeffs2(const unsigned deg) const
126  {
127  const std::pair<long double,long double>& p = recurrenceCoeffs(deg);
128  return std::pair<double,double>(p.first, p.second);
129  }
130 #endif // SWIG
131  };
132 }
133 
134 #include "npstat/nm/StorablePolySeries1D.icc"
135 
136 #endif // NPSTAT_STORABLEPOLYSERIES1D_HH_
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