CMS 3D CMS Logo

ConstantStepOdeSolver.h
Go to the documentation of this file.
1 #ifndef CalibCalorimetry_HcalAlgos_ConstantStepOdeSolver_h_
2 #define CalibCalorimetry_HcalAlgos_ConstantStepOdeSolver_h_
3 
4 #include <vector>
5 #include <iostream>
7 
9 
10 //
11 // ODE solver with a constant time step. The derived classes are supposed
12 // to implement concrete ODE solving algorithms (Runge-Kutta, etc).
13 //
15 {
16 public:
18  : rhs_(nullptr), dt_(0.0), dim_(0), runLen_(0), lastIntegrated_(0) {}
19 
20  inline ConstantStepOdeSolver(const AbsODERHS& rhs) :
21  rhs_(nullptr), dt_(0.0), dim_(0), runLen_(0), lastIntegrated_(0)
22  {
23  rhs_ = rhs.clone();
24  }
25 
26  // The copy constructor and the assignment operator are explicitly provided
29 
30  inline virtual ~ConstantStepOdeSolver() {delete rhs_;}
31 
32  // Access the equation right hand side
33  inline void setRHS(const AbsODERHS& rhs)
34  {
35  delete rhs_;
36  rhs_ = rhs.clone();
37  }
38  inline const AbsODERHS* getRHS() const {return rhs_;}
39  inline AbsODERHS* getRHS() {return rhs_;}
40 
41  // Inspectors (will be valid after at least one "run" call)
42  inline unsigned lastDim() const {return dim_;}
43  inline unsigned lastRunLength() const {return runLen_;}
44  inline double lastDeltaT() const {return dt_;}
45  inline double lastMaxT() const {return runLen_ ? dt_*(runLen_-1U) : 0.0;}
46 
47  inline double getTime(const unsigned idx) const
48  {
49  if (idx >= runLen_) throw cms::Exception(
50  "In ConstantStepOdeSolver::getTime: index out of range");
51  return idx*dt_;
52  }
53 
54  inline double getCoordinate(const unsigned which, const unsigned idx) const
55  {
56  if (which >= dim_ || idx >= runLen_) throw cms::Exception(
57  "In ConstantStepOdeSolver::getCoordinate: index out of range");
58  return historyBuffer_[dim_*idx + which];
59  }
60 
61  // Integrate the node with the given number and get
62  // the value of the integral at the given history point
63  double getIntegrated(unsigned which, unsigned idx) const;
64 
65  // Truncate some coordinate
66  void truncateCoordinate(unsigned which, double minValue, double maxValue);
67 
68  // Linear interpolation methods will be used in case the "cubic"
69  // argument is false, and cubic in case it is true
70  double interpolateCoordinate(unsigned which, double t,
71  bool cubic = false) const;
72 
73  // Interpolate the integrated node
74  double interpolateIntegrated(unsigned which, double t,
75  bool cubic = false) const;
76 
77  // Get the time of the peak position
78  double getPeakTime(unsigned which) const;
79 
80  // Solve the ODE and remember the history
81  void run(const double* initialConditions, unsigned lenConditions,
82  double dt, unsigned nSteps);
83 
84  // Set the history from some external source. The size
85  // of the data array should be at least dim*runLen.
86  void setHistory(double dt, const double* data,
87  unsigned dim, unsigned runLen);
88 
89  // Write out the history
90  void writeHistory(std::ostream& os, double dt, bool cubic = false) const;
91 
92  // Write out the integrated node
93  void writeIntegrated(std::ostream& os, unsigned which,
94  double dt, bool cubic = false) const;
95 
96  // The following method must be overriden by derived classes
97  virtual const char* methodName() const = 0;
98 
99 protected:
101 
102 private:
103  // The following method must be overriden by derived classes
104  virtual void step(double t, double dt,
105  const double* x, unsigned lenX,
106  double* coordIncrement) const = 0;
107 
108  // The following integration corresponds to the cubic
109  // interpolation of the coordinate
110  void integrateCoordinate(const unsigned which);
111 
112  double dt_;
113  unsigned dim_;
114  unsigned runLen_;
115  unsigned lastIntegrated_;
116 
117  std::vector<double> historyBuffer_;
118  std::vector<double> chargeBuffer_;
119 };
120 
121 
122 // Dump the coordinate history as it was collected
123 inline std::ostream& operator<<(std::ostream& os,
124  const ConstantStepOdeSolver& s)
125 {
126  s.writeHistory(os, s.lastDeltaT());
127  return os;
128 }
129 
130 
131 // A few concrete ODE solvers
133 {
134 public:
136 
137  inline explicit EulerOdeSolver(const AbsODERHS& rhs)
138  : ConstantStepOdeSolver(rhs) {}
139 
140  inline const char* methodName() const override {return "Euler";}
141 
142 private:
143  void step(double t, double dt,
144  const double* x, unsigned lenX,
145  double* coordIncrement) const override;
146 };
147 
148 
150 {
151 public:
152  inline RK2() : ConstantStepOdeSolver() {}
153 
154  inline explicit RK2(const AbsODERHS& rhs) : ConstantStepOdeSolver(rhs) {}
155 
156  inline const char* methodName() const override {return "2nd order Runge-Kutta";}
157 
158 private:
159  void step(double t, double dt,
160  const double* x, unsigned lenX,
161  double* coordIncrement) const override;
162 
163  mutable std::vector<double> buf_;
164 };
165 
166 
168 {
169 public:
170  inline RK4() : ConstantStepOdeSolver() {}
171 
172  inline explicit RK4(const AbsODERHS& rhs) : ConstantStepOdeSolver(rhs) {}
173 
174  inline const char* methodName() const override {return "4th order Runge-Kutta";}
175 
176 private:
177  void step(double t, double dt,
178  const double* x, unsigned lenX,
179  double* coordIncrement) const override;
180 
181  mutable std::vector<double> buf_;
182 };
183 
184 #endif // CalibCalorimetry_HcalAlgos_ConstantStepOdeSolver_h_
double getCoordinate(const unsigned which, const unsigned idx) const
float dt
Definition: AMPTWrapper.h:126
RK4(const AbsODERHS &rhs)
double interpolateIntegrated(unsigned which, double t, bool cubic=false) const
std::ostream & operator<<(std::ostream &os, const ConstantStepOdeSolver &s)
double getTime(const unsigned idx) const
#define nullptr
ConstantStepOdeSolver(const AbsODERHS &rhs)
std::vector< double > chargeBuffer_
virtual AbsODERHS * clone() const =0
double getPeakTime(unsigned which) const
double getIntegrated(unsigned which, unsigned idx) const
void writeIntegrated(std::ostream &os, unsigned which, double dt, bool cubic=false) const
const char * methodName() const override
std::vector< double > buf_
EulerOdeSolver(const AbsODERHS &rhs)
void setRHS(const AbsODERHS &rhs)
void integrateCoordinate(const unsigned which)
const AbsODERHS * getRHS() const
ConstantStepOdeSolver & operator=(const ConstantStepOdeSolver &r)
const char * methodName() const override
double interpolateCoordinate(unsigned which, double t, bool cubic=false) const
std::vector< double > historyBuffer_
void run(const double *initialConditions, unsigned lenConditions, double dt, unsigned nSteps)
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
std::vector< double > buf_
void truncateCoordinate(unsigned which, double minValue, double maxValue)
unsigned lastRunLength() const
const char * methodName() const override
virtual const char * methodName() const =0
def which(cmd)
Definition: eostools.py:335
void writeHistory(std::ostream &os, double dt, bool cubic=false) const
void setHistory(double dt, const double *data, unsigned dim, unsigned runLen)
virtual void step(double t, double dt, const double *x, unsigned lenX, double *coordIncrement) const =0
RK2(const AbsODERHS &rhs)