Go to the documentation of this file.
4 #ifndef LSODA_FUNCTIONS_H
5 #define LSODA_FUNCTIONS_H
6 #define ETA 2.2204460492503131e-16
10 if(value !=0)
iopt = 1;
15 if(value !=0)
iopt = 1;
26 if(value!=0)
iopt = 1;
31 if(value!=0)
iopt = 1;
36 return (std::abs(a) < std::abs(b));
40 size_t LSODA::idamax1(
const vector<double> &dx,
const size_t n,
const size_t offset = 0)
43 size_t v = 0, vmax = 0;
45 for (
size_t i = 1; i <=
n; ++i)
47 v = abs(dx[i + offset]);
64 const size_t offset = 0)
66 std::transform(dx.begin() + 1 + offset, dx.end(), dx.begin() + 1 + offset,
67 [&da](
double x) ->
double {
return da * x; });
71 double LSODA::ddot1(
const vector<double> &a,
const vector<double> &b,
72 const size_t n,
const size_t offsetA = 0,
73 const size_t offsetB = 0)
76 for (
size_t i = 1; i <=
n; ++i)
77 sum += a[i + offsetA] * b[i + offsetB];
82 vector<double> &dy,
const size_t n,
const size_t offsetX = 0,
83 const size_t offsetY = 0)
86 for (
size_t i = 1; i <=
n; ++i)
87 dy[i + offsetY] = da * dx[i + offsetX] + dy[i + offsetY];
92 vector<int> &ipvt, vector<double> &b,
const size_t job)
105 for (k = 1; k <=
n; ++k)
107 t =
ddot1(a[k], b, k - 1);
108 b[k] = (b[k] - t) / a[k][k];
113 for (k =
n - 1; k >= 1; k--)
115 b[k] = b[k] +
ddot1(a[k], b,
n - k, k, k);
131 for (k = 1; k <=
n - 1; ++k)
140 daxpy1(t, a[k], b,
n - k, k, k);
145 for (k =
n; k >= 1; k--)
147 b[k] = b[k] / a[k][k];
149 daxpy1(t, a[k], b, k - 1);
154 void LSODA::dgefa1(vector<vector<double>> &a,
const size_t n, vector<int> &ipvt,
157 size_t j = 0, k = 0, i = 0;
163 for (k = 1; k <=
n - 1; ++k)
169 j =
idamax1(a[k],
n - k + 1, k - 1) + k - 1;
197 for (i = k + 1; i <=
n; ++i)
205 daxpy1(t, a[k], a[i],
n - k, k, k);
221 REprintf(
"[lsoda] repeated occurrence of illegal input. run aborted;\n");
222 REprintf(
" apparent infinite loop.\n");
234 for (
size_t i = 1; i <=
n; ++i)
249 double tcrit,
int *istate)
251 for (
size_t i = 1; i <=
n; ++i)
275 assert(neq + 1 == y.size());
277 size_t corflag = 0, orderflag = 0;
278 size_t i = 0, i1 = 0, j = 0, m = 0, ncf = 0;
279 double del = 0.0, delp = 0.0, dsm = 0.0, dup = 0.0, exup = 0.0, r = 0.0,
280 rh = 0.0, rhup = 0.0, told = 0.0;
281 double pdh = 0.0, pnorm = 0.0;
353 for (i = 1; i <= 5; ++i)
356 for (i = 1; i <= 12; ++i)
418 for (
size_t j =
nq; j >= 1; j--)
419 for (
size_t i1 = j; i1 <=
nq; ++i1)
420 for (i = 1; i <=
n; ++i)
421 yh_[i1][i] +=
yh_[i1 + 1][i];
424 correction(neq, y, f, &corflag, pnorm, &del, &delp, &told, &ncf, &rh, &m,
430 rh = max(rh,
hmin / fabs(
h_));
474 for (
size_t j = 1; j <=
l; ++j)
477 for (i = 1; i <=
n; ++i)
486 rh = max(rh,
hmin / fabs(
h_));
502 for (i = 1; i <=
n; ++i)
505 exup = 1. / (double)(
l + 1);
506 rhup = 1. / (1.4 * pow(dup, exup) + 0.0000014);
524 rh = max(rh,
hmin / fabs(
h_));
536 rh = max(rh,
hmin / fabs(
h_));
549 for (
size_t i = 1; i <=
n; ++i)
567 for (j =
nq; j >= 1; j--)
569 for (i1 = j; i1 <=
nq; ++i1)
570 for (i = 1; i <=
n; ++i)
571 yh_[i1][i] -=
yh_[i1 + 1][i];
574 if (fabs(
h_) <=
hmin * 1.00001)
585 if (orderflag == 1 || orderflag == 0)
589 rh = max(rh,
hmin / fabs(
h_));
595 rh = max(rh,
hmin / fabs(
h_));
622 rh = max(
hmin / fabs(
h_), rh);
624 for (i = 1; i <=
n; ++i)
626 (*f)(
tn_, &y[1], &
savf[1], _data);
628 for (i = 1; i <=
n; ++i)
650 for (
size_t i = 1; i <=
n; ++i)
654 for (
size_t i = 1; i <=
n; ++i)
658 for (
size_t i = 1; i <=
n; ++i)
662 for (
size_t i = 1; i <=
n; ++i)
697 if (k < 0 || k > (
int)
nq)
700 REprintf(
"[intdy] k = %d illegal.\n", k);
705 if ((t - tp) * (t -
tn_) > 0.)
710 REprintf(
"[intdy] t = %g illegal. t not in interval tcur - hu to tcur.\n\n",t);
716 for (
size_t jj =
l - k; jj <=
nq; ++jj)
719 for (
size_t i = 1; i <=
n; ++i)
720 dky[i] = c *
yh_[
l][i];
722 for (
int j =
nq - 1; j >= k; j--)
726 for (
int jj = jp1 - k; jj <= j; ++jj)
730 for (
size_t i = 1; i <=
n; ++i)
731 dky[i] = c *
yh_[jp1][i] + s * dky[i];
735 r = pow(
h_, (
double)(-k));
737 for (
size_t i = 1; i <=
n; ++i)
744 int i,
nq, nqm1, nqp1;
745 double agamq, fnq, fnqm1, pc[13], pint, ragq, rqfac, rq1fac, tsign, xpin;
788 for (
nq = 2;
nq <= 12;
nq++)
798 rqfac = rqfac / (double)
nq;
800 fnqm1 = (double)nqm1;
806 for (i =
nq; i >= 2; i--)
807 pc[i] = pc[i - 1] + fnqm1 * pc[i];
808 pc[1] = fnqm1 * pc[1];
815 for (i = 2; i <=
nq; ++i)
818 pint += tsign * pc[i] / (double)i;
819 xpin += tsign * pc[i] / (double)(i + 1);
824 elco[
nq][1] = pint * rq1fac;
826 for (i = 2; i <=
nq; ++i)
827 elco[
nq][i + 1] = rq1fac * pc[i] / (
double)i;
828 agamq = rqfac * xpin;
832 tesco[nqp1][1] = ragq * rqfac / (double)nqp1;
833 tesco[nqm1][3] = ragq;
847 for (
nq = 1;
nq <= 5; ++
nq)
855 for (i =
nq + 1; i >= 2; i--)
856 pc[i] = pc[i - 1] + fnq * pc[i];
861 for (i = 1; i <= nqp1; ++i)
862 elco[
nq][i] = pc[i] / pc[2];
882 *rh = min(*rh,
rmax);
883 *rh = *rh / max(1., fabs(
h_) *
hmxi * *rh);
892 *pdh = max(fabs(
h_) *
pdlast, 0.000001);
893 if ((*rh * *pdh * 1.00001) >=
sm1[
nq])
895 *rh =
sm1[
nq] / *pdh;
900 for (
size_t j = 2; j <=
l; ++j)
903 for (
size_t i = 1; i <=
n; ++i)
915 size_t i = 0, ier = 0, j = 0;
916 double fac = 0.0, hl0 = 0.0, r = 0.0, r0 = 0.0, yj = 0.0;
938 REprintf(
"[prja] miter != 2.\n");
944 r0 = 1000. * fabs(
h_) *
ETA * ((double)
n) * fac;
947 for (j = 1; j <=
n; ++j)
953 (*f)(
tn_, &y[1], &
acor[1], _data);
954 for (i = 1; i <=
n; ++i)
966 for (i = 1; i <=
n; ++i)
986 const vector<double> &w)
989 for (
size_t i = 1; i <=
n; ++i)
990 vm = max(vm, fabs(v[i]) * w[i]);
995 const vector<double> &w)
1006 double an = 0, sum = 0;
1008 for (
size_t i = 1; i <= (size_t)
n; ++i)
1011 for (
size_t j = 1; j <= (size_t)
n; ++j)
1012 sum += fabs(a[i][j]) / w[j];
1013 an = max(an, sum * w[i]);
1025 double *del,
double *delp,
double *told,
size_t *ncf,
1026 double *rh,
size_t *m,
dtype _data)
1028 double rm = 0.0, rate = 0.0, dcon = 0.0;
1041 for (
size_t i = 1; i <=
n; ++i)
1044 (*f)(
tn_, &y[1], &
savf[1], _data);
1058 prja(neq, y, f, _data);
1069 for (
size_t i = 1; i <=
n; ++i)
1078 for (
size_t i = 1; i <=
n; ++i)
1084 for (
size_t i = 1; i <=
n; ++i)
1098 for (
size_t i = 1; i <=
n; ++i)
1104 for (
size_t i = 1; i <=
n; ++i)
1122 if (*del <= 100. * pnorm *
ETA)
1124 if (*m != 0 ||
meth_ != 1)
1129 if (*del <= (1024. * *delp))
1131 rate = max(rate, rm);
1151 if (*m ==
maxcor || (*m >= 2 && *del > 2. * *delp))
1165 for (
size_t i = 1; i <=
n; ++i)
1168 (*f)(
tn_, &y[1], &
savf[1], _data);
1178 (*f)(
tn_, &y[1], &
savf[1], _data);
1189 for (
size_t j =
nq; j >= 1; j--)
1190 for (
size_t i1 = j; i1 <=
nq; ++i1)
1191 for (
size_t i = 1; i <=
n; ++i)
1192 yh_[i1][i] -=
yh_[i1 + 1][i];
1219 REprintf(
"solsy -- miter != 2.\n");
1229 int lm1, lm1p1, lm2, lm2p1, nqm1, nqm2;
1230 double rh1, rh2, rh1it, exm2, dm2, exm1, dm1, alpha, exsm;
1254 if (dsm <= (100. * pnorm *
ETA) ||
pdest == 0.)
1263 exsm = 1. / (double)
l;
1264 rh1 = 1. / (1.2 * pow(dsm, exsm) + 0.0000012);
1267 if ((*pdh * rh1) > 0.00001)
1268 rh1it =
sm1[
nq] / *pdh;
1269 rh1 = min(rh1, rh1it);
1274 exm2 = 1. / (double)lm2;
1277 rh2 = 1. / (1.2 * pow(dm2, exm2) + 0.0000012);
1282 rh2 = 1. / (1.2 * pow(dm2, exsm) + 0.0000012);
1285 if (rh2 <
ratio * rh1)
1311 exsm = 1. / (double)
l;
1316 exm1 = 1. / (double)lm1;
1319 rh1 = 1. / (1.2 * pow(dm1, exm1) + 0.0000012);
1324 rh1 = 1. / (1.2 * pow(dm1, exsm) + 0.0000012);
1330 if ((*pdh * rh1) > 0.00001)
1331 rh1it =
sm1[nqm1] / *pdh;
1332 rh1 = min(rh1, rh1it);
1333 rh2 = 1. / (1.2 * pow(dsm, exsm) + 0.0000012);
1334 if ((rh1 *
ratio) < (5. * rh2))
1336 alpha = max(0.001, rh1);
1337 dm1 *= pow(alpha, exm1);
1338 if (dm1 <= 1000. *
ETA * pnorm)
1359 for (
size_t i = 1; i <=
n; ++i)
1382 double exsm, rhdn, rhsm, ddn, exdn, r;
1386 exsm = 1. / (double)
l;
1387 rhsm = 1. / (1.2 * pow(dsm, exsm) + 0.0000012);
1393 exdn = 1. / (double)
nq;
1394 rhdn = 1. / (1.3 * pow(ddn, exdn) + 0.0000013);
1401 *pdh = max(fabs(
h_) *
pdlast, 0.000001);
1403 *rhup = min(*rhup,
sm1[
l] / *pdh);
1404 rhsm = min(rhsm,
sm1[
nq] / *pdh);
1406 rhdn = min(rhdn,
sm1[
nq - 1] / *pdh);
1420 if (kflag < 0 && *rh > 1.)
1430 if (kflag < 0 && *rh > 1.)
1438 r =
el[
l] / (double)
l;
1441 for (
size_t i = 1; i <=
n; ++i)
1459 if ((*rh * *pdh * 1.00001) <
sm1[newq])
1460 if (
kflag == 0 && *rh < 1.1)
1468 if (
kflag == 0 && *rh < 1.1)
1475 *rh = min(*rh, 0.2);
1498 array<double, 14> ep1;
1501 for (
size_t i = 1; i <=
l; ++i)
1505 conit = 0.5 / (double)(
nq + 2);
void daxpy1(const double da, const vector< double > &dx, vector< double > &dy, const size_t n, const size_t offsetX, const size_t offsetY)
Definition: lsoda_functions.h:81
vector< int > ipvt
Definition: LSODA.h:159
size_t nslp
Definition: LSODA.h:177
void orderswitch(double *rhup, double dsm, double *pdh, double *rh, size_t *orderflag)
Definition: lsoda_functions.h:1378
size_t nje
Definition: LSODA.h:167
int jstart
Definition: LSODA.h:136
size_t l
Definition: LSODA.h:162
void endstoda(void)
Definition: lsoda_functions.h:1356
#define ETA
Definition: lsoda_functions.h:6
void ewset(const vector< double > &ycur)
Definition: lsoda_functions.h:645
int itol_
Definition: LSODA.h:182
int irflag
Definition: LSODA.h:179
double h_
Definition: LSODA.h:171
double sqrteta
Definition: LSODA.h:139
void ixpr_(const int value)
Definition: lsoda_functions.h:24
double el0
Definition: LSODA.h:171
size_t mxords
Definition: LSODA.h:164
size_t nst
Definition: LSODA.h:167
size_t mused
Definition: LSODA.h:164
double pdest
Definition: LSODA.h:178
double hmin
Definition: LSODA.h:172
void intdy(double t, int k, vector< double > &dky, int *iflag)
Definition: lsoda_functions.h:691
void resetcoeff(void)
Definition: lsoda_functions.h:1492
size_t msbp
Definition: LSODA.h:162
double ratio
Definition: LSODA.h:178
void correction(const size_t neq, vector< double > &y, LSODA_ODE_SYSTEM_TYPE f, size_t *corflag, double pnorm, double *del, double *delp, double *told, size_t *ncf, double *rh, size_t *m, dtype _data)
Definition: lsoda_functions.h:1023
double rc
Definition: LSODA.h:172
size_t ialth
Definition: LSODA.h:176
size_t mxordn
Definition: LSODA.h:164
array< double, 13 > sm1
Definition: LSODA.h:145
void mxhnil_(const int value)
Definition: lsoda_functions.h:29
int kflag
Definition: LSODA.h:136
void dgesl1(const vector< vector< double >> &a, const size_t n, vector< int > &ipvt, vector< double > &b, const size_t job)
Definition: lsoda_functions.h:91
void(* LSODA_ODE_SYSTEM_TYPE)(double t, double *y, double *dydt, dtype _data)
Definition: LSODA.h:35
void hmin_(const double value)
Definition: lsoda_functions.h:13
int Maxsteps
Definition: LSODA.h:55
void stoda(const size_t neq, vector< double > &y, LSODA_ODE_SYSTEM_TYPE f, dtype _data)
Definition: lsoda_functions.h:272
array< array< double, 14 >, 13 > elco
Definition: LSODA.h:151
double vmnorm(const size_t n, const vector< double > &v, const vector< double > &w)
Definition: lsoda_functions.h:985
void hmax_(const double value)
Definition: lsoda_functions.h:8
double pdlast
Definition: LSODA.h:178
array< array< double, 4 >, 13 > tesco
Definition: LSODA.h:152
vector< double > acor
Definition: LSODA.h:156
size_t maxord
Definition: LSODA.h:162
double crate
Definition: LSODA.h:174
size_t idamax1(const vector< double > &dx, const size_t n, const size_t offset)
Definition: lsoda_functions.h:40
odeproblem * dtype
Definition: lsoda_functions.h:2
void terminate2(vector< double > &y, double *t)
Definition: lsoda_functions.h:232
void maxsteps_(const int value)
Definition: lsoda_functions.h:18
array< int, 7 > iworks
Definition: LSODA.h:64
size_t ierpj
Definition: LSODA.h:161
vector< double > savf
Definition: LSODA.h:155
vector< double > ewt
Definition: LSODA.h:154
void cfode(int meth_)
Definition: lsoda_functions.h:742
size_t n
Definition: LSODA.h:167
size_t mxncf
Definition: LSODA.h:162
size_t nqu
Definition: LSODA.h:167
void successreturn(vector< double > &y, double *t, int itask, int ihit, double tcrit, int *istate)
Definition: lsoda_functions.h:248
size_t nfe
Definition: LSODA.h:167
void methodswitch(double dsm, double pnorm, double *pdh, double *rh)
Definition: lsoda_functions.h:1227
static bool abs_compare(double a, double b)
Definition: lsoda_functions.h:34
array< double, 14 > el
Definition: LSODA.h:147
size_t lmax
Definition: LSODA.h:176
void dscal1(const double da, vector< double > &dx, const size_t n, const size_t offset)
Definition: lsoda_functions.h:63
void scaleh(double *rh, double *pdh)
Definition: lsoda_functions.h:873
size_t jtyp
Definition: LSODA.h:164
vector< vector< double > > wm_
Definition: LSODA.h:158
size_t meth_
Definition: LSODA.h:165
double fnorm(int n, const vector< vector< double >> &a, const vector< double > &w)
Definition: lsoda_functions.h:994
vector< vector< double > > yh_
Definition: LSODA.h:157
array< double, 6 > cm2
Definition: LSODA.h:149
double tn_
Definition: LSODA.h:172
size_t miter
Definition: LSODA.h:162
void corfailure(double *told, double *rh, size_t *ncf, size_t *corflag)
Definition: lsoda_functions.h:1184
size_t ipup
Definition: LSODA.h:176
double rmax
Definition: LSODA.h:174
double pdnorm
Definition: LSODA.h:173
size_t iersl
Definition: LSODA.h:161
double ccmax
Definition: LSODA.h:171
double conit
Definition: LSODA.h:174
void dgefa1(vector< vector< double >> &a, const size_t n, vector< int > &ipvt, size_t *const info)
Definition: lsoda_functions.h:154
std::vector< double > rtol_
Definition: LSODA.h:183
void terminate(int *istate)
Definition: lsoda_functions.h:215
int iopt
Definition: LSODA.h:50
void prja(const size_t neq, vector< double > &y, LSODA_ODE_SYSTEM_TYPE f, dtype _data)
Definition: lsoda_functions.h:912
array< double, 13 > cm1
Definition: LSODA.h:148
std::vector< double > atol_
Definition: LSODA.h:184
double hu
Definition: LSODA.h:172
size_t jcur
Definition: LSODA.h:161
size_t maxcor
Definition: LSODA.h:162
double hmxi
Definition: LSODA.h:172
array< double, 4 > rworks
Definition: LSODA.h:65
void solsy(vector< double > &y)
Definition: lsoda_functions.h:1213
double hold
Definition: LSODA.h:174
size_t illin
Definition: LSODA.h:161
Definition: odeproblem.h:94
int icount
Definition: LSODA.h:179
int itask
Definition: LSODA.h:52
size_t nq
Definition: LSODA.h:167
double ddot1(const vector< double > &a, const vector< double > &b, const size_t n, const size_t offsetA, const size_t offsetB)
Definition: lsoda_functions.h:71