1 #include <../../nrnconf.h>
16 #include "ida/ida_impl.h"
41 extern booleantype
IDAEwtSet(IDAMem IDA_mem, N_Vector ycur);
45 #define nt_dt nrn_threads->_dt
46 #define nt_t nrn_threads->_t
52 static int res_gvardt(realtype
t, N_Vector y, N_Vector yp, N_Vector delta,
void* rdata);
54 static int minit(IDAMem);
56 static int msetup(IDAMem mem,
64 static int msolve(IDAMem mem, N_Vector b, N_Vector ycur, N_Vector ypcur, N_Vector deltacur);
66 static int mfree(IDAMem);
70 #define thread_t nrn_thread_t
93 static int res_gvardt(realtype
t, N_Vector y, N_Vector yp, N_Vector delta,
void* rdata) {
111 static int msetup(IDAMem mem, N_Vector y, N_Vector yp, N_Vector, N_Vector, N_Vector, N_Vector) {
128 static int msolve(IDAMem mem, N_Vector b, N_Vector w, N_Vector ycur, N_Vector, N_Vector) {
158 IDAFree((IDAMem)
mem_);
171 IDAMem mem = (IDAMem) IDACreate();
175 IDASetRdata(mem,
cv_);
178 mem->ida_linit =
minit;
181 mem->ida_lfree =
mfree;
182 mem->ida_setupNonNull =
false;
211 double norm = N_VWrmsNorm(ida->
delta_, ((IDAMem) (ida->
mem_))->ida_ewt);
212 Printf(
"ida check t=%.15g norm=%g\n",
t, norm);
215 printf(
" %3d %22.15g %22.15g %22.15g\n",
i,
216 N_VGetArrayPointer(ida->
cv_->
y_)[
i],
217 N_VGetArrayPointer(ida->
yp_)[
i],
218 N_VGetArrayPointer(ida->
delta_)[
i]);
228 printf(
"Daspk_init t_=%20.12g t-t_=%g t0_-t_=%g\n",
239 double dtinv = 1. /
dteps_;
270 N_VScale(dtinv,
yp_,
yp_);
287 double norm = N_VWrmsNorm(
parasite_, ((IDAMem)
mem_)->ida_ewt);
292 Printf(
"IDA initialization failure, weighted norm of residual=%g\n", norm);
296 Printf(
"IDA initialization warning, weighted norm of residual=%g\n", norm);
299 Printf(
"IDA initialization warning, weighted norm of residual=%g\n", norm);
302 Printf(
" subtracting (for next 1e-6 ms): f(y', y, %g)*exp(-1e7*(t-%g))\n",
nt_t,
nt_t);
327 IDASetStopTime(
mem_, tstop);
334 if (ier > 0 && t < cv_->t_) {
353 IDASetStopTime(
mem_, tt);
356 Printf(
"DASPK interpolate error\n");
371 printf(
"rwork size = %d\n", iwork_[18-1]);
372 printf(
"iwork size = %d\n", iwork_[17-1]);
373 printf(
"Number of time steps = %d\n", iwork_[11-1]);
374 printf(
"Number of residual evaluations = %d\n", iwork_[12-1]);
375 printf(
"Number of Jac evaluations = %d\n", iwork_[13-1]);
376 printf(
"Number of preconditioner solves = %d\n", iwork_[21-1]);
377 printf(
"Number of nonlinear iterations = %d\n", iwork_[19-1]);
378 printf(
"Number of linear iterations = %d\n", iwork_[20-1]);
379 double avlin = double(iwork_[20-1])/double(iwork_[19-1]);
380 printf(
"Average Krylov subspace dimension = %g\n", avlin);
381 printf(
"nonlinear conv. failures = %d\n", iwork_[15-1]);
382 printf(
"linear conv. failures = %d\n", iwork_[16-1]);
411 for (
i = 0;
i <
n; ++
i) {
433 for (
i = 0;
i <
n; ++
i) {
459 printf(
"Cvode::res enter tt=%g\n", tt);
461 printf(
" %d %g %g %g\n",
i, y[
i], yprime[
i], delta[
i]);
485 printf(
"Cvode::res after ode and gather_ydot into delta\n");
487 printf(
" %d %g %g %g\n",
i, y[
i], yprime[
i], delta[
i]);
503 for (
i = 0;
i <
n; ++
i) {
504 double* cd = ml->
data[
i];
510 cdvm = 1
e-3 * cd[0] * (yprime[
j] - yprime[
j + 1]);
512 delta[
j + 1] += cdvm;
521 cdvm = 1
e-3 * cd[0] * yprime[
j];
536 for (
i = 0;
i <
n; ++
i) {
537 double* cd = ml->
data[
i];
550 delta[
j] -= 1
e-3 * cd[2] * yprime[
j];
556 delta[jj] -= 1
e-3 * cd[2 *
nlayer +
k] * (yprime[jj]);
561 x = 1
e-3 * cd[2 *
nlayer +
k] * (yprime[jj] - yprime[jj + 1]);
574 delta[
i] -= yprime[
i];
584 delta[
i] -= tps[
i] * fac;
589 printf(
"Cvode::res exit res_=%d tt=%20.12g\n", res_, tt);
591 printf(
" %d %g %g %g\n",
i, y[
i], yprime[
i], delta[
i]);
598 e += delta[
i]*delta[
i];
600 printf(
"Cvode::res %d e=%g t=%.15g\n", res_,
e, tt);
634 printf(
"before nrn_solve matrix cj=%g\n", cj);
636 printf(
"before nrn_solve actual_rhs=\n");
638 printf(
"%d %g\n",
i+1, actual_rhs[
i+1]);
645 printf(
"after nrn_solve actual_rhs=\n");
646 for (
i=0;
i < neq_v_; ++
i) {
647 printf(
"%d %g\n",
i+1, actual_rhs[
i+1]);
670 return ((IDAMem)
mem_)->ida_ewt;
674 return ((IDAMem)
mem_)->ida_delta;
void play_continuous_thread(double t, NrnThread *)
int psol(double, double *, double *, double, NrnThread *)
void daspk_scatter_y(N_Vector)
void play_continuous(double t)
void before_after(BAMechList *, NrnThread *)
int res(double, double *, double *, double *, NrnThread *)
void scatter_y(double *, int)
double * n_vector_data(N_Vector, int)
void solvemem(NrnThread *)
void daspk_gather_y(N_Vector)
void scatter_ydot(double *, int)
void gather_ydot(N_Vector)
BAMechList * after_solve_
static int init_try_again_
int interpolate(double tout)
static int first_try_init_failures_
int advance_tn(double tstop)
static int init_failure_style_
static bool eq(double x, double y, double e)
static double eps(double x)
void hoc_execerror(const char *, const char *)
void nrn_multithread_job(void *(*job)(NrnThread *))
void nrndae_dkpsol(double)
void nrndae_dkres(double *, double *, double *)
void nrn_rhs(NrnThread *)
static int res_gvardt(realtype t, N_Vector y, N_Vector yp, N_Vector delta, void *rdata)
static double check(double t, Daspk *ida)
static void daspk_nrn_solve(NrnThread *nt)
void nrn_daspk_init_step(double, double, int)
static N_Vector nvec_delta
static void * msolve_thread(NrnThread *nt)
void nrn_solve(NrnThread *)
static int msolve(IDAMem mem, N_Vector b, N_Vector ycur, N_Vector ypcur, N_Vector deltacur)
void nrn_lhs(NrnThread *)
booleantype IDAEwtSet(IDAMem IDA_mem, N_Vector ycur)
static void * do_ode_thread(NrnThread *nt)
static void * daspk_gather_thread(NrnThread *nt)
static void * daspk_scatter_thread(NrnThread *nt)
static void * res_thread(NrnThread *nt)
static int msetup(IDAMem mem, N_Vector y, N_Vector ydot, N_Vector delta, N_Vector tempv1, N_Vector tempv2, N_Vector tempv3)
int const size_t const size_t n
static philox4x32_key_t k
Represent main neuron object computed by single thread.
_nrn_Fast_Imem * _nrn_fast_imem