NEURON
nrndaspk.cpp
Go to the documentation of this file.
1 #include <../../nrnconf.h>
2 // differential algebraic system solver interface to DDASPK
3 
4 // DDASPK is translated from fortran with f2c. Hence all args are
5 // call by reference
6 
7 #include <stdio.h>
8 #include <math.h>
9 #include <InterViews/resource.h>
10 #include "spmatrix.h"
11 #include "nrnoc2iv.h"
12 #include "cvodeobj.h"
13 #include "nrndaspk.h"
14 #include "netcvode.h"
15 #include "ida/ida.h"
16 #include "ida/ida_impl.h"
17 #include "mymath.h"
18 
19 // the state of the g - d2/dx2 matrix for voltages
20 #define INVALID 0
21 #define NO_CAP 1
22 #define SETUP 2
23 #define FACTORED 3
24 static int solve_state_;
25 
26 // prototypes
27 
28 double Daspk::dteps_;
29 
30 
31 extern void nrndae_dkres(double*, double*, double*);
32 extern void nrndae_dkpsol(double);
33 extern void nrn_rhs(NrnThread*);
34 extern void nrn_lhs(NrnThread*);
35 extern void nrn_solve(NrnThread*);
36 void nrn_daspk_init_step(double, double, int);
37 // this is private in ida.cpp but we want to check if our initialization
38 // is good. Unfortunately ewt is set on the first call to solve which
39 // is too late for us.
40 extern "C" {
41 extern booleantype IDAEwtSet(IDAMem IDA_mem, N_Vector ycur);
42 } // extern "C"
43 
44 // extern double t, dt;
45 #define nt_dt nrn_threads->_dt
46 #define nt_t nrn_threads->_t
47 
48 static void daspk_nrn_solve(NrnThread* nt) {
49  nrn_solve(nt);
50 }
51 
52 static int res_gvardt(realtype t, N_Vector y, N_Vector yp, N_Vector delta, void* rdata);
53 
54 static int minit(IDAMem);
55 
56 static int msetup(IDAMem mem,
57  N_Vector y,
58  N_Vector ydot,
59  N_Vector delta,
60  N_Vector tempv1,
61  N_Vector tempv2,
62  N_Vector tempv3);
63 
64 static int msolve(IDAMem mem, N_Vector b, N_Vector ycur, N_Vector ypcur, N_Vector deltacur);
65 
66 static int mfree(IDAMem);
67 
68 
69 // at least in DARWIN the following is already declared so avoid conflict
70 #define thread_t nrn_thread_t
71 
72 // residual
73 static N_Vector nvec_y;
74 static N_Vector nvec_yp;
75 static N_Vector nvec_delta;
76 static double thread_t;
77 static double thread_cj;
78 static int thread_ier;
79 static Cvode* thread_cv;
80 static void* res_thread(NrnThread* nt) {
81  int i = nt->id;
82  Cvode* cv = thread_cv;
83  int ier = cv->res(thread_t,
84  cv->n_vector_data(nvec_y, i),
85  cv->n_vector_data(nvec_yp, i),
87  nt);
88  if (ier != 0) {
89  thread_ier = ier;
90  }
91  return 0;
92 }
93 static int res_gvardt(realtype t, N_Vector y, N_Vector yp, N_Vector delta, void* rdata) {
94  thread_cv = (Cvode*) rdata;
95  nvec_y = y;
96  nvec_yp = yp;
97  nvec_delta = delta;
98  thread_t = t;
99  thread_ier = 0;
101  return thread_ier;
102 }
103 
104 // linear solver specific allocation and initialization
105 static int minit(IDAMem) {
106  return IDA_SUCCESS;
107 }
108 
109 // linear solver preparation for subsequent calls to msolve
110 // approximation to jacobian. Everything necessary for solving P*x = b
111 static int msetup(IDAMem mem, N_Vector y, N_Vector yp, N_Vector, N_Vector, N_Vector, N_Vector) {
112  Cvode* cv = (Cvode*) mem->ida_rdata;
113  ++cv->jac_calls_;
114  return 0;
115 }
116 
117 /* solve P*x = b */
118 static void* msolve_thread(NrnThread* nt) {
119  int i = nt->id;
120  Cvode* cv = thread_cv;
121  int ier = cv->psol(
123  if (ier != 0) {
124  thread_ier = ier;
125  }
126  return 0;
127 }
128 static int msolve(IDAMem mem, N_Vector b, N_Vector w, N_Vector ycur, N_Vector, N_Vector) {
129  thread_cv = (Cvode*) mem->ida_rdata;
130  thread_t = mem->ida_tn;
131  nvec_y = ycur;
132  nvec_yp = b;
133  thread_cj = mem->ida_cj;
135  return thread_ier;
136 }
137 
138 static int mfree(IDAMem) {
139  return IDA_SUCCESS;
140 }
141 
142 Daspk::Daspk(Cvode* cv, int neq) {
143  // printf("Daspk::Daspk\n");
144  cv_ = cv;
145  yp_ = cv->nvnew(neq);
146  delta_ = cv->nvnew(neq);
147  parasite_ = cv->nvnew(neq);
148  use_parasite_ = false;
149  spmat_ = nil;
150  mem_ = nil;
151 }
152 
154  // printf("Daspk::~Daspk\n");
155  N_VDestroy(delta_);
156  N_VDestroy(yp_);
157  if (mem_) {
158  IDAFree((IDAMem) mem_);
159  }
160 }
161 
163  int ier;
164  if (mem_) {
165  ier = IDAReInit(
166  mem_, res_gvardt, cv_->t_, cv_->y_, yp_, IDA_SV, &cv_->ncv_->rtol_, cv_->atolnvec_);
167  if (ier < 0) {
168  hoc_execerror("IDAReInit error", 0);
169  }
170  } else {
171  IDAMem mem = (IDAMem) IDACreate();
172  if (!mem) {
173  hoc_execerror("IDAMalloc error", 0);
174  }
175  IDASetRdata(mem, cv_);
176  ier = IDAMalloc(
177  mem, res_gvardt, cv_->t_, cv_->y_, yp_, IDA_SV, &cv_->ncv_->rtol_, cv_->atolnvec_);
178  mem->ida_linit = minit;
179  mem->ida_lsetup = msetup;
180  mem->ida_lsolve = msolve;
181  mem->ida_lfree = mfree;
182  mem->ida_setupNonNull = false;
183  mem_ = mem;
184  }
185 }
186 
187 void Daspk::info() {}
188 
189 
190 // last two bits, 0 error, 1 warning, 2 apply parasitic
191 // if init_failure_style & 010, then use the original method
195 
196 static void* do_ode_thread(NrnThread* nt) {
197  int i;
198  Cvode* cv = thread_cv;
199  nt->_t = cv->t_;
200  cv->do_ode(nt);
201  CvodeThreadData& z = cv->ctd_[nt->id];
202  double* yp = cv->n_vector_data(nvec_yp, nt->id);
203  for (i = z.neq_v_; i < z.nvsize_; ++i) {
204  yp[i] = *(z.pvdot_[i]);
205  }
206  return 0;
207 }
208 
209 static double check(double t, Daspk* ida) {
210  res_gvardt(t, ida->cv_->y_, ida->yp_, ida->delta_, ida->cv_);
211  double norm = N_VWrmsNorm(ida->delta_, ((IDAMem) (ida->mem_))->ida_ewt);
212  Printf("ida check t=%.15g norm=%g\n", t, norm);
213 #if 0
214  for (int i=0; i < ida->cv_->neq_; ++i) {
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]);
219  }
220 #endif
221  return norm;
222 }
223 
224 int Daspk::init() {
225  extern double t;
226  int i;
227 #if 0
228 printf("Daspk_init t_=%20.12g t-t_=%g t0_-t_=%g\n",
229 cv_->t_, t-cv_->t_, cv_->t0_-cv_->t_);
230 #endif
231  N_VConst(0., yp_);
232 
233  // the new initial condition is based on a dteps_ step backward euler
234  // linear solution with respect to the old state in order to
235  // start the following initial condition calculation with a "valid"
236  // (in a linear system sense) initial state.
237 
238  double tt = cv_->t_;
239  double dtinv = 1. / dteps_;
240  if (init_failure_style_ & 010) {
241  cv_->play_continuous(tt);
242  nrn_daspk_init_step(tt, dteps_, 1);
243  nrn_daspk_init_step(tt, dteps_, 1);
245  cv_->play_continuous(tt);
246  nrn_daspk_init_step(tt, dteps_, 1);
248  N_VLinearSum(dtinv, cv_->y_, -dtinv, yp_, yp_);
249  } else {
250 #if 0
251  cv_->play_continuous(tt);
252  nrn_daspk_init_step(tt, dteps_, 1);
254  tt = cv_->t_ + dteps_;
255  cv_->play_continuous(tt);
256  nrn_daspk_init_step(tt, dteps_, 1);
258  N_VLinearSum(dtinv, yp_, -dtinv, cv_->y_, yp_);
260 #else
261  cv_->play_continuous(tt);
262  nrn_daspk_init_step(tt, dteps_, 1); // first approx to y (and maybe good enough)
263  nrn_daspk_init_step(tt, dteps_, 1); // 2nd approx to y (trouble with 2sramp.hoc)
264 
266  tt = cv_->t_ + dteps_;
267  cv_->play_continuous(tt);
268  nrn_daspk_init_step(tt, dteps_, 0); // rhs contains delta y (for v, vext, linmod
269  cv_->gather_ydot(yp_);
270  N_VScale(dtinv, yp_, yp_);
271 #endif
272  }
273  thread_cv = cv_;
274  nvec_yp = yp_;
276  ida_init();
277  t = cv_->t_;
278 #if 1
279  // test
280  // printf("test\n");
281  if (!IDAEwtSet((IDAMem) mem_, cv_->y_)) {
282  hoc_execerror("Bad Ida error weight vector", 0);
283  }
284  use_parasite_ = false;
285  // check(cv_->t_, this);
287  double norm = N_VWrmsNorm(parasite_, ((IDAMem) mem_)->ida_ewt);
288  // printf("norm=%g at t=%g\n", norm, t);
289  if (norm > 1.) {
290  switch (init_failure_style_ & 03) {
291  case 0:
292  Printf("IDA initialization failure, weighted norm of residual=%g\n", norm);
293  return IDA_ERR_FAIL;
294  break;
295  case 1:
296  Printf("IDA initialization warning, weighted norm of residual=%g\n", norm);
297  break;
298  case 2:
299  Printf("IDA initialization warning, weighted norm of residual=%g\n", norm);
300  use_parasite_ = true;
301  t_parasite_ = nt_t;
302  Printf(" subtracting (for next 1e-6 ms): f(y', y, %g)*exp(-1e7*(t-%g))\n", nt_t, nt_t);
303  break;
304  }
305 #if 0
306 for (i=0; i < cv_->neq_; ++i) {
307  printf(" %d %g %g %g %g\n", i, nt_t, N_VGetArrayPointer(cv_->y_)[i], N_VGetArrayPointer(yp_)[i], N_VGetArrayPointer(delta_)[i]);
308 }
309 #endif
310  if (init_try_again_ < 0) {
312  init_try_again_ += 1;
313  int err = init();
314  init_try_again_ = 0;
315  return err;
316  }
317  return 0;
318  }
319 #endif
320 
321  return 0;
322 }
323 
324 int Daspk::advance_tn(double tstop) {
325  // printf("Daspk::advance_tn(%g)\n", tstop);
326  double tn = cv_->tn_;
327  IDASetStopTime(mem_, tstop);
328  int ier = IDASolve(mem_, tstop, &cv_->t_, cv_->y_, yp_, IDA_ONE_STEP_TSTOP);
329  if (ier < 0) {
330  // printf("DASPK advance_tn error\n");
331  return ier;
332  }
333 #if 0
334  if (ier > 0 && t < cv_->t_) {
335  // interpolation to tstop does not call res. So we have to.
336  cv_->res(cv_->t_, N_VGetArrayPointer(cv_->y_), N_VGetArrayPointer(yp_), N_VGetArrayPointer(delta_));
338  }
339 #else
340  // this is very bad, performance-wise. However ida modifies its states
341  // after a call to fun with the proper t.
342  res_gvardt(cv_->t_, cv_->y_, yp_, delta_, cv_);
343 #endif
344  cv_->t0_ = tn;
345  cv_->tn_ = cv_->t_;
346  // printf("Daspk::advance_tn complete.\n");
347  return ier;
348 }
349 
350 int Daspk::interpolate(double tt) {
351  // printf("Daspk::interpolate %.15g\n", tt);
352  assert(tt >= cv_->t0_ && tt <= cv_->tn_);
353  IDASetStopTime(mem_, tt);
354  int ier = IDASolve(mem_, tt, &cv_->t_, cv_->y_, yp_, IDA_NORMAL);
355  if (ier < 0) {
356  Printf("DASPK interpolate error\n");
357  return ier;
358  }
360  // interpolation does not call res. So we have to.
361  res_gvardt(cv_->t_, cv_->y_, yp_, delta_, cv_);
362  // if(MyMath::eq(t, cv_->t_, NetCvode::eps(cv_->t_))) {
363  // printf("t=%.15g t_=%.15g\n", t, cv_->t_);
364  //}
365  // assert(MyMath::eq(t, cv_->t_, NetCvode::eps(cv_->t_)));
366  return ier;
367 }
368 
370 #if 0
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]);
383 #endif
385  Printf(" %d First try Initialization failures\n", first_try_init_failures_);
386  }
387 }
388 
389 static void* daspk_scatter_thread(NrnThread* nt) {
391  return 0;
392 }
393 void Cvode::daspk_scatter_y(N_Vector y) {
394  thread_cv = this;
395  nvec_y = y;
397 }
398 void Cvode::daspk_scatter_y(double* y, int tid) {
399  // the dependent variables in daspk are vi,vx,etc
400  // whereas in the node structure we need vm, vx
401  // note that a corresponding transformation for gather_ydot is
402  // not needed since the matrix solve is already with respect to vi,vx
403  // in all cases. (i.e. the solution vector is in the right hand side
404  // and refers to vi, vx.
405  scatter_y(y, tid);
406  // transform the vm+vext to vm
407  CvodeThreadData& z = ctd_[tid];
408  if (z.cmlext_) {
409  Memb_list* ml = z.cmlext_->ml;
410  int i, n = ml->nodecount;
411  for (i = 0; i < n; ++i) {
412  Node* nd = ml->nodelist[i];
413  NODEV(nd) -= nd->extnode->v[0];
414  }
415  }
416 }
417 static void* daspk_gather_thread(NrnThread* nt) {
419  return 0;
420 }
421 void Cvode::daspk_gather_y(N_Vector y) {
422  thread_cv = this;
423  nvec_y = y;
425 }
426 void Cvode::daspk_gather_y(double* y, int tid) {
427  gather_y(y, tid);
428  // transform vm to vm+vext
429  CvodeThreadData& z = ctd_[tid];
430  if (z.cmlext_) {
431  Memb_list* ml = z.cmlext_->ml;
432  int i, n = ml->nodecount;
433  for (i = 0; i < n; ++i) {
434  Node* nd = ml->nodelist[i];
435  int j = nd->eqn_index_;
436  y[j - 1] += y[j];
437  }
438  }
439 }
440 
441 // for res and psol the equations for c*yp = f(y) are
442 // cast in the form G(t,y,yp) = f(y) - c*yp
443 // So res calculates delta = f(y) - c*yp
444 // and psol solves (c*cj - df/dy)*x = -b
445 // Note that since cvode uses J = 1 - gam*df/dy and
446 // ida uses J = df/dy - cj*df/dyp that is the origin of the use of -b in our
447 // psol and also why all the non-voltage odes are scaled by dt at the
448 // end of it.
449 
450 int Cvode::res(double tt, double* y, double* yprime, double* delta, NrnThread* nt) {
451  CvodeThreadData& z = ctd_[nt->id];
452  ++f_calls_;
453  static int res_;
454  int i;
455  nt->_t = tt;
456  res_++;
457 
458 #if 0
459 printf("Cvode::res enter tt=%g\n", tt);
460 for (i=0; i < z.nvsize_; ++i) {
461  printf(" %d %g %g %g\n", i, y[i], yprime[i], delta[i]);
462 }
463 #endif
464  nt->_vcv = this; // some models may need to know this
465  daspk_scatter_y(y, nt->id); // vi, vext, channel states, linmod non-node y.
466  // rhs of cy' = f(y)
467  play_continuous_thread(tt, nt);
468  nrn_rhs(nt);
469  do_ode(nt);
470  // accumulate into delta
471  gather_ydot(delta, nt->id);
472 
473  // now calculate -c*yp. i.e.
474  // cm*vm' + c_linmod*vi' internal current balance
475  // cx*vx' + c_linmod*vx' external current balance
476  // c_linmod*y' non-node linmod states
477  // y' mechanism states
478 
479  // this can be accumulated into delta in several stages
480  // -cm*vm'and -cx*vx for current balance equation delta's
481  // -c_linmod*yp (but note that the node yp yp(vm)+yp(vx))
482  // subtract yp from mechanism state delta's
483 
484 #if 0
485 printf("Cvode::res after ode and gather_ydot into delta\n");
486 for (i=0; i < z.nvsize_; ++i) {
487  printf(" %d %g %g %g\n", i, y[i], yprime[i], delta[i]);
488 }
489 #endif
490  // the cap nodes : see nrnoc/capac.cpp for location of cm, etc.
491  // these are not in same order as for cvode but are in
492  // spmatrix order mixed with nocap nodes and extracellular
493  // therefore we use the Node.eqn_index to calculate the delta index.
494  // assert(use_sparse13 == true && nlayer <= 1);
495  assert(use_sparse13 == true);
496  if (z.cmlcap_) {
497  Memb_list* ml = z.cmlcap_->ml;
498  int n = ml->nodecount;
499  double* p = NULL;
500  if (nt->_nrn_fast_imem) {
502  }
503  for (i = 0; i < n; ++i) {
504  double* cd = ml->data[i];
505  Node* nd = ml->nodelist[i];
506  int j = nd->eqn_index_ - 1;
507  Extnode* nde = nd->extnode;
508  double cdvm;
509  if (nde) {
510  cdvm = 1e-3 * cd[0] * (yprime[j] - yprime[j + 1]);
511  delta[j] -= cdvm;
512  delta[j + 1] += cdvm;
513  // i_cap
514  cd[1] = cdvm;
515 #if I_MEMBRANE
516  // add i_cap to i_ion which is in sav_g
517  // this will be copied to i_membrane below
518  nde->param[3 + 3 * nlayer] += cdvm;
519 #endif
520  } else {
521  cdvm = 1e-3 * cd[0] * yprime[j];
522  delta[j] -= cdvm;
523  cd[1] = cdvm;
524  }
525  if (p) {
526  int i = nd->v_node_index;
527  p[i] += cdvm;
528  p[i] *= NODEAREA(nd) * 0.01;
529  }
530  }
531  }
532  // See nrnoc/excelln.cpp for location of cx.
533  if (z.cmlext_) {
534  Memb_list* ml = z.cmlext_->ml;
535  int n = ml->nodecount;
536  for (i = 0; i < n; ++i) {
537  double* cd = ml->data[i];
538  Node* nd = ml->nodelist[i];
539  int j = nd->eqn_index_;
540 #if EXTRACELLULAR
541 #if I_MEMBRANE
542  // i_membrane = sav_rhs --- even for zero area nodes
543  cd[1 + 3 * nlayer] = cd[3 + 3 * nlayer];
544 #endif /*I_MEMBRANE*/
545  if (nlayer == 1) {
546  // only works for one layer
547  // otherwise loop over layer,
548  // xc is (pd + 2*(nlayer))[layer]
549  // and deal with yprime[i+layer]-yprime[i+layer+1]
550  delta[j] -= 1e-3 * cd[2] * yprime[j];
551  } else {
552  int k, jj;
553  double x;
554  k = nlayer - 1;
555  jj = j + k;
556  delta[jj] -= 1e-3 * cd[2 * nlayer + k] * (yprime[jj]);
557  for (k = nlayer - 2; k >= 0; --k) {
558  // k=0 refers to stuff between layer 0 and 1
559  // j is for vext[0]
560  jj = j + k;
561  x = 1e-3 * cd[2 * nlayer + k] * (yprime[jj] - yprime[jj + 1]);
562  delta[jj] -= x;
563  delta[jj + 1] += x; // last one in iteration is nlayer-1
564  }
565  }
566 #endif /*EXTRACELLULAR*/
567  }
568  }
569 
570  nrndae_dkres(y, yprime, delta);
571 
572  // the ode's
573  for (i = z.neq_v_; i < z.nvsize_; ++i) {
574  delta[i] -= yprime[i];
575  }
576 
577  for (i = 0; i < z.nvsize_; ++i) {
578  delta[i] *= -1.;
579  }
580  if (daspk_->use_parasite_ && tt - daspk_->t_parasite_ < 1e-6) {
581  double fac = exp(1e7 * (daspk_->t_parasite_ - tt));
582  double* tps = n_vector_data(daspk_->parasite_, nt->id);
583  for (i = 0; i < z.nvsize_; ++i) {
584  delta[i] -= tps[i] * fac;
585  }
586  }
587  before_after(z.after_solve_, nt);
588 #if 0
589 printf("Cvode::res exit res_=%d tt=%20.12g\n", res_, tt);
590 for (i=0; i < z.nvsize_; ++i) {
591  printf(" %d %g %g %g\n", i, y[i], yprime[i], delta[i]);
592 }
593 #endif
594  nt->_vcv = 0;
595 #if 0
596 double e = 0;
597 for (i=0; i < z.nvsize_; ++i) {
598  e += delta[i]*delta[i];
599 }
600 printf("Cvode::res %d e=%g t=%.15g\n", res_, e, tt);
601 #endif
602  return 0;
603 }
604 
605 int Cvode::psol(double tt, double* y, double* b, double cj, NrnThread* _nt) {
606  CvodeThreadData& z = ctd_[_nt->id];
607  ++mxb_calls_;
608  int i;
609  _nt->_t = tt;
610 
611 #if 0
612 printf("Cvode::psol tt=%g solvestate=%d \n", tt, solve_state_);
613 for (i=0; i < z.nvsize_; ++i) {
614 printf(" %g", b[i]);
615 }
616 printf("\n");
617 #endif
618 
619  _nt->cj = cj;
620  _nt->_dt = 1. / cj;
621 
622  _nt->_vcv = this;
623  daspk_scatter_y(y, _nt->id); // I'm not sure this is necessary.
624  if (solve_state_ == INVALID) {
625  nrn_lhs(_nt); // designed to setup M*[dvm+dvext, dvext, dy] = ...
627  }
628  if (solve_state_ == SETUP) {
629  // if using sparse 13 then should factor
631  }
632  scatter_ydot(b, _nt->id);
633 #if 0
634 printf("before nrn_solve matrix cj=%g\n", cj);
635 spPrint(sp13mat_, 1,1,1);
636 printf("before nrn_solve actual_rhs=\n");
637 for (i=0; i < z.neq_v_; ++i) {
638  printf("%d %g\n", i+1, actual_rhs[i+1]);
639 }
640 #endif
641  daspk_nrn_solve(_nt); // not the cvode one
642 #if 0
643 //printf("after nrn_solve matrix\n");
644 //spPrint(sp13mat_, 1,1,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]);
648 }
649 #endif
650  solve_state_ = INVALID; // but not if using sparse13
651  solvemem(_nt);
652  gather_ydot(b, _nt->id);
653  // the ode's of the form m' = (minf - m)/mtau in model descriptions compute
654  // b = b/(1 + dt*mtau) since cvode required J = 1 - gam*df/dy
655  // so we need to scale those results by 1/cj.
656  for (i = z.neq_v_; i < z.nvsize_; ++i) {
657  b[i] *= _nt->_dt;
658  }
659 #if 0
660 for (i=0; i < z.nvsize_; ++i) {
661 printf(" %g", b[i]);
662 }
663 printf("\n");
664 #endif
665  _nt->_vcv = 0;
666  return 0;
667 }
668 
669 N_Vector Daspk::ewtvec() {
670  return ((IDAMem) mem_)->ida_ewt;
671 }
672 
673 N_Vector Daspk::acorvec() {
674  return ((IDAMem) mem_)->ida_delta;
675 }
#define nil
Definition: enter-scope.h:36
Memb_list * ml
Definition: cvodeobj.h:28
Definition: cvodeobj.h:76
int neq_
Definition: cvodeobj.h:219
double t_
Definition: cvodeobj.h:105
int mxb_calls_
Definition: cvodeobj.h:111
void play_continuous_thread(double t, NrnThread *)
Definition: occvode.cpp:1082
N_Vector atolnvec_
Definition: cvodeobj.h:204
int jac_calls_
Definition: cvodeobj.h:111
int psol(double, double *, double *, double, NrnThread *)
Definition: nrndaspk.cpp:605
CvodeThreadData * ctd_
Definition: cvodeobj.h:214
void daspk_scatter_y(N_Vector)
Definition: nrndaspk.cpp:393
void gather_y(N_Vector)
Definition: occvode.cpp:490
void play_continuous(double t)
Definition: occvode.cpp:1068
N_Vector y_
Definition: cvodeobj.h:203
Daspk * daspk_
Definition: cvodeobj.h:181
NetCvode * ncv_
Definition: cvodeobj.h:218
void before_after(BAMechList *, NrnThread *)
Definition: occvode.cpp:778
int res(double, double *, double *, double *, NrnThread *)
Definition: nrndaspk.cpp:450
void scatter_y(double *, int)
Definition: occvode.cpp:464
double * n_vector_data(N_Vector, int)
Definition: occvode.cpp:451
void solvemem(NrnThread *)
Definition: occvode.cpp:642
void daspk_gather_y(N_Vector)
Definition: nrndaspk.cpp:421
double t0_
Definition: cvodeobj.h:105
double tn_
Definition: cvodeobj.h:105
void scatter_ydot(double *, int)
Definition: occvode.cpp:508
N_Vector nvnew(long)
Definition: cvodeobj.cpp:811
void gather_ydot(N_Vector)
Definition: occvode.cpp:521
int f_calls_
Definition: cvodeobj.h:111
void do_ode(NrnThread *)
Definition: occvode.cpp:901
double ** pvdot_
Definition: cvodeobj.h:66
CvMembList * cmlext_
Definition: cvodeobj.h:54
BAMechList * after_solve_
Definition: cvodeobj.h:57
CvMembList * cmlcap_
Definition: cvodeobj.h:53
Definition: nrndaspk.h:11
virtual ~Daspk()
Definition: nrndaspk.cpp:153
void statistics()
Definition: nrndaspk.cpp:369
double t_parasite_
Definition: nrndaspk.h:32
bool use_parasite_
Definition: nrndaspk.h:33
void ida_init()
Definition: nrndaspk.cpp:162
int init()
Definition: nrndaspk.cpp:224
static double dteps_
Definition: nrndaspk.h:36
static int init_try_again_
Definition: nrndaspk.h:37
int interpolate(double tout)
Definition: nrndaspk.cpp:350
static int first_try_init_failures_
Definition: nrndaspk.h:38
void info()
Definition: nrndaspk.cpp:187
N_Vector yp_
Definition: nrndaspk.h:29
N_Vector delta_
Definition: nrndaspk.h:30
void * mem_
Definition: nrndaspk.h:27
Cvode * cv_
Definition: nrndaspk.h:28
Daspk(Cvode *, int neq)
Definition: nrndaspk.cpp:142
int advance_tn(double tstop)
Definition: nrndaspk.cpp:324
N_Vector acorvec()
Definition: nrndaspk.cpp:673
char * spmat_
Definition: nrndaspk.h:34
N_Vector ewtvec()
Definition: nrndaspk.cpp:669
static int init_failure_style_
Definition: nrndaspk.h:35
N_Vector parasite_
Definition: nrndaspk.h:31
static bool eq(double x, double y, double e)
Definition: mymath.h:66
static double eps(double x)
Definition: netcvode.h:139
double rtol_
Definition: netcvode.h:161
#define spPrint
Definition: cspredef.h:33
double t
Definition: cvodeobj.cpp:59
void hoc_execerror(const char *, const char *)
Definition: hoc.cpp:754
#define assert(ex)
Definition: hocassrt.h:32
#define i
Definition: md1redef.h:12
#define Printf
Definition: model.h:237
exp
Definition: extdef.h:5
void nrn_multithread_job(void *(*job)(NrnThread *))
Definition: multicore.cpp:1136
#define printf
Definition: mwprefix.h:26
void nrndae_dkpsol(double)
void nrndae_dkres(double *, double *, double *)
Definition: nrndae.cpp:92
void nrn_rhs(NrnThread *)
Definition: treeset.cpp:357
#define FACTORED
Definition: nrndaspk.cpp:23
static int res_gvardt(realtype t, N_Vector y, N_Vector yp, N_Vector delta, void *rdata)
Definition: nrndaspk.cpp:93
static N_Vector nvec_yp
Definition: nrndaspk.cpp:74
static double check(double t, Daspk *ida)
Definition: nrndaspk.cpp:209
static void daspk_nrn_solve(NrnThread *nt)
Definition: nrndaspk.cpp:48
#define thread_t
Definition: nrndaspk.cpp:70
void nrn_daspk_init_step(double, double, int)
Definition: fadvance.cpp:335
static N_Vector nvec_delta
Definition: nrndaspk.cpp:75
static void * msolve_thread(NrnThread *nt)
Definition: nrndaspk.cpp:118
static int solve_state_
Definition: nrndaspk.cpp:24
static double thread_cj
Definition: nrndaspk.cpp:77
void nrn_solve(NrnThread *)
Definition: solve.cpp:339
static int thread_ier
Definition: nrndaspk.cpp:78
static int msolve(IDAMem mem, N_Vector b, N_Vector ycur, N_Vector ypcur, N_Vector deltacur)
static N_Vector nvec_y
Definition: nrndaspk.cpp:73
void nrn_lhs(NrnThread *)
Definition: treeset.cpp:491
static int mfree(IDAMem)
Definition: nrndaspk.cpp:138
booleantype IDAEwtSet(IDAMem IDA_mem, N_Vector ycur)
#define SETUP
Definition: nrndaspk.cpp:22
static void * do_ode_thread(NrnThread *nt)
Definition: nrndaspk.cpp:196
#define INVALID
Definition: nrndaspk.cpp:20
#define nt_t
Definition: nrndaspk.cpp:46
static void * daspk_gather_thread(NrnThread *nt)
Definition: nrndaspk.cpp:417
static void * daspk_scatter_thread(NrnThread *nt)
Definition: nrndaspk.cpp:389
static void * res_thread(NrnThread *nt)
Definition: nrndaspk.cpp:80
static int msetup(IDAMem mem, N_Vector y, N_Vector ydot, N_Vector delta, N_Vector tempv1, N_Vector tempv2, N_Vector tempv3)
Definition: nrndaspk.cpp:111
static int minit(IDAMem)
Definition: nrndaspk.cpp:105
static Cvode * thread_cv
Definition: nrndaspk.cpp:79
int const size_t const size_t n
Definition: nrngsl.h:11
size_t p
size_t j
static philox4x32_key_t k
Definition: nrnran123.cpp:11
#define e
Definition: passive0.cpp:22
int use_sparse13
Definition: treeset.cpp:71
#define NODEV(n)
Definition: section.h:115
#define NODEAREA(n)
Definition: section.h:116
#define nlayer
Definition: section.h:188
#define NULL
Definition: sptree.h:16
double * _nrn_sav_rhs
Definition: multicore.h:46
double * param
Definition: section.h:190
double * v
Definition: section.h:196
int nodecount
Definition: nrnoc_ml.h:18
Node ** nodelist
Definition: nrnoc_ml.h:5
double ** data
Definition: nrnoc_ml.h:14
Definition: section.h:133
struct Extnode * extnode
Definition: section.h:161
int eqn_index_
Definition: section.h:149
int v_node_index
Definition: section.h:175
Represent main neuron object computed by single thread.
Definition: multicore.h:58
double _dt
Definition: multicore.h:60
_nrn_Fast_Imem * _nrn_fast_imem
Definition: multicore.h:82
int id
Definition: multicore.h:66
double cj
Definition: multicore.h:61
void * _vcv
Definition: multicore.h:83
double _t
Definition: multicore.h:59