Hex  1.0
Hydrogen-electron collision solver
 All Data Structures Namespaces Files Functions Variables Typedefs Macros Pages
nodeintegrate.h
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *\
2  * *
3  * / / / / __ \ \ / / *
4  * / /__ / / / _ \ \ \/ / *
5  * / ___ / | |/_/ / /\ \ *
6  * / / / / \_\ / / \ \ *
7  * *
8  * Jakub Benda (c) 2014 *
9  * Charles University in Prague *
10  * *
11 \* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12 
13 #ifndef HEX_DIOPHANTINE
14 #define HEX_DIOPHANTINE
15 
16 #include <iostream>
17 #include <string>
18 
19 #include <gsl/gsl_sf.h>
20 
21 #include "misc.h"
22 #include "special.h"
23 
24 template <class Functor, class Integrator>
26 {
27  private:
28 
29  Integrator Q_;
30  Functor f_;
31 
32  double result_;
33  bool ok_;
34  std::string status_;
35 
36  double epsabs_;
37  double epsrel_;
38  int limit_;
39 
40 
41  public:
42 
43  NodeIntegrator (Functor f)
44  : Q_(f), f_(f), result_(0), ok_(true), status_(),
45  epsabs_(1e-8), epsrel_(1e-6), limit_(1000) {}
46 
47  double result () const { return result_; }
48  bool ok () const { return ok_; }
49  std::string const & status () const { return status_; }
50 
51  double epsabs () const { return epsabs_; }
52  void setEpsAbs (double eps) { epsabs_ = eps; }
53 
54  double epsrel () const { return epsrel_; }
55  void setEpsRel (double eps) { epsrel_ = eps; }
56 
57  int limit () const { return limit_; }
58  void setLimit (double n) { limit_ = n; }
59 
60  virtual double nthNode (int n) const = 0;
61  virtual double turningPoint () const { return 0; }
62 
63  bool integrate (double a, double b)
64  {
65  // initialize
66  ok_ = true;
67  status_ = "";
68  result_ = 0;
69 
70  // end of previous integration parcel
71  double prevR = 0;
72 
73  // start by integrating in the classically forbidden region (if any)
74  double rt = this->turningPoint();
75  if (rt != 0 and a < rt)
76  {
77  double r0 = a;
78  double r1 = std::min(rt,b);
79 
80  double quotient = 0.75; // FIXME variable ?
81  for (int samples = 8; ; samples *= 2)
82  {
83  rArray grid = geomspace(r0,r1,samples,quotient);
84  rArray eval = grid.transform(f_);
85  double estimate = special::integral::trapz(grid, eval);
86 
87  // check iteration limit
88  if ((int)log2(samples) == limit_)
89  {
90  throw exception
91  (
92  "Cannot integrate classically forbidden region - reached maximal number of iterations."
93  );
94  }
95 
96  // check convergence
97  if (std::abs(estimate - result_) < epsrel_ * std::abs(estimate))
98  {
99  result_ = estimate;
100  break;
101  }
102 
103  // check hopelessness
104  if (estimate == 0 and result_ == 0)
105  {
106  // this is probably hopeless...
107  break;
108  }
109 
110  // update estimate and go to the next iteration
111  result_ = estimate;
112  }
113 
114  prevR = r1;
115  }
116 
117  // check if the classically forbidden region was the only one to integrate
118  if (prevR == b)
119  return ok_;
120 
121  // for all integration parcels (nodes of the Bessel function)
122  for (int inode = 1; inode < limit_; inode++)
123  {
124  // get integration bounds
125  double rmin = prevR;
126  double rmax = nthNode(inode);
127 
128  // skip intervals that are below the lower limit
129  if (rmax < a)
130  continue;
131 
132  // skip intervals that are above the upper limit
133  if (rmin > b)
134  break;
135 
136  // shrink integration interval, if necessary
137  rmin = std::max (a, rmin);
138  rmax = std::min (rmax, b);
139 
140  // integrate and check success
141  if (not Q_.integrate(rmin,rmax))
142  {
143  ok_ = false;
144  status_ = format("Node integrator failed on <%g,%g>: %s", rmin, rmax, Q_.status().c_str());
145  return ok_;
146  }
147 
148  // update result
149  result_ += Q_.result();
150  prevR = rmax;
151 
152  // check convergence
153  if (std::abs(Q_.result()) < epsabs_ * std::abs(rmax - rmin) or std::abs(Q_.result()) < epsrel_ * std::abs(result_))
154  break;
155  }
156 
157  return ok_;
158  }
159 };
160 
161 template <class Functor, class Integrator>
162 class BesselNodeIntegrator : public NodeIntegrator<Functor,Integrator>
163 {
164  public:
165 
166  BesselNodeIntegrator (Functor f, double k, int l)
167  : NodeIntegrator<Functor,Integrator>(f), k_(k), l_(l) {}
168 
169  virtual double nthNode (int n) const
170  {
171  gsl_sf_result res;
172  int err = gsl_sf_bessel_zero_Jnu_e(l_+0.5,n,&res);
173  if (err != GSL_SUCCESS)
174  {
175  throw exception
176  (
177  "Cannot find %d-th root of the Bessel function j[%d](%g r) -- %s.",
178  n, l_, k_, gsl_strerror(err)
179  );
180  }
181  return res.val / k_;
182  }
183 
184  virtual double turningPoint () const
185  {
186  return std::sqrt(l_*(l_+1)) / k_;
187  };
188 
189  private:
190 
191  double k_;
192  int l_;
193 };
194 
195 template <class Functor, class Integrator>
196 class CoulombNodeIntegrator : public NodeIntegrator<Functor,Integrator>
197 {
198  public:
199 
200  CoulombNodeIntegrator (Functor f, double k, int l)
201  : NodeIntegrator<Functor,Integrator>(f), k_(k), l_(l), nzeros_(0), zeros_(nullptr) {}
202 
204  {
205  if (zeros_ != nullptr)
206  delete [] zeros_;
207  }
208 
209  virtual double nthNode (int n) const
210  {
211  if (n <= 0)
212  return 0;
213 
214  // compute more zeros if necessary
215  if (n > nzeros_)
216  {
217  // release old data
218  if (zeros_ != nullptr)
219  delete [] zeros_;
220 
221  // compute new number of nodes
222  nzeros_ = 2 * std::max(n, nzeros_);
223 
224  // compute the zeros
225  zeros_ = new double [nzeros_];
226  special::coulomb_zeros(-1/k_, l_, nzeros_, zeros_);
227  }
228 
229  // return the requested zero
230  return zeros_[n-1] / k_;
231  }
232 
233  virtual double turningPoint () const
234  {
235  return (std::sqrt(1 + l_*(l_+1)*k_*k_) - 1) / (k_*k_);
236  };
237 
238  private:
239 
240  double k_;
241  int l_;
242 
243  mutable int nzeros_;
244  mutable double * zeros_;
245 };
246 
247 template <class Functor, class Integrator>
248 class FixedNodeIntegrator : public NodeIntegrator<Functor,Integrator>
249 {
250  public:
251 
252  FixedNodeIntegrator (Functor f, const rArrayView zeros, double rt = 0)
253  : NodeIntegrator<Functor,Integrator>(f), zeros_(zeros), rt_(rt) {}
254 
255  virtual double nthNode (int n) const
256  {
257  if (n <= 0)
258  return 0;
259 
260  if (n <= (int)zeros_.size())
261  return zeros_[n-1];
262 
263  return special::constant::Inf;
264  }
265 
266  virtual double turningPoint () const
267  {
268  return rt_;
269  };
270 
271  private:
272 
273  double k_;
274  int l_;
275 
276  const rArray zeros_;
277  double rt_;
278 };
279 
280 #endif
double epsrel() const
Definition: nodeintegrate.h:54
Array view.
Definition: arrays.h:186
Definition: nodeintegrate.h:248
std::string const & status() const
Definition: nodeintegrate.h:49
virtual double nthNode(int n) const
Definition: nodeintegrate.h:209
CoulombNodeIntegrator(Functor f, double k, int l)
Definition: nodeintegrate.h:200
T trapz(double h, NumberArray< T > const &y)
Uniform trapezoidal integration.
Definition: special.h:78
A comfortable number array class.
Definition: arrays.h:171
T min(const ArrayView< T > a)
Minimal element of array.
Definition: arrays.h:1663
virtual double turningPoint() const
Definition: nodeintegrate.h:233
FixedNodeIntegrator(Functor f, const rArrayView zeros, double rt=0)
Definition: nodeintegrate.h:252
int limit() const
Definition: nodeintegrate.h:57
void setEpsRel(double eps)
Definition: nodeintegrate.h:55
virtual double turningPoint() const
Definition: nodeintegrate.h:184
const double e
Definition: special.h:44
void eval(TFunctor f, TArray grid, TArray &vals)
Definition: arrays.h:1844
double epsabs() const
Definition: nodeintegrate.h:51
~CoulombNodeIntegrator()
Definition: nodeintegrate.h:203
auto transform(Functor f) const -> NumberArray< decltype(f(T(0)))>
Apply a user transformation.
Definition: arrays.h:990
class NodeIntegrator exception
const double Inf
Definition: special.h:67
Definition: nodeintegrate.h:25
char const * format(Params...p)
printf-like formatting.
Definition: misc.h:233
bool integrate(double a, double b)
Definition: nodeintegrate.h:63
double result() const
Definition: nodeintegrate.h:47
int coulomb_zeros(double eta, int L, int nzeros, double *zeros, double epsrel=1e-8)
Get zeros of the Coulomb wave function .
Definition: special.cpp:60
void setLimit(double n)
Definition: nodeintegrate.h:58
virtual double nthNode(int n) const
Definition: nodeintegrate.h:255
void setEpsAbs(double eps)
Definition: nodeintegrate.h:52
bool ok() const
Definition: nodeintegrate.h:48
virtual double nthNode(int n) const
Definition: nodeintegrate.h:169
NumberArray< T > geomspace(T x0, T x1, size_t samples, double q)
Generate geometric grid.
Definition: arrays.h:1505
BesselNodeIntegrator(Functor f, double k, int l)
Definition: nodeintegrate.h:166
Definition: nodeintegrate.h:162
NumberArray< T > sqrt(NumberArray< T > const &A)
Return per-element square root.
Definition: arrays.h:1711
rArray abs(const cArrayView u)
Definition: arrays.cpp:19
virtual double turningPoint() const
Definition: nodeintegrate.h:61
virtual double nthNode(int n) const =0
virtual double turningPoint() const
Definition: nodeintegrate.h:266
T max(const ArrayView< T > a)
Maximal element of array.
Definition: arrays.h:1673
NodeIntegrator(Functor f)
Definition: nodeintegrate.h:43
size_t size() const
Item count.
Definition: arrays.h:673
Definition: nodeintegrate.h:196