ampsci
High-precision calculations for one- and two-valence atomic systems
AdamsMoulton.hpp
1#pragma once
2#include <algorithm>
3#include <array>
4#include <cassert>
5#include <complex>
6#include <tuple>
7#include <type_traits>
8
9//! Contains classes and functions which use general N-step Adams Moulton method
10//! to solve systems of 2x2 ODEs, up to N=12.
11namespace AdamsMoulton {
12
13//==============================================================================
14
15/*!
16 @brief Pure-virtual struct defining the derivative matrix for a 2x2 ODE system.
17 @details
18 Used by ODESolver2D to define the derivative matrix D and optional inhomogeneous
19 term S. Derive from this and implement a(t), b(t), c(t), d(t) to define the ODE.
20
21 The system of ODEs is:
22
23 \f[ \frac{dF(t)}{dt} = D(t) F(t) + S(t) \f]
24
25 where:
26
27 \f[
28 F(t) = \begin{pmatrix} f(t)\\ g(t) \end{pmatrix}, \quad
29 D(t) = \begin{pmatrix} a(t) & b(t)\\ c(t) & d(t) \end{pmatrix}, \quad
30 S(t) = \begin{pmatrix} s_f(t)\\ s_g(t) \end{pmatrix}.
31 \f]
32
33 D (and optionally S) must be provided by implementing this struct.
34 The four functions a, b, c, d must be overridden to define the ODE system.
35 Sf and Sg default to zero if not overridden.
36
37 Template parameter T is the type of the argument t (usually double or
38 complex<double>, but may be an index type such as std::size_t if the matrix
39 is known only at discrete grid points). Template parameter Y is the return
40 type (usually double, but may be float or complex<double>).
41
42 @note In benchmarks, deriving from a struct was significantly faster than
43 using std::function, slightly faster than function pointers, and
44 comparable to a direct implementation.
45*/
46template <typename T = double, typename Y = double>
48 //! a,b,c,d are derivative matrix functions; all must be user implemented
49 virtual Y a(T t) const = 0;
50 virtual Y b(T t) const = 0;
51 virtual Y c(T t) const = 0;
52 virtual Y d(T t) const = 0;
53 //! Sf and Sg are optional inhomogenous terms
54 virtual Y Sf(T) const { return Y(0); };
55 virtual Y Sg(T) const { return Y(0); };
56 virtual ~DerivativeMatrix() = default;
57};
58
59//==============================================================================
60
61// User-defined type-trait: Checks whether T is a std::complex type
62template <typename T>
63struct is_complex : std::false_type {};
64// User-defined type-trait: Checks whether T is a std::complex type
65template <typename T>
66struct is_complex<std::complex<T>> : std::true_type {};
67/*!
68 @brief Type trait: true if T is std::complex<U> for some U, false otherwise.
69 @details
70 Examples:
71 ```cpp
72 static_assert(!is_complex_v<double>);
73 static_assert(!is_complex_v<float>);
74 static_assert(is_complex_v<std::complex<double>>);
75 static_assert(is_complex_v<std::complex<float>>);
76 static_assert(is_complex<std::complex<float>>::value);
77 ```
78*/
79template <typename T>
80constexpr bool is_complex_v = is_complex<T>::value;
81
82//==============================================================================
83
84/*!
85 @brief Inner product of two std::arrays: sum_i a_i * b_i.
86 @details
87
88 \f[ \text{inner\_product}(a,\, b) = \sum_{i=0}^{N-1} a_i \, b_i \f]
89
90 where \f$ N = \min(\text{a.size()}, \text{b.size()}) \f$.
91 The array types may differ (T and U), but U must be convertible to T.
92 The return type is T (same as the first array).
93*/
94template <typename T, typename U, std::size_t N, std::size_t M>
95constexpr T inner_product(const std::array<T, N> &a,
96 const std::array<U, M> &b) {
97 static_assert(std::is_convertible_v<U, T>,
98 "In inner_product, type of second array (U) must be "
99 "convertable to that of dirst (T)");
100 constexpr std::size_t Size = std::min(N, M);
101
102 if constexpr (Size == 0) {
103 return T{0};
104 } else if constexpr (!std::is_same_v<T, U> && is_complex_v<T>) {
105 // This is to avoid float conversion warning in case that U=double,
106 // T=complex<float>; want to case b to float, then to complex<float>
107 T res = a[0] * static_cast<typename T::value_type>(b[0]);
108 for (std::size_t i = 1; i < Size; ++i) {
109 res += a[i] * static_cast<typename T::value_type>(b[i]);
110 }
111 return res;
112 } else {
113 T res = a[0] * static_cast<T>(b[0]);
114 for (std::size_t i = 1; i < Size; ++i) {
115 res += a[i] * static_cast<T>(b[i]);
116 }
117 return res;
118 }
119}
120
121//==============================================================================
122
123namespace helper {
124
125// Simple struct for storing "Raw" Adams ("B") coefficients.
126/*
127Stored as integers, with 'denominator' factored out. \n
128Converted to double ("A" coefs) once, at compile time (see below). \n
129Adams coefficients, a_k, defined such that: \n
130 \f[ F_{n+K} = F_{n+K-1} + dx * \sum_{k=0}^K a_k y_{n+k} \f]
131where:
132 \f[ y = d(F)/dr \f]
133Note: the 'order' of the coefs is reversed compared to some sources.
134The final coefficient is separated, such that: \n
135 \f[ a_k = b_k / denom \f]
136for k = {0,1,...,K-1} \n
137and
138 \f[ a_K = b_K / denom \f]
139*/
140template <std::size_t K>
141struct AdamsB {
142 long denom;
143 std::array<long, K> bk;
144 long bK;
145};
146
147// List of Adams coefficients data
148/*
149Note: there is (of course) no 0-step Adams method.
150The entry at [0] is invalid, and will produce 0.
151It is included so that the index of this list matches the order of the method.
152Program will not compile (static_asser) is [0] is requested.
153Note: assumes that the kth element corresponds to k-order AM method.
154*/
155static constexpr auto ADAMS_data = std::tuple{
156 AdamsB<0>{1, {}, 0}, // invalid entry, but want index to match order
157 AdamsB<1>{2, {1}, 1},
158 AdamsB<2>{12, {-1, 8}, 5},
159 AdamsB<3>{24, {1, -5, 19}, 9},
160 AdamsB<4>{720, {-19, 106, -264, 646}, 251},
161 AdamsB<5>{1440, {27, -173, 482, -798, 1427}, 475},
162 AdamsB<6>{60480, {-863, 6312, -20211, 37504, -46461, 65112}, 19087},
163 AdamsB<7>{
164 120960, {1375, -11351, 41499, -88547, 123133, -121797, 139849}, 36799},
165 AdamsB<8>{
166 3628800,
167 {-33953, 312874, -1291214, 3146338, -5033120, 5595358, -4604594, 4467094},
168 1070017},
169 AdamsB<9>{7257600,
170 {57281, -583435, 2687864, -7394032, 13510082, -17283646, 16002320,
171 -11271304, 9449717},
172 2082753},
173 AdamsB<10>{479001600,
174 {-3250433, 36284876, -184776195, 567450984, -1170597042,
175 1710774528, -1823311566, 1446205080, -890175549, 656185652},
176 134211265},
177 AdamsB<11>{958003200,
178 {5675265, -68928781, 384709327, -1305971115, 3007739418,
179 -4963166514, 6043521486, -5519460582, 3828828885, -2092490673,
180 1374799219},
181 262747265},
182 AdamsB<12>{2615348736000,
183 {-13695779093, 179842822566, -1092096992268, 4063327863170,
184 -10344711794985, 19058185652796, -26204344465152, 27345870698436,
185 -21847538039895, 13465774256510, -6616420957428, 3917551216986},
186 703604254357}};
187
188} // namespace helper
189
190//==============================================================================
191
192//! Stores maximum K (order of AM method) for which we have coefficients
193//! implemented.
194static constexpr std::size_t K_max =
195 std::tuple_size_v<decltype(helper::ADAMS_data)> - 1;
196
197//==============================================================================
198
199/*!
200 @brief Holds the K+1 Adams-Moulton coefficients for the K-step AM method.
201 @details
202 The Adams coefficients a_k are defined such that:
203
204 \f[ F_{n+K} = F_{n+K-1} + dx \sum_{k=0}^{K} a_k y_{n+k}, \quad y \equiv \frac{dF}{dr} \f]
205
206 The order of the coefficients is reversed compared to some sources.
207 The final coefficient a_K is stored separately:
208
209 \f[ a_k = b_k / \text{denom}, \quad k = 0, 1, \ldots, K-1 \f]
210 \f[ a_K = b_K / \text{denom} \f]
211
212 All coefficients are stored as doubles regardless of other template parameters.
213*/
214template <std::size_t K, typename = std::enable_if_t<(K > 0)>,
215 typename = std::enable_if_t<(K <= K_max)>>
216struct AM_Coefs {
217
218private:
219 // Forms the (double) ak coefficients, from the raw (int) bk ones
220 static constexpr std::array<double, K> make_ak() {
221 const auto &am = std::get<K>(helper::ADAMS_data);
222 static_assert(
223 am.bk.size() == K,
224 "Kth Entry in ADAMS_data must correspond to K-order AM method");
225 std::array<double, K> tak{};
226 for (std::size_t i = 0; i < K; ++i) {
227 tak.at(i) = double(am.bk.at(i)) / double(am.denom);
228 }
229 return tak;
230 }
231
232 // Forms the final (double) aK coefficient, from the raw (int) bK one
233 static constexpr double make_aK() {
234 const auto &am = std::get<K>(helper::ADAMS_data);
235 return (double(am.bK) / double(am.denom));
236 }
237
238public:
239 //! First K coefficients: ak for k={0,1,...,K-1}
240 static constexpr std::array<double, K> ak{make_ak()};
241 //! Final aK coefficients: ak for k=K
242 static constexpr double aK{make_aK()};
243};
244
245//==============================================================================
246/*!
247 @brief Solves a 2x2 system of ODEs using a K-step Adams-Moulton method.
248 @details
249 The system of ODEs is defined such that:
250
251\f[ \frac{dF(t)}{dt} = D(t) * F(t) + S(t) \f]
252
253Where F is a 2D set of functions:
254
255\f[
256 F(t) = \begin{pmatrix}
257 f(t)\\
258 g(t)
259 \end{pmatrix},
260\f]
261
262D is the 2x2 "derivative matrix":
263
264\f[
265 D(t) = \begin{pmatrix}
266 a(t) & b(t)\\
267 c(t) & d(t)
268 \end{pmatrix},
269\f]
270
271and S(t) is the (optional) 2D inhomogenous term:
272
273\f[
274 S(t) = \begin{pmatrix}
275 s_f(t)\\
276 s_g(t)
277 \end{pmatrix}.
278\f]
279
280See struct `DerivativeMatrix` - which is a pure virtual struct that must be
281implmented by the user in order to define the ODE.
282
283The step-size, dt, must remain constant (since it must remain consistant
284between the K+1 and previous K points). It may be positive or negative,
285however (or complex). It's perfectly possible to have a non-uniform grid - this
286just introduces a Jacobian into the Derivative matrix; dt must still be
287constant.
288
289The template parameter, T, is the type of the argument of the Derivative
290Matrix (i.e., type of `t`). This is often `double` or `complex<double>`, but may
291also be an index type (e.g., std::size_t) if the derivative matrix is only known
292numerically at certain grid points/stored in an array.
293
294The template parameter, Y, is the type of the function value F(t), and the
295type of dt, and the return value of the Derivative Matrix. This is often
296`double`, but may also be another floating-point type, or std::complex.
297
298The first K points of the function F, and derivative dF/dt, must be known.
299You may directly access the f,g (function) and df,dg (derivative) arrays, to
300set these points.
301
302Alternatively, you may use the provided function
303 ```cpp
304 void solve_initial_K(T t0, Y f0, Y g0);
305 ```
306which automatically sets the first K values for F (and dF), given a single
307initial value for F, f0=f(t0), fg=g(t0), by using successive N-step AM
308methods, for N={1,2,...,K-1}.
309
310For now, just a 2x2 system. In theory, simple to generalise to an N*N system,
311though requires a matrix inversion.
312
313\par
314
315**Example:** Bessel's differential equation
316
317\f[
318 x^2 \frac{d^2y}{dx^2} + x \frac{dy}{dx} + (x^2-n^2)y = 0
319\f]
320
321With y(0) = 1.0, the solutions are the Bessel functions, y(x) = J_n(x)
322
323This can be re-written as a pair of coupled first-order ODEs:
324
325\f[
326 \begin{aligned}
327 \frac{dy}{dx} &= p \\
328 \frac{dp}{dx} &= \left[\left(\frac{n}{x}\right)^2 - 1\right]y - \frac{1}{x}p
329 \end{aligned}
330\f]
331
332Putting this into the notation/form required for the solver we have:
333
334\f[
335 F(x) = \begin{pmatrix}
336 y(x)\\
337 p(x)
338 \end{pmatrix}
339\f]
340
341with the "derivative matrix":
342
343\f[
344 D(x) = \begin{pmatrix}
345 0 & 1 \\
346 \left(\frac{n}{x}\right)^2 - 1 & \frac{-1}{x}
347 \end{pmatrix}
348\f]
349
350i.e.,
351
352for n=1:
353
354```cpp
355 struct BesselDerivative : AdamsMoulton::DerivativeMatrix<double, double> {
356 double a(double) const final { return 0.0; }
357 double b(double) const final { return 1.0; }
358 double c(double t) const final { return 1.0/t/t - 1.0; }
359 double d(double t) const final { return -1.0 / t; }
360 };
361```
362
363Or, more generally (for example):
364
365```cpp
366 struct BesselDerivative : AdamsMoulton::DerivativeMatrix<double, double> {
367 int n;
368 BesselDerivative(int tn) : n(tn) {}
369 double a(double) const final { return 0.0; }
370 double b(double) const final { return 1.0; }
371 double c(double t) const final { return std::pow(n/t,2) - 1.0; }
372 double d(double t) const final { return -1.0 / t; }
373 };
374```
375
376Minimal example: -- see full examples included elsewhere
377
378```
379 // Construct the Derivative matrix (BesselDerivative defined above) with n=0
380 int n = 0;
381 BesselDerivative D{n};
382
383 // Set the step size:
384 double dt = 0.01;
385
386 // Construct the Solver, using K=6-step method:
387 AdamsMoulton::ODESolver2D<6> ode{dt, &D};
388
389 // Since 1/t appears in D, we cannot start at zero. Instead, begin at small t
390 double t0 = 1.0e-6;
391
392 // Set initial points:
393 // Note: these are *approximate*, since they are technically f(0.0)
394 double f0 = 1.0;
395 double g0 = 0.0;
396
397 // Use automatic solver for first K points:
398 ode.solve_initial_K(t0, f0, g0);
399
400 // Drive forwards another 100 steps
401 for (int i = 0; i < 100; ++i) {
402 ode.drive();
403 std::cout << ode.last_t() << " " << ode.last_f() << '\n';
404 }
405```
406
407*/
408template <std::size_t K, typename T = double, typename Y = double>
410 static_assert(K > 0, "Order (K) for Adams method must be K>0");
411 static_assert(K <= K_max,
412 "Order (K) requested for Adams method too "
413 "large. Adams coefficients are implemented up to K_max-1 only");
414 static_assert(
415 is_complex_v<Y> || std::is_floating_point_v<Y>,
416 "Template parameter Y (function values and dt) must be floating point "
417 "or complex");
418 static_assert(is_complex_v<Y> || std::is_floating_point_v<Y> ||
419 std::is_integral_v<T>,
420 "Template parameter T (derivative matrix argument) must be "
421 "floating point, complex, or integral");
422
423private:
424 // Stores the AM coefficients
425 static constexpr AM_Coefs<K> am{};
426 // step size:
427 Y m_dt;
428 // previous 't' value
429 // T m_t{T(0)}; // Should be invalid (nan), but no nan for int
430 // Pointer to the derivative matrix
431 const DerivativeMatrix<T, Y> *m_D;
432
433public:
434 /*!
435 @brief Arrays storing the previous K values of f and g.
436 @details
437 Stored in chronological order regardless of the sign of dt (i.e. whether
438 driving forwards or backwards). f[0] is the oldest value, f[K-1] is newest.
439 */
440 std::array<Y, K> f{}, g{};
441
442 //! Arrays to store the previous K values of derivatives, df and dg
443 std::array<Y, K> df{}, dg{};
444
445 //! Array to store the previous K values of t: f.at(i) = f(t.at(i))
446 std::array<T, K> t{};
447
448 Y S_scale{1.0};
449
450public:
451 /*!
452 @brief Constructs the ODE solver with a given step size and derivative matrix.
453 @details
454 The step-size dt may be positive (drive forwards), negative (drive backwards),
455 or complex. A raw pointer to D is stored internally and must not be null; it
456 must outlive the ODESolver2D instance.
457 @param dt Constant step size.
458 @param D Pointer to the derivative matrix. Must not be null.
459 */
460 ODESolver2D(Y dt, const DerivativeMatrix<T, Y> *D) : m_dt(dt), m_D(D) {
461 assert(dt != Y{0.0} && "Cannot have zero step-size in ODESolver2D");
462 assert(D != nullptr && "Cannot have null Derivative Matrix in ODESolver2D");
463 }
464
465 //! Returns the AM order (number of steps), K
466 constexpr std::size_t K_steps() const { return K; }
467
468 //! Returns most recent f value. Can also access f array directly
469 Y last_f() { return f.back(); }
470 //! Returns most recent g value. Can also access g array directly
471 Y last_g() { return g.back(); }
472 //! Returns most recent t value; last_f() := f(last_t())
473 T last_t() { return t.back(); }
474 //! Returns the step size
475 Y dt() { return m_dt; }
476
477 //! Returns derivative, df/dt(t), given f(t),g(t),t
478 Y dfdt(Y ft, Y gt, T tt) const {
479 return m_D->a(tt) * ft + m_D->b(tt) * gt + S_scale * m_D->Sf(tt);
480 }
481 //! Returns derivative, dg/dt(t), given f(t),g(t),t
482 Y dgdt(Y ft, Y gt, T tt) const {
483 return m_D->c(tt) * ft + m_D->d(tt) * gt + S_scale * m_D->Sg(tt);
484 }
485
486 /*!
487 @brief Drives the ODE system one step to t_next, given the K previous values.
488 @details
489 Assumes the system has already been solved for the K previous values
490 {t-K*dt, ..., t-dt}. The value t_next should satisfy t_next = last_t + dt.
491
492 t is passed explicitly to avoid accumulation of floating-point errors: the
493 10,000th grid point may not equal t0 + 10000*dt exactly, particularly on
494 non-linear grids. The no-argument overload drive() avoids this but should be
495 used with care.
496
497 The type T of t_next must match the type expected by DerivativeMatrix; usually
498 T=double for an analytic derivative, or T=std::size_t when the derivative is
499 stored on a discrete grid.
500 @param t_next The target t value for the new step.
501 */
502 void drive(T t_next) {
503
504 // assert (t_next = Approx[next_t(last_t())])
505
506 // Use AM method to determine new values, given previous K values:
507 const auto sf = f.back() + m_dt * (inner_product(df, am.ak) +
508 am.aK * S_scale * m_D->Sf(t_next));
509 const auto sg = g.back() + m_dt * (inner_product(dg, am.ak) +
510 am.aK * S_scale * m_D->Sg(t_next));
511
512 const auto a = m_D->a(t_next);
513 const auto b = m_D->b(t_next);
514 const auto c = m_D->c(t_next);
515 const auto d = m_D->d(t_next);
516
517 const auto a0 = m_dt * static_cast<Y>(am.aK);
518 const auto a02 = a0 * a0;
519 const auto det_inv =
520 Y{1.0} / (Y{1.0} - (a02 * (b * c - a * d) + a0 * (a + d)));
521 const auto fi = (sf - a0 * (d * sf - b * sg)) * det_inv;
522 const auto gi = (sg - a0 * (-c * sf + a * sg)) * det_inv;
523
524 // Shift values along. nb: rotate({1,2,3,4}) = {2,3,4,1}
525 // We keep track of previous K values in order to determine next value
526 std::rotate(f.begin(), f.begin() + 1, f.end());
527 std::rotate(g.begin(), g.begin() + 1, g.end());
528 std::rotate(df.begin(), df.begin() + 1, df.end());
529 std::rotate(dg.begin(), dg.begin() + 1, dg.end());
530 std::rotate(t.begin(), t.begin() + 1, t.end());
531
532 // Sets new values:
533 t.back() = t_next;
534 f.back() = fi;
535 g.back() = gi;
536 df.back() = dfdt(fi, gi, t_next);
537 dg.back() = dgdt(fi, gi, t_next);
538 }
539
540 //! Overload of drive(T t) for 'default' case, where next t is defined as
541 //! last_t + dt (for arithmetic/complex types), or last_t++/last_t-- for
542 //! integral types (grid index).
543 void drive() { drive(next_t(last_t())); }
544
545 /*!
546 @brief Sets the first K values of F (and dF) given a single initial condition.
547 @details
548 Uses successive N-step AM methods for N = {1, 2, ..., K-1}:
549 - F[0] is set from the supplied initial values.
550 - F[1] is determined from F[0] using a 1-step AM method.
551 - F[2] is determined from F[0], F[1] using a 2-step AM method.
552 - ...
553 - F[K-1] is determined from F[0], ..., F[K-2] using a (K-1)-step AM method.
554
555 @param t0 Initial value of t.
556 @param f0 Initial value f(t0).
557 @param g0 Initial value g(t0).
558 */
559 void solve_initial_K(T t0, Y f0, Y g0) {
560 t.at(0) = t0;
561 f.at(0) = f0;
562 g.at(0) = g0;
563 df.at(0) = dfdt(f0, g0, t0);
564 dg.at(0) = dgdt(f0, g0, t0);
565 first_k_i<1>(next_t(t0));
566 }
567
568private:
569 // only used in solve_initial_K(). Use this, because it works for double t,
570 // where next value is t+dt, and for integral t, where next value is t++ is
571 // driving forward (dt>0), or t-- if driving backwards (dt<0)
572 T next_t(T last_t) {
573 if constexpr (std::is_integral_v<T> && is_complex_v<Y>) {
574 return (m_dt.real() > 0.0) ? last_t + 1 : last_t - 1;
575 } else if constexpr (std::is_integral_v<T>) {
576 return (m_dt > 0.0) ? last_t + 1 : last_t - 1;
577 } else {
578 return last_t + m_dt;
579 }
580 }
581
582 // Used recursively by solve_initial_K() to find first K points
583 template <std::size_t ik>
584 void first_k_i(T t_next) {
585 if constexpr (ik >= K) {
586 (void)t_next; // suppress unused variable warning on old g++ versions
587 return;
588 } else {
589 constexpr AM_Coefs<ik> ai{};
590 // nb: ai.ak is smaller than df; inner_product still works
591 const auto sf = f.at(ik - 1) + m_dt * (inner_product(df, ai.ak) +
592 am.aK * S_scale * m_D->Sf(t_next));
593 const auto sg = g.at(ik - 1) + m_dt * (inner_product(dg, ai.ak) +
594 am.aK * S_scale * m_D->Sg(t_next));
595 const auto a0 = m_dt * static_cast<Y>(ai.aK);
596 const auto a02 = a0 * a0;
597 const auto a = m_D->a(t_next);
598 const auto b = m_D->b(t_next);
599 const auto c = m_D->c(t_next);
600 const auto d = m_D->d(t_next);
601 const auto det_inv =
602 Y{1.0} / (Y{1.0} - (a02 * (b * c - a * d) + a0 * (a + d)));
603 const auto fi = (sf - a0 * (d * sf - b * sg)) * det_inv;
604 const auto gi = (sg - a0 * (-c * sf + a * sg)) * det_inv;
605 // Sets new values:
606 t.at(ik) = t_next;
607 f.at(ik) = fi;
608 g.at(ik) = gi;
609 df.at(ik) = dfdt(fi, gi, t_next);
610 dg.at(ik) = dgdt(fi, gi, t_next);
611 // call recursively
612 first_k_i<ik + 1>(next_t(t_next));
613 }
614 }
615};
616} // namespace AdamsMoulton
Solves a 2x2 system of ODEs using a K-step Adams-Moulton method.
Definition AdamsMoulton.hpp:409
T last_t()
Returns most recent t value; last_f() := f(last_t())
Definition AdamsMoulton.hpp:473
Y last_g()
Returns most recent g value. Can also access g array directly.
Definition AdamsMoulton.hpp:471
std::array< T, K > t
Array to store the previous K values of t: f.at(i) = f(t.at(i))
Definition AdamsMoulton.hpp:446
constexpr std::size_t K_steps() const
Returns the AM order (number of steps), K.
Definition AdamsMoulton.hpp:466
void drive(T t_next)
Drives the ODE system one step to t_next, given the K previous values.
Definition AdamsMoulton.hpp:502
Y last_f()
Returns most recent f value. Can also access f array directly.
Definition AdamsMoulton.hpp:469
void solve_initial_K(T t0, Y f0, Y g0)
Sets the first K values of F (and dF) given a single initial condition.
Definition AdamsMoulton.hpp:559
Y dgdt(Y ft, Y gt, T tt) const
Returns derivative, dg/dt(t), given f(t),g(t),t.
Definition AdamsMoulton.hpp:482
Y dfdt(Y ft, Y gt, T tt) const
Returns derivative, df/dt(t), given f(t),g(t),t.
Definition AdamsMoulton.hpp:478
void drive()
Overload of drive(T t) for 'default' case, where next t is defined as last_t + dt (for arithmetic/com...
Definition AdamsMoulton.hpp:543
ODESolver2D(Y dt, const DerivativeMatrix< T, Y > *D)
Constructs the ODE solver with a given step size and derivative matrix.
Definition AdamsMoulton.hpp:460
std::array< Y, K > df
Arrays to store the previous K values of derivatives, df and dg.
Definition AdamsMoulton.hpp:443
std::array< Y, K > f
Arrays storing the previous K values of f and g.
Definition AdamsMoulton.hpp:440
Y dt()
Returns the step size.
Definition AdamsMoulton.hpp:475
Contains classes and functions which use general N-step Adams Moulton method to solve systems of 2x2 ...
Definition AdamsMoulton.hpp:11
constexpr bool is_complex_v
Type trait: true if T is std::complex<U> for some U, false otherwise.
Definition AdamsMoulton.hpp:80
constexpr T inner_product(const std::array< T, N > &a, const std::array< U, M > &b)
Inner product of two std::arrays: sum_i a_i * b_i.
Definition AdamsMoulton.hpp:95
constexpr double c
speed of light in a.u. (=1/alpha)
Definition PhysConst_constants.hpp:63
Holds the K+1 Adams-Moulton coefficients for the K-step AM method.
Definition AdamsMoulton.hpp:216
static constexpr std::array< double, K > ak
First K coefficients: ak for k={0,1,...,K-1}.
Definition AdamsMoulton.hpp:240
static constexpr double aK
Final aK coefficients: ak for k=K.
Definition AdamsMoulton.hpp:242
Pure-virtual struct defining the derivative matrix for a 2x2 ODE system.
Definition AdamsMoulton.hpp:47
virtual Y a(T t) const =0
a,b,c,d are derivative matrix functions; all must be user implemented
virtual Y Sf(T) const
Sf and Sg are optional inhomogenous terms.
Definition AdamsMoulton.hpp:54