46template <
typename T =
double,
typename Y =
double>
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;
54 virtual Y
Sf(T)
const {
return Y(0); };
55 virtual Y Sg(T)
const {
return Y(0); };
56 virtual ~DerivativeMatrix() =
default;
63struct is_complex : std::false_type {};
66struct is_complex<std::complex<T>> : std::true_type {};
94template <
typename T,
typename U, std::
size_t N, std::
size_t M>
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);
102 if constexpr (Size == 0) {
104 }
else if constexpr (!std::is_same_v<T, U> && is_complex_v<T>) {
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]);
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]);
140template <std::
size_t K>
143 std::array<long, K> bk;
155static constexpr auto ADAMS_data = std::tuple{
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},
164 120960, {1375, -11351, 41499, -88547, 123133, -121797, 139849}, 36799},
167 {-33953, 312874, -1291214, 3146338, -5033120, 5595358, -4604594, 4467094},
170 {57281, -583435, 2687864, -7394032, 13510082, -17283646, 16002320,
173 AdamsB<10>{479001600,
174 {-3250433, 36284876, -184776195, 567450984, -1170597042,
175 1710774528, -1823311566, 1446205080, -890175549, 656185652},
177 AdamsB<11>{958003200,
178 {5675265, -68928781, 384709327, -1305971115, 3007739418,
179 -4963166514, 6043521486, -5519460582, 3828828885, -2092490673,
182 AdamsB<12>{2615348736000,
183 {-13695779093, 179842822566, -1092096992268, 4063327863170,
184 -10344711794985, 19058185652796, -26204344465152, 27345870698436,
185 -21847538039895, 13465774256510, -6616420957428, 3917551216986},
194static constexpr std::size_t K_max =
195 std::tuple_size_v<
decltype(helper::ADAMS_data)> - 1;
214template <std::
size_t K,
typename = std::enable_if_t<(K > 0)>,
215 typename = std::enable_if_t<(K <= K_max)>>
220 static constexpr std::array<double, K> make_ak() {
221 const auto &am = std::get<K>(helper::ADAMS_data);
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);
233 static constexpr double make_aK() {
234 const auto &am = std::get<K>(helper::ADAMS_data);
235 return (
double(am.bK) /
double(am.denom));
240 static constexpr std::array<double, K>
ak{make_ak()};
242 static constexpr double aK{make_aK()};
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");
415 is_complex_v<Y> || std::is_floating_point_v<Y>,
416 "Template parameter Y (function values and dt) must be floating point "
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");
440 std::array<Y, K>
f{}, g{};
443 std::array<Y, K>
df{}, dg{};
446 std::array<T, K>
t{};
461 assert(
dt != Y{0.0} &&
"Cannot have zero step-size in ODESolver2D");
462 assert(D !=
nullptr &&
"Cannot have null Derivative Matrix in ODESolver2D");
466 constexpr std::size_t
K_steps()
const {
return K; }
475 Y
dt() {
return m_dt; }
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);
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);
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));
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);
517 const auto a0 = m_dt *
static_cast<Y
>(am.aK);
518 const auto a02 = a0 * a0;
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;
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());
536 df.back() =
dfdt(fi, gi, t_next);
537 dg.back() =
dgdt(fi, gi, t_next);
563 df.at(0) =
dfdt(f0, g0, t0);
564 dg.at(0) =
dgdt(f0, g0, t0);
565 first_k_i<1>(next_t(t0));
573 if constexpr (std::is_integral_v<T> && is_complex_v<Y>) {
575 }
else if constexpr (std::is_integral_v<T>) {
583 template <std::
size_t ik>
584 void first_k_i(T t_next) {
585 if constexpr (ik >= K) {
589 constexpr AM_Coefs<ik> ai{};
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);
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;
609 df.at(ik) =
dfdt(fi, gi, t_next);
610 dg.at(ik) =
dgdt(fi, gi, t_next);
612 first_k_i<ik + 1>(next_t(t_next));
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