1 /*
2  * -----------------------------------------------------------------
3  * $Revision$
4  * $Date$
5  * -----------------------------------------------------------------
6  * Programmer(s): Alan C. Hindmarsh and Radu Serban @ LLNL
7  * -----------------------------------------------------------------
8  * SUNDIALS Copyright Start
9  * Copyright (c) 2002-2020, Lawrence Livermore National Security
10  * and Southern Methodist University.
11  * All rights reserved.
12  *
13  * See the top-level LICENSE and NOTICE files for details.
14  *
15  * SPDX-License-Identifier: BSD-3-Clause
16  * SUNDIALS Copyright End
17  * -----------------------------------------------------------------
18  * This is the implementation file for the optional input and output
19  * functions for the CVODES solver.
20  * -----------------------------------------------------------------
21  */
22 
23 #include <stdio.h>
24 #include <stdlib.h>
25 
26 #include "cvodes_impl.h"
27 
28 #include <sundials/sundials_math.h>
29 #include <sundials/sundials_types.h>
30 
31 #define ZERO   RCONST(0.0)
32 #define HALF   RCONST(0.5)
33 #define ONE    RCONST(1.0)
34 #define TWOPT5 RCONST(2.5)
35 
36 /*
37  * =================================================================
38  * CVODES optional input functions
39  * =================================================================
40  */
41 
42 /*
43  * CVodeSetErrHandlerFn
44  *
45  * Specifies the error handler function
46  */
47 
CVodeSetErrHandlerFn(void * cvode_mem,CVErrHandlerFn ehfun,void * eh_data)48 int CVodeSetErrHandlerFn(void *cvode_mem, CVErrHandlerFn ehfun, void *eh_data)
49 {
50   CVodeMem cv_mem;
51 
52   if (cvode_mem==NULL) {
53     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetErrHandlerFn", MSGCV_NO_MEM);
54     return(CV_MEM_NULL);
55   }
56 
57   cv_mem = (CVodeMem) cvode_mem;
58 
59   cv_mem->cv_ehfun = ehfun;
60   cv_mem->cv_eh_data = eh_data;
61 
62   return(CV_SUCCESS);
63 }
64 
65 /*
66  * CVodeSetErrFile
67  *
68  * Specifies the FILE pointer for output (NULL means no messages)
69  */
70 
CVodeSetErrFile(void * cvode_mem,FILE * errfp)71 int CVodeSetErrFile(void *cvode_mem, FILE *errfp)
72 {
73   CVodeMem cv_mem;
74 
75   if (cvode_mem==NULL) {
76     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetErrFile", MSGCV_NO_MEM);
77     return(CV_MEM_NULL);
78   }
79 
80   cv_mem = (CVodeMem) cvode_mem;
81 
82   cv_mem->cv_errfp = errfp;
83 
84   return(CV_SUCCESS);
85 }
86 
87 /*
88  * CVodeSetUserData
89  *
90  * Specifies the user data pointer for f
91  */
92 
CVodeSetUserData(void * cvode_mem,void * user_data)93 int CVodeSetUserData(void *cvode_mem, void *user_data)
94 {
95   CVodeMem cv_mem;
96 
97   if (cvode_mem==NULL) {
98     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetUserData", MSGCV_NO_MEM);
99     return(CV_MEM_NULL);
100   }
101 
102   cv_mem = (CVodeMem) cvode_mem;
103 
104   cv_mem->cv_user_data = user_data;
105 
106   return(CV_SUCCESS);
107 }
108 
109 /*
110  * CVodeSetMaxOrd
111  *
112  * Specifies the maximum method order
113  */
114 
CVodeSetMaxOrd(void * cvode_mem,int maxord)115 int CVodeSetMaxOrd(void *cvode_mem, int maxord)
116 {
117   CVodeMem cv_mem;
118   int qmax_alloc;
119 
120   if (cvode_mem==NULL) {
121     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetMaxOrd", MSGCV_NO_MEM);
122     return(CV_MEM_NULL);
123   }
124 
125   cv_mem = (CVodeMem) cvode_mem;
126 
127   if (maxord <= 0) {
128     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetMaxOrd", MSGCV_NEG_MAXORD);
129     return(CV_ILL_INPUT);
130   }
131 
132   /* Cannot increase maximum order beyond the value that
133      was used when allocating memory */
134   qmax_alloc = cv_mem->cv_qmax_alloc;
135   qmax_alloc = SUNMIN(qmax_alloc, cv_mem->cv_qmax_allocQ);
136   qmax_alloc = SUNMIN(qmax_alloc, cv_mem->cv_qmax_allocS);
137 
138   if (maxord > qmax_alloc) {
139     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetMaxOrd", MSGCV_BAD_MAXORD);
140     return(CV_ILL_INPUT);
141   }
142 
143   cv_mem->cv_qmax = maxord;
144 
145   return(CV_SUCCESS);
146 }
147 
148 /*
149  * CVodeSetMaxNumSteps
150  *
151  * Specifies the maximum number of integration steps
152  */
153 
CVodeSetMaxNumSteps(void * cvode_mem,long int mxsteps)154 int CVodeSetMaxNumSteps(void *cvode_mem, long int mxsteps)
155 {
156   CVodeMem cv_mem;
157 
158   if (cvode_mem==NULL) {
159     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetMaxNumSteps", MSGCV_NO_MEM);
160     return(CV_MEM_NULL);
161   }
162 
163   cv_mem = (CVodeMem) cvode_mem;
164 
165   /* Passing mxsteps=0 sets the default. Passing mxsteps<0 disables the test. */
166   if (mxsteps == 0)
167     cv_mem->cv_mxstep = MXSTEP_DEFAULT;
168   else
169     cv_mem->cv_mxstep = mxsteps;
170 
171   return(CV_SUCCESS);
172 }
173 
174 /*
175  * CVodeSetMaxHnilWarns
176  *
177  * Specifies the maximum number of warnings for small h
178  */
179 
CVodeSetMaxHnilWarns(void * cvode_mem,int mxhnil)180 int CVodeSetMaxHnilWarns(void *cvode_mem, int mxhnil)
181 {
182   CVodeMem cv_mem;
183 
184   if (cvode_mem==NULL) {
185     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetMaxHnilWarns", MSGCV_NO_MEM);
186     return(CV_MEM_NULL);
187   }
188 
189   cv_mem = (CVodeMem) cvode_mem;
190 
191   cv_mem->cv_mxhnil = mxhnil;
192 
193   return(CV_SUCCESS);
194 }
195 
196 /*
197  *CVodeSetStabLimDet
198  *
199  * Turns on/off the stability limit detection algorithm
200  */
201 
CVodeSetStabLimDet(void * cvode_mem,booleantype sldet)202 int CVodeSetStabLimDet(void *cvode_mem, booleantype sldet)
203 {
204   CVodeMem cv_mem;
205 
206   if (cvode_mem==NULL) {
207     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetStabLimDet", MSGCV_NO_MEM);
208     return(CV_MEM_NULL);
209   }
210 
211   cv_mem = (CVodeMem) cvode_mem;
212 
213   if( sldet && (cv_mem->cv_lmm != CV_BDF) ) {
214     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetStabLimDet", MSGCV_SET_SLDET);
215     return(CV_ILL_INPUT);
216   }
217 
218   cv_mem->cv_sldeton = sldet;
219 
220   return(CV_SUCCESS);
221 }
222 
223 /*
224  * CVodeSetInitStep
225  *
226  * Specifies the initial step size
227  */
228 
CVodeSetInitStep(void * cvode_mem,realtype hin)229 int CVodeSetInitStep(void *cvode_mem, realtype hin)
230 {
231   CVodeMem cv_mem;
232 
233   if (cvode_mem==NULL) {
234     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetInitStep", MSGCV_NO_MEM);
235     return(CV_MEM_NULL);
236   }
237 
238   cv_mem = (CVodeMem) cvode_mem;
239 
240   cv_mem->cv_hin = hin;
241 
242   return(CV_SUCCESS);
243 }
244 
245 /*
246  * CVodeSetMinStep
247  *
248  * Specifies the minimum step size
249  */
250 
CVodeSetMinStep(void * cvode_mem,realtype hmin)251 int CVodeSetMinStep(void *cvode_mem, realtype hmin)
252 {
253   CVodeMem cv_mem;
254 
255   if (cvode_mem==NULL) {
256     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetMinStep", MSGCV_NO_MEM);
257     return(CV_MEM_NULL);
258   }
259 
260   cv_mem = (CVodeMem) cvode_mem;
261 
262   if (hmin<0) {
263     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetMinStep", MSGCV_NEG_HMIN);
264     return(CV_ILL_INPUT);
265   }
266 
267   /* Passing 0 sets hmin = zero */
268   if (hmin == ZERO) {
269     cv_mem->cv_hmin = HMIN_DEFAULT;
270     return(CV_SUCCESS);
271   }
272 
273   if (hmin * cv_mem->cv_hmax_inv > ONE) {
274     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetMinStep", MSGCV_BAD_HMIN_HMAX);
275     return(CV_ILL_INPUT);
276   }
277 
278   cv_mem->cv_hmin = hmin;
279 
280   return(CV_SUCCESS);
281 }
282 
283 /*
284  * CVodeSetMaxStep
285  *
286  * Specifies the maximum step size
287  */
288 
CVodeSetMaxStep(void * cvode_mem,realtype hmax)289 int CVodeSetMaxStep(void *cvode_mem, realtype hmax)
290 {
291   realtype hmax_inv;
292   CVodeMem cv_mem;
293 
294   if (cvode_mem==NULL) {
295     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetMaxStep", MSGCV_NO_MEM);
296     return (CV_MEM_NULL);
297   }
298 
299   cv_mem = (CVodeMem) cvode_mem;
300 
301   if (hmax < 0) {
302     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetMaxStep", MSGCV_NEG_HMAX);
303     return(CV_ILL_INPUT);
304   }
305 
306   /* Passing 0 sets hmax = infinity */
307   if (hmax == ZERO) {
308     cv_mem->cv_hmax_inv = HMAX_INV_DEFAULT;
309     return(CV_SUCCESS);
310   }
311 
312   hmax_inv = ONE/hmax;
313   if (hmax_inv * cv_mem->cv_hmin > ONE) {
314     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetMaxStep", MSGCV_BAD_HMIN_HMAX);
315     return(CV_ILL_INPUT);
316   }
317 
318   cv_mem->cv_hmax_inv = hmax_inv;
319 
320   return(CV_SUCCESS);
321 }
322 
323 /*
324  * CVodeSetStopTime
325  *
326  * Specifies the time beyond which the integration is not to proceed.
327  */
328 
CVodeSetStopTime(void * cvode_mem,realtype tstop)329 int CVodeSetStopTime(void *cvode_mem, realtype tstop)
330 {
331   CVodeMem cv_mem;
332 
333   if (cvode_mem==NULL) {
334     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetStopTime", MSGCV_NO_MEM);
335     return (CV_MEM_NULL);
336   }
337   cv_mem = (CVodeMem) cvode_mem;
338 
339   /* If CVode was called at least once, test if tstop is legal
340    * (i.e. if it was not already passed).
341    * If CVodeSetStopTime is called before the first call to CVode,
342    * tstop will be checked in CVode. */
343   if (cv_mem->cv_nst > 0) {
344 
345     if ( (tstop - cv_mem->cv_tn) * cv_mem->cv_h < ZERO ) {
346       cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetStopTime", MSGCV_BAD_TSTOP, tstop, cv_mem->cv_tn);
347       return(CV_ILL_INPUT);
348     }
349 
350   }
351 
352   cv_mem->cv_tstop = tstop;
353   cv_mem->cv_tstopset = SUNTRUE;
354 
355   return(CV_SUCCESS);
356 }
357 
358 /*
359  * CVodeSetMaxErrTestFails
360  *
361  * Specifies the maximum number of error test failures during one
362  * step try.
363  */
364 
CVodeSetMaxErrTestFails(void * cvode_mem,int maxnef)365 int CVodeSetMaxErrTestFails(void *cvode_mem, int maxnef)
366 {
367   CVodeMem cv_mem;
368 
369   if (cvode_mem==NULL) {
370     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetMaxErrTestFails", MSGCV_NO_MEM);
371     return (CV_MEM_NULL);
372   }
373 
374   cv_mem = (CVodeMem) cvode_mem;
375 
376   cv_mem->cv_maxnef = maxnef;
377 
378   return(CV_SUCCESS);
379 }
380 
381 /*
382  * CVodeSetMaxConvFails
383  *
384  * Specifies the maximum number of nonlinear convergence failures
385  * during one step try.
386  */
387 
CVodeSetMaxConvFails(void * cvode_mem,int maxncf)388 int CVodeSetMaxConvFails(void *cvode_mem, int maxncf)
389 {
390   CVodeMem cv_mem;
391 
392   if (cvode_mem==NULL) {
393     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetMaxConvFails", MSGCV_NO_MEM);
394     return (CV_MEM_NULL);
395   }
396 
397   cv_mem = (CVodeMem) cvode_mem;
398 
399   cv_mem->cv_maxncf = maxncf;
400 
401   return(CV_SUCCESS);
402 }
403 
404 /*
405  * CVodeSetMaxNonlinIters
406  *
407  * Specifies the maximum number of nonlinear iterations during
408  * one solve.
409  */
410 
CVodeSetMaxNonlinIters(void * cvode_mem,int maxcor)411 int CVodeSetMaxNonlinIters(void *cvode_mem, int maxcor)
412 {
413   CVodeMem cv_mem;
414   booleantype sensi_sim;
415 
416   if (cvode_mem==NULL) {
417     cvProcessError(NULL, CV_MEM_NULL, "CVODES",
418                    "CVodeSetMaxNonlinIters", MSGCV_NO_MEM);
419     return (CV_MEM_NULL);
420   }
421 
422   cv_mem = (CVodeMem) cvode_mem;
423 
424   /* Are we computing sensitivities with the simultaneous approach? */
425   sensi_sim = (cv_mem->cv_sensi && (cv_mem->cv_ism==CV_SIMULTANEOUS));
426 
427   if (sensi_sim) {
428 
429     /* check that the NLS is non-NULL */
430     if (cv_mem->NLSsim == NULL) {
431       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
432                       "CVodeSetMaxNonlinIters", MSGCV_MEM_FAIL);
433       return(CV_MEM_FAIL);
434     }
435 
436     return(SUNNonlinSolSetMaxIters(cv_mem->NLSsim, maxcor));
437 
438   } else {
439 
440     /* check that the NLS is non-NULL */
441     if (cv_mem->NLS == NULL) {
442       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
443                       "CVodeSetMaxNonlinIters", MSGCV_MEM_FAIL);
444       return(CV_MEM_FAIL);
445     }
446 
447     return(SUNNonlinSolSetMaxIters(cv_mem->NLS, maxcor));
448   }
449 
450   return(CV_SUCCESS);
451 }
452 
453 /*
454  * CVodeSetNonlinConvCoef
455  *
456  * Specifies the coeficient in the nonlinear solver convergence
457  * test
458  */
459 
CVodeSetNonlinConvCoef(void * cvode_mem,realtype nlscoef)460 int CVodeSetNonlinConvCoef(void *cvode_mem, realtype nlscoef)
461 {
462   CVodeMem cv_mem;
463 
464   if (cvode_mem==NULL) {
465     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetNonlinConvCoef", MSGCV_NO_MEM);
466     return(CV_MEM_NULL);
467   }
468 
469   cv_mem = (CVodeMem) cvode_mem;
470 
471   cv_mem->cv_nlscoef = nlscoef;
472 
473   return(CV_SUCCESS);
474 }
475 
476 /*
477  * CVodeSetRootDirection
478  *
479  * Specifies the direction of zero-crossings to be monitored.
480  * The default is to monitor both crossings.
481  */
482 
CVodeSetRootDirection(void * cvode_mem,int * rootdir)483 int CVodeSetRootDirection(void *cvode_mem, int *rootdir)
484 {
485   CVodeMem cv_mem;
486   int i, nrt;
487 
488   if (cvode_mem==NULL) {
489     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetRootDirection", MSGCV_NO_MEM);
490     return(CV_MEM_NULL);
491   }
492 
493   cv_mem = (CVodeMem) cvode_mem;
494 
495   nrt = cv_mem->cv_nrtfn;
496   if (nrt==0) {
497     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetRootDirection", MSGCV_NO_ROOT);
498     return(CV_ILL_INPUT);
499   }
500 
501   for(i=0; i<nrt; i++) cv_mem->cv_rootdir[i] = rootdir[i];
502 
503   return(CV_SUCCESS);
504 }
505 
506 
507 /*
508  * CVodeSetNoInactiveRootWarn
509  *
510  * Disables issuing a warning if some root function appears
511  * to be identically zero at the beginning of the integration
512  */
513 
CVodeSetNoInactiveRootWarn(void * cvode_mem)514 int CVodeSetNoInactiveRootWarn(void *cvode_mem)
515 {
516   CVodeMem cv_mem;
517 
518   if (cvode_mem==NULL) {
519     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetNoInactiveRootWarn", MSGCV_NO_MEM);
520     return(CV_MEM_NULL);
521   }
522 
523   cv_mem = (CVodeMem) cvode_mem;
524 
525   cv_mem->cv_mxgnull = 0;
526 
527   return(CV_SUCCESS);
528 }
529 
530 /*
531  * CVodeSetConstraints
532  *
533  * Setup for constraint handling feature
534  */
535 
CVodeSetConstraints(void * cvode_mem,N_Vector constraints)536 int CVodeSetConstraints(void *cvode_mem, N_Vector constraints)
537 {
538   CVodeMem cv_mem;
539   realtype temptest;
540 
541   if (cvode_mem==NULL) {
542     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetConstraints", MSGCV_NO_MEM);
543     return(CV_MEM_NULL);
544   }
545 
546   cv_mem = (CVodeMem) cvode_mem;
547 
548   /* If there are no constraints, destroy data structures */
549   if (constraints == NULL) {
550     if (cv_mem->cv_constraintsMallocDone) {
551       N_VDestroy(cv_mem->cv_constraints);
552       cv_mem->cv_lrw -= cv_mem->cv_lrw1;
553       cv_mem->cv_liw -= cv_mem->cv_liw1;
554     }
555     cv_mem->cv_constraintsMallocDone = SUNFALSE;
556     cv_mem->cv_constraintsSet = SUNFALSE;
557     return(CV_SUCCESS);
558   }
559 
560   /* Test if required vector ops. are defined */
561 
562   if (constraints->ops->nvdiv         == NULL ||
563       constraints->ops->nvmaxnorm     == NULL ||
564       constraints->ops->nvcompare     == NULL ||
565       constraints->ops->nvconstrmask  == NULL ||
566       constraints->ops->nvminquotient == NULL) {
567     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetConstraints", MSGCV_BAD_NVECTOR);
568     return(CV_ILL_INPUT);
569   }
570 
571   /* Check the constraints vector */
572   temptest = N_VMaxNorm(constraints);
573   if ((temptest > TWOPT5) || (temptest < HALF)) {
574     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetConstraints", MSGCV_BAD_CONSTR);
575     return(CV_ILL_INPUT);
576   }
577 
578   if ( !(cv_mem->cv_constraintsMallocDone) ) {
579     cv_mem->cv_constraints = N_VClone(constraints);
580     cv_mem->cv_lrw += cv_mem->cv_lrw1;
581     cv_mem->cv_liw += cv_mem->cv_liw1;
582     cv_mem->cv_constraintsMallocDone = SUNTRUE;
583   }
584 
585   /* Load the constraints vector */
586   N_VScale(ONE, constraints, cv_mem->cv_constraints);
587 
588   cv_mem->cv_constraintsSet = SUNTRUE;
589 
590   return(CV_SUCCESS);
591 }
592 
593 /*
594  * =================================================================
595  * Quadrature optional input functions
596  * =================================================================
597  */
598 
CVodeSetQuadErrCon(void * cvode_mem,booleantype errconQ)599 int CVodeSetQuadErrCon(void *cvode_mem, booleantype errconQ)
600 {
601   CVodeMem cv_mem;
602 
603   if (cvode_mem==NULL) {
604     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetQuadErrCon", MSGCV_NO_MEM);
605     return(CV_MEM_NULL);
606   }
607   cv_mem = (CVodeMem) cvode_mem;
608 
609   cv_mem->cv_errconQ = errconQ;
610 
611   return(CV_SUCCESS);
612 }
613 
614 /*
615  * =================================================================
616  * FSA optional input functions
617  * =================================================================
618  */
619 
CVodeSetSensDQMethod(void * cvode_mem,int DQtype,realtype DQrhomax)620 int CVodeSetSensDQMethod(void *cvode_mem, int DQtype, realtype DQrhomax)
621 {
622   CVodeMem cv_mem;
623 
624   if (cvode_mem==NULL) {
625     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetSensDQMethod", MSGCV_NO_MEM);
626     return(CV_MEM_NULL);
627   }
628 
629   cv_mem = (CVodeMem) cvode_mem;
630 
631   if ( (DQtype != CV_CENTERED) && (DQtype != CV_FORWARD) ) {
632     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetSensDQMethod", MSGCV_BAD_DQTYPE);
633     return(CV_ILL_INPUT);
634   }
635 
636   if (DQrhomax < ZERO ) {
637     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetSensDQMethod", MSGCV_BAD_DQRHO);
638     return(CV_ILL_INPUT);
639   }
640 
641   cv_mem->cv_DQtype = DQtype;
642   cv_mem->cv_DQrhomax = DQrhomax;
643 
644   return(CV_SUCCESS);
645 }
646 
647 /*-----------------------------------------------------------------*/
648 
CVodeSetSensErrCon(void * cvode_mem,booleantype errconS)649 int CVodeSetSensErrCon(void *cvode_mem, booleantype errconS)
650 {
651   CVodeMem cv_mem;
652 
653   if (cvode_mem==NULL) {
654     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetSensErrCon", MSGCV_NO_MEM);
655     return(CV_MEM_NULL);
656   }
657   cv_mem = (CVodeMem) cvode_mem;
658 
659   cv_mem->cv_errconS = errconS;
660 
661   return(CV_SUCCESS);
662 }
663 
664 /*-----------------------------------------------------------------*/
665 
CVodeSetSensMaxNonlinIters(void * cvode_mem,int maxcorS)666 int CVodeSetSensMaxNonlinIters(void *cvode_mem, int maxcorS)
667 {
668   CVodeMem cv_mem;
669   booleantype sensi_stg;
670 
671   if (cvode_mem==NULL) {
672     cvProcessError(NULL, CV_MEM_NULL, "CVODES",
673                    "CVodeSetSensMaxNonlinIters", MSGCV_NO_MEM);
674     return (CV_MEM_NULL);
675   }
676 
677   cv_mem = (CVodeMem) cvode_mem;
678 
679   /* Are we computing sensitivities with a staggered approach? */
680   sensi_stg = (cv_mem->cv_sensi && (cv_mem->cv_ism==CV_STAGGERED));
681 
682   if (sensi_stg) {
683 
684     /* check that the NLS is non-NULL */
685     if (cv_mem->NLSstg == NULL) {
686       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
687                       "CVodeSetSensMaxNonlinIters", MSGCV_MEM_FAIL);
688       return(CV_MEM_FAIL);
689     }
690 
691     return(SUNNonlinSolSetMaxIters(cv_mem->NLSstg, maxcorS));
692 
693   } else {
694 
695     /* check that the NLS is non-NULL */
696     if (cv_mem->NLSstg1 == NULL) {
697       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
698                       "CVodeSetMaxNonlinIters", MSGCV_MEM_FAIL);
699       return(CV_MEM_FAIL);
700     }
701 
702     return(SUNNonlinSolSetMaxIters(cv_mem->NLSstg1, maxcorS));
703   }
704 
705   return(CV_SUCCESS);
706 }
707 
708 /*-----------------------------------------------------------------*/
709 
CVodeSetSensParams(void * cvode_mem,realtype * p,realtype * pbar,int * plist)710 int CVodeSetSensParams(void *cvode_mem, realtype *p, realtype *pbar, int *plist)
711 {
712   CVodeMem cv_mem;
713   int is, Ns;
714 
715   if (cvode_mem==NULL) {
716     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetSensParams", MSGCV_NO_MEM);
717     return(CV_MEM_NULL);
718   }
719 
720   cv_mem = (CVodeMem) cvode_mem;
721 
722   /* Was sensitivity initialized? */
723 
724   if (cv_mem->cv_SensMallocDone == SUNFALSE) {
725     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeSetSensParams", MSGCV_NO_SENSI);
726     return(CV_NO_SENS);
727   }
728 
729   Ns = cv_mem->cv_Ns;
730 
731   /* Parameters */
732 
733   cv_mem->cv_p = p;
734 
735   /* pbar */
736 
737   if (pbar != NULL)
738     for (is=0; is<Ns; is++) {
739       if (pbar[is] == ZERO) {
740         cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetSensParams", MSGCV_BAD_PBAR);
741         return(CV_ILL_INPUT);
742       }
743       cv_mem->cv_pbar[is] = SUNRabs(pbar[is]);
744     }
745   else
746     for (is=0; is<Ns; is++)
747       cv_mem->cv_pbar[is] = ONE;
748 
749   /* plist */
750 
751   if (plist != NULL)
752     for (is=0; is<Ns; is++) {
753       if ( plist[is] < 0 ) {
754         cvProcessError(cv_mem, CV_ILL_INPUT, "CVODES", "CVodeSetSensParams", MSGCV_BAD_PLIST);
755         return(CV_ILL_INPUT);
756       }
757       cv_mem->cv_plist[is] = plist[is];
758     }
759   else
760     for (is=0; is<Ns; is++)
761       cv_mem->cv_plist[is] = is;
762 
763   return(CV_SUCCESS);
764 }
765 
766 /*-----------------------------------------------------------------*/
767 
CVodeSetQuadSensErrCon(void * cvode_mem,booleantype errconQS)768 int CVodeSetQuadSensErrCon(void *cvode_mem, booleantype errconQS)
769 {
770   CVodeMem cv_mem;
771 
772   if (cvode_mem==NULL) {
773     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeSetQuadSensErrCon", MSGCV_NO_MEM);
774     return(CV_MEM_NULL);
775   }
776   cv_mem = (CVodeMem) cvode_mem;
777 
778   /* Was sensitivity initialized? */
779 
780   if (cv_mem->cv_SensMallocDone == SUNFALSE) {
781     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeSetQuadSensTolerances", MSGCV_NO_SENSI);
782     return(CV_NO_SENS);
783   }
784 
785   /* Ckeck if quadrature sensitivity was initialized? */
786 
787   if (cv_mem->cv_QuadSensMallocDone == SUNFALSE) {
788     cvProcessError(cv_mem, CV_NO_QUADSENS, "CVODES", "CVodeSetQuadSensErrCon", MSGCV_NO_QUADSENSI);
789     return(CV_NO_QUAD);
790   }
791 
792   cv_mem->cv_errconQS = errconQS;
793 
794   return(CV_SUCCESS);
795 }
796 
797 /*
798  * =================================================================
799  * CVODES optional output functions
800  * =================================================================
801  */
802 
803 /*
804  * CVodeGetNumSteps
805  *
806  * Returns the current number of integration steps
807  */
808 
CVodeGetNumSteps(void * cvode_mem,long int * nsteps)809 int CVodeGetNumSteps(void *cvode_mem, long int *nsteps)
810 {
811   CVodeMem cv_mem;
812 
813   if (cvode_mem==NULL) {
814     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetNumSteps", MSGCV_NO_MEM);
815     return(CV_MEM_NULL);
816   }
817 
818   cv_mem = (CVodeMem) cvode_mem;
819 
820   *nsteps = cv_mem->cv_nst;
821 
822   return(CV_SUCCESS);
823 }
824 
825 /*
826  * CVodeGetNumRhsEvals
827  *
828  * Returns the current number of calls to f
829  */
830 
CVodeGetNumRhsEvals(void * cvode_mem,long int * nfevals)831 int CVodeGetNumRhsEvals(void *cvode_mem, long int *nfevals)
832 {
833   CVodeMem cv_mem;
834 
835   if (cvode_mem==NULL) {
836     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetNumRhsEvals", MSGCV_NO_MEM);
837     return(CV_MEM_NULL);
838   }
839 
840   cv_mem = (CVodeMem) cvode_mem;
841 
842   *nfevals = cv_mem->cv_nfe;
843 
844   return(CV_SUCCESS);
845 }
846 
847 /*
848  * CVodeGetNumLinSolvSetups
849  *
850  * Returns the current number of calls to the linear solver setup routine
851  */
852 
CVodeGetNumLinSolvSetups(void * cvode_mem,long int * nlinsetups)853 int CVodeGetNumLinSolvSetups(void *cvode_mem, long int *nlinsetups)
854 {
855   CVodeMem cv_mem;
856 
857   if (cvode_mem==NULL) {
858     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetNumLinSolvSetups", MSGCV_NO_MEM);
859     return(CV_MEM_NULL);
860   }
861 
862   cv_mem = (CVodeMem) cvode_mem;
863 
864   *nlinsetups = cv_mem->cv_nsetups;
865 
866   return(CV_SUCCESS);
867 }
868 
869 /*
870  * CVodeGetNumErrTestFails
871  *
872  * Returns the current number of error test failures
873  */
874 
CVodeGetNumErrTestFails(void * cvode_mem,long int * netfails)875 int CVodeGetNumErrTestFails(void *cvode_mem, long int *netfails)
876 {
877   CVodeMem cv_mem;
878 
879   if (cvode_mem==NULL) {
880     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetNumErrTestFails", MSGCV_NO_MEM);
881     return(CV_MEM_NULL);
882   }
883 
884   cv_mem = (CVodeMem) cvode_mem;
885 
886   *netfails = cv_mem->cv_netf;
887 
888   return(CV_SUCCESS);
889 }
890 
891 /*
892  * CVodeGetLastOrder
893  *
894  * Returns the order on the last succesful step
895  */
896 
CVodeGetLastOrder(void * cvode_mem,int * qlast)897 int CVodeGetLastOrder(void *cvode_mem, int *qlast)
898 {
899   CVodeMem cv_mem;
900 
901   if (cvode_mem==NULL) {
902     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetLastOrder", MSGCV_NO_MEM);
903     return(CV_MEM_NULL);
904   }
905 
906   cv_mem = (CVodeMem) cvode_mem;
907 
908   *qlast = cv_mem->cv_qu;
909 
910   return(CV_SUCCESS);
911 }
912 
913 /*
914  * CVodeGetCurrentOrder
915  *
916  * Returns the order to be attempted on the next step
917  */
918 
CVodeGetCurrentOrder(void * cvode_mem,int * qcur)919 int CVodeGetCurrentOrder(void *cvode_mem, int *qcur)
920 {
921   CVodeMem cv_mem;
922 
923   if (cvode_mem==NULL) {
924     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetCurrentOrder", MSGCV_NO_MEM);
925     return(CV_MEM_NULL);
926   }
927 
928   cv_mem = (CVodeMem) cvode_mem;
929 
930   *qcur = cv_mem->cv_next_q;
931 
932   return(CV_SUCCESS);
933 }
934 
935 /*
936  * CVodeGetCurrentGamma
937  *
938  * Returns the value of gamma for the current step.
939  */
940 
CVodeGetCurrentGamma(void * cvode_mem,realtype * gamma)941 int CVodeGetCurrentGamma(void *cvode_mem, realtype *gamma)
942 {
943   CVodeMem cv_mem;
944 
945   if (cvode_mem==NULL) {
946     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetCurrentGamma", MSGCV_NO_MEM);
947     return(CV_MEM_NULL);
948   }
949 
950   cv_mem = (CVodeMem) cvode_mem;
951 
952   *gamma = cv_mem->cv_gamma;
953 
954   return(CV_SUCCESS);
955 }
956 
957 /*
958  * CVodeGetNumStabLimOrderReds
959  *
960  * Returns the number of order reductions triggered by the stability
961  * limit detection algorithm
962  */
963 
CVodeGetNumStabLimOrderReds(void * cvode_mem,long int * nslred)964 int CVodeGetNumStabLimOrderReds(void *cvode_mem, long int *nslred)
965 {
966   CVodeMem cv_mem;
967 
968   if (cvode_mem==NULL) {
969     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetNumStabLimOrderReds", MSGCV_NO_MEM);
970     return(CV_MEM_NULL);
971   }
972 
973   cv_mem = (CVodeMem) cvode_mem;
974 
975   if (cv_mem->cv_sldeton==SUNFALSE)
976     *nslred = 0;
977   else
978     *nslred = cv_mem->cv_nor;
979 
980   return(CV_SUCCESS);
981 }
982 
983 /*
984  * CVodeGetActualInitStep
985  *
986  * Returns the step size used on the first step
987  */
988 
CVodeGetActualInitStep(void * cvode_mem,realtype * hinused)989 int CVodeGetActualInitStep(void *cvode_mem, realtype *hinused)
990 {
991   CVodeMem cv_mem;
992 
993   if (cvode_mem==NULL) {
994     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetActualInitStep", MSGCV_NO_MEM);
995     return(CV_MEM_NULL);
996   }
997 
998   cv_mem = (CVodeMem) cvode_mem;
999 
1000   *hinused = cv_mem->cv_h0u;
1001 
1002   return(CV_SUCCESS);
1003 }
1004 
1005 /*
1006  * CVodeGetLastStep
1007  *
1008  * Returns the step size used on the last successful step
1009  */
1010 
CVodeGetLastStep(void * cvode_mem,realtype * hlast)1011 int CVodeGetLastStep(void *cvode_mem, realtype *hlast)
1012 {
1013   CVodeMem cv_mem;
1014 
1015   if (cvode_mem==NULL) {
1016     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetLastStep", MSGCV_NO_MEM);
1017     return(CV_MEM_NULL);
1018   }
1019 
1020   cv_mem = (CVodeMem) cvode_mem;
1021 
1022   *hlast = cv_mem->cv_hu;
1023 
1024   return(CV_SUCCESS);
1025 }
1026 
1027 /*
1028  * CVodeGetCurrentStep
1029  *
1030  * Returns the step size to be attempted on the next step
1031  */
1032 
CVodeGetCurrentStep(void * cvode_mem,realtype * hcur)1033 int CVodeGetCurrentStep(void *cvode_mem, realtype *hcur)
1034 {
1035   CVodeMem cv_mem;
1036 
1037   if (cvode_mem==NULL) {
1038     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetCurrentStep", MSGCV_NO_MEM);
1039     return(CV_MEM_NULL);
1040   }
1041 
1042   cv_mem = (CVodeMem) cvode_mem;
1043 
1044   *hcur = cv_mem->cv_next_h;
1045 
1046   return(CV_SUCCESS);
1047 }
1048 
1049 /*
1050  * CVodeGetCurrentState
1051  *
1052  * Returns the current state vector
1053  */
1054 
CVodeGetCurrentState(void * cvode_mem,N_Vector * y)1055 int CVodeGetCurrentState(void *cvode_mem, N_Vector *y)
1056 {
1057   CVodeMem cv_mem;
1058 
1059   if (cvode_mem==NULL) {
1060     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetCurrentState", MSGCV_NO_MEM);
1061     return(CV_MEM_NULL);
1062   }
1063 
1064   cv_mem = (CVodeMem) cvode_mem;
1065 
1066   *y = cv_mem->cv_y;
1067 
1068   return(CV_SUCCESS);
1069 }
1070 
1071 /*
1072  * CVodeGetCurrentStateSens
1073  *
1074  * Returns the current sensitivity state vector array
1075  */
1076 
CVodeGetCurrentStateSens(void * cvode_mem,N_Vector ** yS)1077 int CVodeGetCurrentStateSens(void *cvode_mem, N_Vector **yS)
1078 {
1079   CVodeMem cv_mem;
1080 
1081   if (cvode_mem==NULL) {
1082     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetCurrentStateSens", MSGCV_NO_MEM);
1083     return(CV_MEM_NULL);
1084   }
1085 
1086   cv_mem = (CVodeMem) cvode_mem;
1087 
1088   *yS = cv_mem->cv_yS;
1089 
1090   return(CV_SUCCESS);
1091 }
1092 
1093 /*
1094  * CVodeGetCurrentSensSolveIndex
1095  *
1096  * Returns the current index of the sensitivity solve when using
1097  * the staggered1 nonlinear solver.
1098  */
1099 
CVodeGetCurrentSensSolveIndex(void * cvode_mem,int * index)1100 int CVodeGetCurrentSensSolveIndex(void *cvode_mem, int *index)
1101 {
1102   CVodeMem cv_mem;
1103 
1104   if (cvode_mem==NULL) {
1105     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetCurrentSensSolveIndex", MSGCV_NO_MEM);
1106     return(CV_MEM_NULL);
1107   }
1108 
1109   cv_mem = (CVodeMem) cvode_mem;
1110 
1111   *index = cv_mem->sens_solve_idx;
1112 
1113   return(CV_SUCCESS);
1114 }
1115 
1116 /*
1117  * CVodeGetCurrentTime
1118  *
1119  * Returns the current value of the independent variable
1120  */
1121 
CVodeGetCurrentTime(void * cvode_mem,realtype * tcur)1122 int CVodeGetCurrentTime(void *cvode_mem, realtype *tcur)
1123 {
1124   CVodeMem cv_mem;
1125 
1126   if (cvode_mem==NULL) {
1127     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetCurrentTime", MSGCV_NO_MEM);
1128     return(CV_MEM_NULL);
1129   }
1130 
1131   cv_mem = (CVodeMem) cvode_mem;
1132 
1133   *tcur = cv_mem->cv_tn;
1134 
1135   return(CV_SUCCESS);
1136 }
1137 
1138 /*
1139  * CVodeGetTolScaleFactor
1140  *
1141  * Returns a suggested factor for scaling tolerances
1142  */
1143 
CVodeGetTolScaleFactor(void * cvode_mem,realtype * tolsfact)1144 int CVodeGetTolScaleFactor(void *cvode_mem, realtype *tolsfact)
1145 {
1146   CVodeMem cv_mem;
1147 
1148   if (cvode_mem==NULL) {
1149     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetTolScaleFactor", MSGCV_NO_MEM);
1150     return(CV_MEM_NULL);
1151   }
1152 
1153   cv_mem = (CVodeMem) cvode_mem;
1154 
1155   *tolsfact = cv_mem->cv_tolsf;
1156 
1157   return(CV_SUCCESS);
1158 }
1159 
1160 /*
1161  * CVodeGetErrWeights
1162  *
1163  * This routine returns the current weight vector.
1164  */
1165 
CVodeGetErrWeights(void * cvode_mem,N_Vector eweight)1166 int CVodeGetErrWeights(void *cvode_mem, N_Vector eweight)
1167 {
1168   CVodeMem cv_mem;
1169 
1170   if (cvode_mem==NULL) {
1171     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetErrWeights", MSGCV_NO_MEM);
1172     return(CV_MEM_NULL);
1173   }
1174 
1175   cv_mem = (CVodeMem) cvode_mem;
1176 
1177   N_VScale(ONE, cv_mem->cv_ewt, eweight);
1178 
1179   return(CV_SUCCESS);
1180 }
1181 
1182 /*
1183  * CVodeGetEstLocalErrors
1184  *
1185  * Returns an estimate of the local error
1186  */
1187 
CVodeGetEstLocalErrors(void * cvode_mem,N_Vector ele)1188 int CVodeGetEstLocalErrors(void *cvode_mem, N_Vector ele)
1189 {
1190   CVodeMem cv_mem;
1191 
1192   if (cvode_mem==NULL) {
1193     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetEstLocalErrors", MSGCV_NO_MEM);
1194     return(CV_MEM_NULL);
1195   }
1196 
1197   cv_mem = (CVodeMem) cvode_mem;
1198 
1199   N_VScale(ONE, cv_mem->cv_acor, ele);
1200 
1201   return(CV_SUCCESS);
1202 }
1203 
1204 /*
1205  * CVodeGetWorkSpace
1206  *
1207  * Returns integrator work space requirements
1208  */
1209 
CVodeGetWorkSpace(void * cvode_mem,long int * lenrw,long int * leniw)1210 int CVodeGetWorkSpace(void *cvode_mem, long int *lenrw, long int *leniw)
1211 {
1212   CVodeMem cv_mem;
1213 
1214   if (cvode_mem==NULL) {
1215     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetWorkSpace", MSGCV_NO_MEM);
1216     return(CV_MEM_NULL);
1217   }
1218 
1219   cv_mem = (CVodeMem) cvode_mem;
1220 
1221   *leniw = cv_mem->cv_liw;
1222   *lenrw = cv_mem->cv_lrw;
1223 
1224   return(CV_SUCCESS);
1225 }
1226 
1227 /*
1228  * CVodeGetIntegratorStats
1229  *
1230  * Returns integrator statistics
1231  */
1232 
CVodeGetIntegratorStats(void * cvode_mem,long int * nsteps,long int * nfevals,long int * nlinsetups,long int * netfails,int * qlast,int * qcur,realtype * hinused,realtype * hlast,realtype * hcur,realtype * tcur)1233 int CVodeGetIntegratorStats(void *cvode_mem, long int *nsteps, long int *nfevals,
1234                             long int *nlinsetups, long int *netfails, int *qlast,
1235                             int *qcur, realtype *hinused, realtype *hlast,
1236                             realtype *hcur, realtype *tcur)
1237 {
1238   CVodeMem cv_mem;
1239 
1240   if (cvode_mem==NULL) {
1241     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetIntegratorStats", MSGCV_NO_MEM);
1242     return(CV_MEM_NULL);
1243   }
1244 
1245   cv_mem = (CVodeMem) cvode_mem;
1246 
1247   *nsteps = cv_mem->cv_nst;
1248   *nfevals = cv_mem->cv_nfe;
1249   *nlinsetups = cv_mem->cv_nsetups;
1250   *netfails = cv_mem->cv_netf;
1251   *qlast = cv_mem->cv_qu;
1252   *qcur = cv_mem->cv_next_q;
1253   *hinused = cv_mem->cv_h0u;
1254   *hlast = cv_mem->cv_hu;
1255   *hcur = cv_mem->cv_next_h;
1256   *tcur = cv_mem->cv_tn;
1257 
1258   return(CV_SUCCESS);
1259 }
1260 
1261 /*
1262  * CVodeGetNumGEvals
1263  *
1264  * Returns the current number of calls to g (for rootfinding)
1265  */
1266 
CVodeGetNumGEvals(void * cvode_mem,long int * ngevals)1267 int CVodeGetNumGEvals(void *cvode_mem, long int *ngevals)
1268 {
1269   CVodeMem cv_mem;
1270 
1271   if (cvode_mem==NULL) {
1272     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetNumGEvals", MSGCV_NO_MEM);
1273     return(CV_MEM_NULL);
1274   }
1275 
1276   cv_mem = (CVodeMem) cvode_mem;
1277 
1278   *ngevals = cv_mem->cv_nge;
1279 
1280   return(CV_SUCCESS);
1281 }
1282 
1283 /*
1284  * CVodeGetRootInfo
1285  *
1286  * Returns pointer to array rootsfound showing roots found
1287  */
1288 
CVodeGetRootInfo(void * cvode_mem,int * rootsfound)1289 int CVodeGetRootInfo(void *cvode_mem, int *rootsfound)
1290 {
1291   CVodeMem cv_mem;
1292   int i, nrt;
1293 
1294   if (cvode_mem==NULL) {
1295     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetRootInfo", MSGCV_NO_MEM);
1296     return(CV_MEM_NULL);
1297   }
1298 
1299   cv_mem = (CVodeMem) cvode_mem;
1300 
1301   nrt = cv_mem->cv_nrtfn;
1302 
1303   for (i=0; i<nrt; i++) rootsfound[i] = cv_mem->cv_iroots[i];
1304 
1305   return(CV_SUCCESS);
1306 }
1307 
1308 
1309 /*
1310  * CVodeGetNumNonlinSolvIters
1311  *
1312  * Returns the current number of iterations in the nonlinear solver
1313  */
1314 
CVodeGetNumNonlinSolvIters(void * cvode_mem,long int * nniters)1315 int CVodeGetNumNonlinSolvIters(void *cvode_mem, long int *nniters)
1316 {
1317   CVodeMem cv_mem;
1318   booleantype sensi_sim;
1319 
1320   if (cvode_mem==NULL) {
1321     cvProcessError(NULL, CV_MEM_NULL, "CVODES",
1322                    "CVodeGetNumNonlinSolvIters", MSGCV_NO_MEM);
1323     return(CV_MEM_NULL);
1324   }
1325 
1326   cv_mem = (CVodeMem) cvode_mem;
1327 
1328   /* are we computing sensitivities with the simultaneous approach? */
1329   sensi_sim = (cv_mem->cv_sensi && (cv_mem->cv_ism==CV_SIMULTANEOUS));
1330 
1331   /* get number of iterations from the NLS */
1332   if (sensi_sim) {
1333 
1334     /* check that the NLS is non-NULL */
1335     if (cv_mem->NLSsim == NULL) {
1336       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
1337                      "CVodeGetNumNonlinSolvIters", MSGCV_MEM_FAIL);
1338       return(CV_MEM_FAIL);
1339     }
1340 
1341     return(SUNNonlinSolGetNumIters(cv_mem->NLSsim, nniters));
1342 
1343   } else {
1344 
1345     /* check that the NLS is non-NULL */
1346     if (cv_mem->NLS == NULL) {
1347       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
1348                      "CVodeGetNumNonlinSolvIters", MSGCV_MEM_FAIL);
1349       return(CV_MEM_FAIL);
1350     }
1351 
1352     return(SUNNonlinSolGetNumIters(cv_mem->NLS, nniters));
1353   }
1354 
1355   return(CV_SUCCESS);
1356 }
1357 
1358 /*
1359  * CVodeGetNumNonlinSolvConvFails
1360  *
1361  * Returns the current number of convergence failures in the
1362  * nonlinear solver
1363  */
1364 
CVodeGetNumNonlinSolvConvFails(void * cvode_mem,long int * nncfails)1365 int CVodeGetNumNonlinSolvConvFails(void *cvode_mem, long int *nncfails)
1366 {
1367   CVodeMem cv_mem;
1368 
1369   if (cvode_mem==NULL) {
1370     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetNumNonlinSolvConvFails", MSGCV_NO_MEM);
1371     return(CV_MEM_NULL);
1372   }
1373 
1374   cv_mem = (CVodeMem) cvode_mem;
1375 
1376   *nncfails = cv_mem->cv_ncfn;
1377 
1378   return(CV_SUCCESS);
1379 }
1380 
1381 /*
1382  * CVodeGetNonlinSolvStats
1383  *
1384  * Returns nonlinear solver statistics
1385  */
1386 
CVodeGetNonlinSolvStats(void * cvode_mem,long int * nniters,long int * nncfails)1387 int CVodeGetNonlinSolvStats(void *cvode_mem, long int *nniters,
1388                             long int *nncfails)
1389 {
1390   CVodeMem cv_mem;
1391   booleantype sensi_sim;
1392 
1393   if (cvode_mem==NULL) {
1394     cvProcessError(NULL, CV_MEM_NULL, "CVODES",
1395                    "CVodeGetNonlinSolvStats", MSGCV_NO_MEM);
1396     return(CV_MEM_NULL);
1397   }
1398 
1399   cv_mem = (CVodeMem) cvode_mem;
1400 
1401   *nncfails = cv_mem->cv_ncfn;
1402 
1403   /* are we computing sensitivities with the simultaneous approach? */
1404   sensi_sim = (cv_mem->cv_sensi && (cv_mem->cv_ism==CV_SIMULTANEOUS));
1405 
1406   /* get number of iterations from the NLS */
1407   if (sensi_sim) {
1408 
1409     /* check that the NLS is non-NULL */
1410     if (cv_mem->NLSsim == NULL) {
1411       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
1412                      "CVodeGetNumNonlinSolvIters", MSGCV_MEM_FAIL);
1413       return(CV_MEM_FAIL);
1414     }
1415 
1416     return(SUNNonlinSolGetNumIters(cv_mem->NLSsim, nniters));
1417 
1418   } else {
1419 
1420     /* check that the NLS is non-NULL */
1421     if (cv_mem->NLS == NULL) {
1422       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
1423                      "CVodeGetNumNonlinSolvIters", MSGCV_MEM_FAIL);
1424       return(CV_MEM_FAIL);
1425     }
1426 
1427     return(SUNNonlinSolGetNumIters(cv_mem->NLS, nniters));
1428   }
1429 
1430   return(CV_SUCCESS);
1431 }
1432 
1433 
1434 /*
1435  * =================================================================
1436  * Quadrature optional output functions
1437  * =================================================================
1438  */
1439 
1440 /*-----------------------------------------------------------------*/
1441 
CVodeGetQuadNumRhsEvals(void * cvode_mem,long int * nfQevals)1442 int CVodeGetQuadNumRhsEvals(void *cvode_mem, long int *nfQevals)
1443 {
1444   CVodeMem cv_mem;
1445 
1446   if (cvode_mem==NULL) {
1447     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetQuadNumRhsEvals", MSGCV_NO_MEM);
1448     return(CV_MEM_NULL);
1449   }
1450 
1451   cv_mem = (CVodeMem) cvode_mem;
1452 
1453   if (cv_mem->cv_quadr==SUNFALSE) {
1454     cvProcessError(cv_mem, CV_NO_QUAD, "CVODES", "CVodeGetQuadNumRhsEvals", MSGCV_NO_QUAD);
1455     return(CV_NO_QUAD);
1456   }
1457 
1458   *nfQevals = cv_mem->cv_nfQe;
1459 
1460   return(CV_SUCCESS);
1461 }
1462 
1463 /*-----------------------------------------------------------------*/
1464 
CVodeGetQuadNumErrTestFails(void * cvode_mem,long int * nQetfails)1465 int CVodeGetQuadNumErrTestFails(void *cvode_mem, long int *nQetfails)
1466 {
1467   CVodeMem cv_mem;
1468 
1469   if (cvode_mem==NULL) {
1470     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetQuadNumErrTestFails", MSGCV_NO_MEM);
1471     return(CV_MEM_NULL);
1472   }
1473 
1474   cv_mem = (CVodeMem) cvode_mem;
1475 
1476   if (cv_mem->cv_quadr==SUNFALSE) {
1477     cvProcessError(cv_mem, CV_NO_QUAD, "CVODES", "CVodeGetQuadNumErrTestFails", MSGCV_NO_QUAD);
1478     return(CV_NO_QUAD);
1479   }
1480 
1481   *nQetfails = cv_mem->cv_netfQ;
1482 
1483   return(CV_SUCCESS);
1484 }
1485 
1486 /*-----------------------------------------------------------------*/
1487 
CVodeGetQuadErrWeights(void * cvode_mem,N_Vector eQweight)1488 int CVodeGetQuadErrWeights(void *cvode_mem, N_Vector eQweight)
1489 {
1490   CVodeMem cv_mem;
1491 
1492   if (cvode_mem==NULL) {
1493     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetQuadErrWeights", MSGCV_NO_MEM);
1494     return(CV_MEM_NULL);
1495   }
1496 
1497   cv_mem = (CVodeMem) cvode_mem;
1498 
1499   if (cv_mem->cv_quadr==SUNFALSE) {
1500     cvProcessError(cv_mem, CV_NO_QUAD, "CVODES", "CVodeGetQuadErrWeights", MSGCV_NO_QUAD);
1501     return(CV_NO_QUAD);
1502   }
1503 
1504   if(cv_mem->cv_errconQ) N_VScale(ONE, cv_mem->cv_ewtQ, eQweight);
1505 
1506   return(CV_SUCCESS);
1507 }
1508 
1509 /*-----------------------------------------------------------------*/
1510 
CVodeGetQuadStats(void * cvode_mem,long int * nfQevals,long int * nQetfails)1511 int CVodeGetQuadStats(void *cvode_mem, long int *nfQevals, long int *nQetfails)
1512 {
1513   CVodeMem cv_mem;
1514 
1515   if (cvode_mem==NULL) {
1516     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetQuadStats", MSGCV_NO_MEM);
1517     return(CV_MEM_NULL);
1518   }
1519 
1520   cv_mem = (CVodeMem) cvode_mem;
1521 
1522   if (cv_mem->cv_quadr==SUNFALSE) {
1523     cvProcessError(cv_mem, CV_NO_QUAD, "CVODES", "CVodeGetQuadStats", MSGCV_NO_QUAD);
1524     return(CV_NO_QUAD);
1525   }
1526 
1527   *nfQevals = cv_mem->cv_nfQe;
1528   *nQetfails = cv_mem->cv_netfQ;
1529 
1530   return(CV_SUCCESS);
1531 }
1532 
1533 /*
1534  * =================================================================
1535  * Quadrature FSA optional output functions
1536  * =================================================================
1537  */
1538 
1539 /*-----------------------------------------------------------------*/
1540 
CVodeGetQuadSensNumRhsEvals(void * cvode_mem,long int * nfQSevals)1541 int CVodeGetQuadSensNumRhsEvals(void *cvode_mem, long int *nfQSevals)
1542 {
1543   CVodeMem cv_mem;
1544 
1545   if (cvode_mem==NULL) {
1546     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetQuadSensNumRhsEvals", MSGCV_NO_MEM);
1547     return(CV_MEM_NULL);
1548   }
1549 
1550   cv_mem = (CVodeMem) cvode_mem;
1551 
1552   if (cv_mem->cv_quadr_sensi == SUNFALSE) {
1553     cvProcessError(cv_mem, CV_NO_QUADSENS, "CVODES", "CVodeGetQuadSensNumRhsEvals", MSGCV_NO_QUADSENSI);
1554     return(CV_NO_QUADSENS);
1555   }
1556 
1557   *nfQSevals = cv_mem->cv_nfQSe;
1558 
1559   return(CV_SUCCESS);
1560 }
1561 
1562 /*-----------------------------------------------------------------*/
1563 
CVodeGetQuadSensNumErrTestFails(void * cvode_mem,long int * nQSetfails)1564 int CVodeGetQuadSensNumErrTestFails(void *cvode_mem, long int *nQSetfails)
1565 {
1566   CVodeMem cv_mem;
1567 
1568   if (cvode_mem==NULL) {
1569     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetQuadSensNumErrTestFails", MSGCV_NO_MEM);
1570     return(CV_MEM_NULL);
1571   }
1572 
1573   cv_mem = (CVodeMem) cvode_mem;
1574 
1575   if (cv_mem->cv_quadr_sensi == SUNFALSE) {
1576     cvProcessError(cv_mem, CV_NO_QUADSENS, "CVODES", "CVodeGetQuadSensNumErrTestFails", MSGCV_NO_QUADSENSI);
1577     return(CV_NO_QUADSENS);
1578   }
1579 
1580   *nQSetfails = cv_mem->cv_netfQS;
1581 
1582   return(CV_SUCCESS);
1583 }
1584 
1585 /*-----------------------------------------------------------------*/
1586 
CVodeGetQuadSensErrWeights(void * cvode_mem,N_Vector * eQSweight)1587 int CVodeGetQuadSensErrWeights(void *cvode_mem, N_Vector *eQSweight)
1588 {
1589   CVodeMem cv_mem;
1590   int is, Ns;
1591 
1592   if (cvode_mem==NULL) {
1593     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetQuadSensErrWeights", MSGCV_NO_MEM);
1594     return(CV_MEM_NULL);
1595   }
1596 
1597   cv_mem = (CVodeMem) cvode_mem;
1598 
1599   if (cv_mem->cv_quadr_sensi == SUNFALSE) {
1600     cvProcessError(cv_mem, CV_NO_QUADSENS, "CVODES", "CVodeGetQuadSensErrWeights", MSGCV_NO_QUADSENSI);
1601     return(CV_NO_QUADSENS);
1602   }
1603   Ns = cv_mem->cv_Ns;
1604 
1605   if (cv_mem->cv_errconQS)
1606     for (is=0; is<Ns; is++)
1607       N_VScale(ONE, cv_mem->cv_ewtQS[is], eQSweight[is]);
1608 
1609   return(CV_SUCCESS);
1610 }
1611 
1612 /*-----------------------------------------------------------------*/
1613 
CVodeGetQuadSensStats(void * cvode_mem,long int * nfQSevals,long int * nQSetfails)1614 int CVodeGetQuadSensStats(void *cvode_mem, long int *nfQSevals, long int *nQSetfails)
1615 {
1616   CVodeMem cv_mem;
1617 
1618   if (cvode_mem==NULL) {
1619     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetQuadSensStats", MSGCV_NO_MEM);
1620     return(CV_MEM_NULL);
1621   }
1622 
1623   cv_mem = (CVodeMem) cvode_mem;
1624 
1625   if (cv_mem->cv_quadr_sensi == SUNFALSE) {
1626     cvProcessError(cv_mem, CV_NO_QUADSENS, "CVODES", "CVodeGetQuadSensStats", MSGCV_NO_QUADSENSI);
1627     return(CV_NO_QUADSENS);
1628   }
1629 
1630   *nfQSevals = cv_mem->cv_nfQSe;
1631   *nQSetfails = cv_mem->cv_netfQS;
1632 
1633   return(CV_SUCCESS);
1634 }
1635 
1636 
1637 /*
1638  * =================================================================
1639  * FSA optional output functions
1640  * =================================================================
1641  */
1642 
1643 /*-----------------------------------------------------------------*/
1644 
CVodeGetSensNumRhsEvals(void * cvode_mem,long int * nfSevals)1645 int CVodeGetSensNumRhsEvals(void *cvode_mem, long int *nfSevals)
1646 {
1647   CVodeMem cv_mem;
1648 
1649   if (cvode_mem==NULL) {
1650     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetSensNumRhsEvals", MSGCV_NO_MEM);
1651     return(CV_MEM_NULL);
1652   }
1653 
1654   cv_mem = (CVodeMem) cvode_mem;
1655 
1656   if (cv_mem->cv_sensi==SUNFALSE) {
1657     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetSensNumRhsEvals", MSGCV_NO_SENSI);
1658     return(CV_NO_SENS);
1659   }
1660 
1661   *nfSevals = cv_mem->cv_nfSe;
1662 
1663   return(CV_SUCCESS);
1664 }
1665 
1666 /*-----------------------------------------------------------------*/
1667 
CVodeGetNumRhsEvalsSens(void * cvode_mem,long int * nfevalsS)1668 int CVodeGetNumRhsEvalsSens(void *cvode_mem, long int *nfevalsS)
1669 {
1670   CVodeMem cv_mem;
1671 
1672   if (cvode_mem==NULL) {
1673     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetNumRhsEvalsSens", MSGCV_NO_MEM);
1674     return(CV_MEM_NULL);
1675   }
1676 
1677   cv_mem = (CVodeMem) cvode_mem;
1678 
1679   if (cv_mem->cv_sensi==SUNFALSE) {
1680     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetNumRhsEvalsSens", MSGCV_NO_SENSI);
1681     return(CV_NO_SENS);
1682   }
1683 
1684   *nfevalsS = cv_mem->cv_nfeS;
1685 
1686   return(CV_SUCCESS);
1687 }
1688 
1689 /*-----------------------------------------------------------------*/
1690 
CVodeGetSensNumErrTestFails(void * cvode_mem,long int * nSetfails)1691 int CVodeGetSensNumErrTestFails(void *cvode_mem, long int *nSetfails)
1692 {
1693   CVodeMem cv_mem;
1694 
1695   if (cvode_mem==NULL) {
1696     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetSensNumErrTestFails", MSGCV_NO_MEM);
1697     return(CV_MEM_NULL);
1698   }
1699 
1700   cv_mem = (CVodeMem) cvode_mem;
1701 
1702   if (cv_mem->cv_sensi==SUNFALSE) {
1703     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetSensNumErrTestFails", MSGCV_NO_SENSI);
1704     return(CV_NO_SENS);
1705   }
1706 
1707   *nSetfails = cv_mem->cv_netfS;
1708 
1709   return(CV_SUCCESS);
1710 }
1711 
1712 /*-----------------------------------------------------------------*/
1713 
CVodeGetSensNumLinSolvSetups(void * cvode_mem,long int * nlinsetupsS)1714 int CVodeGetSensNumLinSolvSetups(void *cvode_mem, long int *nlinsetupsS)
1715 {
1716   CVodeMem cv_mem;
1717 
1718   if (cvode_mem==NULL) {
1719     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetSensNumLinSolvSetups", MSGCV_NO_MEM);
1720     return(CV_MEM_NULL);
1721   }
1722 
1723   cv_mem = (CVodeMem) cvode_mem;
1724 
1725   if (cv_mem->cv_sensi==SUNFALSE) {
1726     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetSensNumLinSolvSetups", MSGCV_NO_SENSI);
1727     return(CV_NO_SENS);
1728   }
1729 
1730   *nlinsetupsS = cv_mem->cv_nsetupsS;
1731 
1732   return(CV_SUCCESS);
1733 }
1734 
1735 /*-----------------------------------------------------------------*/
1736 
CVodeGetSensErrWeights(void * cvode_mem,N_Vector * eSweight)1737 int CVodeGetSensErrWeights(void *cvode_mem, N_Vector *eSweight)
1738 {
1739   CVodeMem cv_mem;
1740   int is, Ns;
1741 
1742   if (cvode_mem==NULL) {
1743     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetSensErrWeights", MSGCV_NO_MEM);
1744     return(CV_MEM_NULL);
1745   }
1746 
1747   cv_mem = (CVodeMem) cvode_mem;
1748 
1749   if (cv_mem->cv_sensi==SUNFALSE) {
1750     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetSensErrWeights", MSGCV_NO_SENSI);
1751     return(CV_NO_SENS);
1752   }
1753 
1754   Ns = cv_mem->cv_Ns;
1755 
1756   for (is=0; is<Ns; is++)
1757     N_VScale(ONE, cv_mem->cv_ewtS[is], eSweight[is]);
1758 
1759   return(CV_SUCCESS);
1760 }
1761 
1762 /*-----------------------------------------------------------------*/
1763 
CVodeGetSensStats(void * cvode_mem,long int * nfSevals,long int * nfevalsS,long int * nSetfails,long int * nlinsetupsS)1764 int CVodeGetSensStats(void *cvode_mem, long int *nfSevals, long int *nfevalsS,
1765                       long int *nSetfails, long int *nlinsetupsS)
1766 {
1767   CVodeMem cv_mem;
1768 
1769   if (cvode_mem==NULL) {
1770     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetSensStats", MSGCV_NO_MEM);
1771     return(CV_MEM_NULL);
1772   }
1773 
1774   cv_mem = (CVodeMem) cvode_mem;
1775 
1776   if (cv_mem->cv_sensi==SUNFALSE) {
1777     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetSensStats", MSGCV_NO_SENSI);
1778     return(CV_NO_SENS);
1779   }
1780 
1781   *nfSevals = cv_mem->cv_nfSe;
1782   *nfevalsS = cv_mem->cv_nfeS;
1783   *nSetfails = cv_mem->cv_netfS;
1784   *nlinsetupsS = cv_mem->cv_nsetupsS;
1785 
1786   return(CV_SUCCESS);
1787 }
1788 
1789 /*-----------------------------------------------------------------*/
1790 
CVodeGetSensNumNonlinSolvIters(void * cvode_mem,long int * nSniters)1791 int CVodeGetSensNumNonlinSolvIters(void *cvode_mem, long int *nSniters)
1792 {
1793   CVodeMem cv_mem;
1794   booleantype sensi_stg;
1795 
1796   if (cvode_mem==NULL) {
1797     cvProcessError(NULL, CV_MEM_NULL, "CVODES",
1798                    "CVodeGetSensNumNonlinSolvIters", MSGCV_NO_MEM);
1799     return(CV_MEM_NULL);
1800   }
1801 
1802   cv_mem = (CVodeMem) cvode_mem;
1803 
1804   if (cv_mem->cv_sensi==SUNFALSE) {
1805     cvProcessError(cv_mem, CV_NO_SENS, "CVODES",
1806                    "CVodeGetSensNumNonlinSolvIters", MSGCV_NO_SENSI);
1807     return(CV_NO_SENS);
1808   }
1809 
1810   /* Are we computing sensitivities with a staggered approach? */
1811   sensi_stg = (cv_mem->cv_sensi && (cv_mem->cv_ism==CV_STAGGERED));
1812 
1813   if (sensi_stg) {
1814 
1815     /* check that the NLS is non-NULL */
1816     if (cv_mem->NLSstg == NULL) {
1817       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
1818                       "CVodeGetSensNumNonlinSolvIters", MSGCV_MEM_FAIL);
1819       return(CV_MEM_FAIL);
1820     }
1821 
1822     return(SUNNonlinSolGetNumIters(cv_mem->NLSstg, nSniters));
1823 
1824   } else {
1825 
1826     /* check that the NLS is non-NULL */
1827     if (cv_mem->NLSstg1 == NULL) {
1828       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
1829                       "CVodeGetSensNumNonlinSolvIters", MSGCV_MEM_FAIL);
1830       return(CV_MEM_FAIL);
1831     }
1832 
1833     return(SUNNonlinSolGetNumIters(cv_mem->NLSstg1, nSniters));
1834   }
1835 
1836 }
1837 
1838 /*-----------------------------------------------------------------*/
1839 
CVodeGetSensNumNonlinSolvConvFails(void * cvode_mem,long int * nSncfails)1840 int CVodeGetSensNumNonlinSolvConvFails(void *cvode_mem, long int *nSncfails)
1841 {
1842   CVodeMem cv_mem;
1843 
1844   if (cvode_mem==NULL) {
1845     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetSensNumNonlinSolvConvFails", MSGCV_NO_MEM);
1846     return(CV_MEM_NULL);
1847   }
1848 
1849   cv_mem = (CVodeMem) cvode_mem;
1850 
1851   if (cv_mem->cv_sensi==SUNFALSE) {
1852     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetSensNumNonlinSolvConvFails", MSGCV_NO_SENSI);
1853     return(CV_NO_SENS);
1854   }
1855 
1856   *nSncfails = cv_mem->cv_ncfnS;
1857 
1858   return(CV_SUCCESS);
1859 }
1860 
1861 /*-----------------------------------------------------------------*/
1862 
CVodeGetStgrSensNumNonlinSolvIters(void * cvode_mem,long int * nSTGR1niters)1863 int CVodeGetStgrSensNumNonlinSolvIters(void *cvode_mem, long int *nSTGR1niters)
1864 {
1865   CVodeMem cv_mem;
1866   int is, Ns;
1867 
1868   if (cvode_mem==NULL) {
1869     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetStgrSensNumNonlinSolvIters", MSGCV_NO_MEM);
1870     return(CV_MEM_NULL);
1871   }
1872 
1873   cv_mem = (CVodeMem) cvode_mem;
1874 
1875   Ns = cv_mem->cv_Ns;
1876 
1877   if (cv_mem->cv_sensi==SUNFALSE) {
1878     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetStgrSensNumNonlinSolvIters", MSGCV_NO_SENSI);
1879     return(CV_NO_SENS);
1880   }
1881 
1882   if(cv_mem->cv_ism==CV_STAGGERED1)
1883     for(is=0; is<Ns; is++) nSTGR1niters[is] = cv_mem->cv_nniS1[is];
1884 
1885   return(CV_SUCCESS);
1886 }
1887 
1888 /*-----------------------------------------------------------------*/
1889 
CVodeGetStgrSensNumNonlinSolvConvFails(void * cvode_mem,long int * nSTGR1ncfails)1890 int CVodeGetStgrSensNumNonlinSolvConvFails(void *cvode_mem, long int *nSTGR1ncfails)
1891 {
1892   CVodeMem cv_mem;
1893   int is, Ns;
1894 
1895   if (cvode_mem==NULL) {
1896     cvProcessError(NULL, CV_MEM_NULL, "CVODES", "CVodeGetStgrSensNumNonlinSolvConvFails", MSGCV_NO_MEM);
1897     return(CV_MEM_NULL);
1898   }
1899 
1900   cv_mem = (CVodeMem) cvode_mem;
1901 
1902   Ns = cv_mem->cv_Ns;
1903 
1904   if (cv_mem->cv_sensi==SUNFALSE) {
1905     cvProcessError(cv_mem, CV_NO_SENS, "CVODES", "CVodeGetStgrSensNumNonlinSolvConvFails", MSGCV_NO_SENSI);
1906     return(CV_NO_SENS);
1907   }
1908 
1909   if(cv_mem->cv_ism==CV_STAGGERED1)
1910     for(is=0; is<Ns; is++) nSTGR1ncfails[is] = cv_mem->cv_ncfnS1[is];
1911 
1912   return(CV_SUCCESS);
1913 }
1914 
1915 /*-----------------------------------------------------------------*/
1916 
CVodeGetSensNonlinSolvStats(void * cvode_mem,long int * nSniters,long int * nSncfails)1917 int CVodeGetSensNonlinSolvStats(void *cvode_mem, long int *nSniters,
1918                                 long int *nSncfails)
1919 {
1920   CVodeMem cv_mem;
1921   booleantype sensi_stg;
1922 
1923   if (cvode_mem==NULL) {
1924     cvProcessError(NULL, CV_MEM_NULL, "CVODES",
1925                    "CVodeGetSensNonlinSolvstats", MSGCV_NO_MEM);
1926     return(CV_MEM_NULL);
1927   }
1928 
1929   cv_mem = (CVodeMem) cvode_mem;
1930 
1931   if (cv_mem->cv_sensi==SUNFALSE) {
1932     cvProcessError(cv_mem, CV_NO_SENS, "CVODES",
1933                    "CVodeGetSensNonlinSolvStats", MSGCV_NO_SENSI);
1934     return(CV_NO_SENS);
1935   }
1936 
1937   *nSncfails = cv_mem->cv_ncfnS;
1938 
1939   /* Are we computing sensitivities with a staggered approach? */
1940   sensi_stg = (cv_mem->cv_sensi && (cv_mem->cv_ism==CV_STAGGERED));
1941 
1942   if (sensi_stg) {
1943 
1944     /* check that the NLS is non-NULL */
1945     if (cv_mem->NLSstg == NULL) {
1946       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
1947                       "CVodeGetSensNumNonlinSolvStats", MSGCV_MEM_FAIL);
1948       return(CV_MEM_FAIL);
1949     }
1950 
1951     return(SUNNonlinSolGetNumIters(cv_mem->NLSstg, nSniters));
1952 
1953   } else {
1954 
1955     /* check that the NLS is non-NULL */
1956     if (cv_mem->NLSstg1 == NULL) {
1957       cvProcessError(NULL, CV_MEM_FAIL, "CVODES",
1958                       "CVodeGetSensNumNonlinSolvStats", MSGCV_MEM_FAIL);
1959       return(CV_MEM_FAIL);
1960     }
1961 
1962     return(SUNNonlinSolGetNumIters(cv_mem->NLSstg1, nSniters));
1963   }
1964 
1965   return(CV_SUCCESS);
1966 }
1967 
1968 /*-----------------------------------------------------------------*/
1969 
CVodeGetReturnFlagName(long int flag)1970 char *CVodeGetReturnFlagName(long int flag)
1971 {
1972   char *name;
1973 
1974   name = (char *)malloc(24*sizeof(char));
1975 
1976   switch(flag) {
1977   case CV_SUCCESS:
1978     sprintf(name,"CV_SUCCESS");
1979     break;
1980   case CV_TSTOP_RETURN:
1981     sprintf(name,"CV_TSTOP_RETURN");
1982     break;
1983   case CV_ROOT_RETURN:
1984     sprintf(name,"CV_ROOT_RETURN");
1985     break;
1986   case CV_TOO_MUCH_WORK:
1987     sprintf(name,"CV_TOO_MUCH_WORK");
1988     break;
1989   case CV_TOO_MUCH_ACC:
1990     sprintf(name,"CV_TOO_MUCH_ACC");
1991     break;
1992   case CV_ERR_FAILURE:
1993     sprintf(name,"CV_ERR_FAILURE");
1994     break;
1995   case CV_CONV_FAILURE:
1996     sprintf(name,"CV_CONV_FAILURE");
1997     break;
1998   case CV_LINIT_FAIL:
1999     sprintf(name,"CV_LINIT_FAIL");
2000     break;
2001   case CV_LSETUP_FAIL:
2002     sprintf(name,"CV_LSETUP_FAIL");
2003     break;
2004   case CV_LSOLVE_FAIL:
2005     sprintf(name,"CV_LSOLVE_FAIL");
2006     break;
2007   case CV_RHSFUNC_FAIL:
2008     sprintf(name,"CV_RHSFUNC_FAIL");
2009     break;
2010   case CV_FIRST_RHSFUNC_ERR:
2011     sprintf(name,"CV_FIRST_RHSFUNC_ERR");
2012     break;
2013   case CV_REPTD_RHSFUNC_ERR:
2014     sprintf(name,"CV_REPTD_RHSFUNC_ERR");
2015     break;
2016   case CV_UNREC_RHSFUNC_ERR:
2017     sprintf(name,"CV_UNREC_RHSFUNC_ERR");
2018     break;
2019   case CV_RTFUNC_FAIL:
2020     sprintf(name,"CV_RTFUNC_FAIL");
2021     break;
2022   case CV_MEM_FAIL:
2023     sprintf(name,"CV_MEM_FAIL");
2024     break;
2025   case CV_MEM_NULL:
2026     sprintf(name,"CV_MEM_NULL");
2027     break;
2028   case CV_ILL_INPUT:
2029     sprintf(name,"CV_ILL_INPUT");
2030     break;
2031   case CV_NO_MALLOC:
2032     sprintf(name,"CV_NO_MALLOC");
2033     break;
2034   case CV_BAD_K:
2035     sprintf(name,"CV_BAD_K");
2036     break;
2037   case CV_BAD_T:
2038     sprintf(name,"CV_BAD_T");
2039     break;
2040   case CV_BAD_DKY:
2041     sprintf(name,"CV_BAD_DKY");
2042     break;
2043   case CV_NO_QUAD:
2044     sprintf(name,"CV_NO_QUAD");
2045     break;
2046   case CV_QRHSFUNC_FAIL:
2047     sprintf(name,"CV_QRHSFUNC_FAIL");
2048     break;
2049   case CV_FIRST_QRHSFUNC_ERR:
2050     sprintf(name,"CV_FIRST_QRHSFUNC_ERR");
2051     break;
2052   case CV_REPTD_QRHSFUNC_ERR:
2053     sprintf(name,"CV_REPTD_QRHSFUNC_ERR");
2054     break;
2055   case CV_UNREC_QRHSFUNC_ERR:
2056     sprintf(name,"CV_UNREC_QRHSFUNC_ERR");
2057     break;
2058   case CV_BAD_IS:
2059     sprintf(name,"CV_BAD_IS");
2060     break;
2061   case CV_NO_SENS:
2062     sprintf(name,"CV_NO_SENS");
2063     break;
2064   case CV_SRHSFUNC_FAIL:
2065     sprintf(name,"CV_SRHSFUNC_FAIL");
2066     break;
2067   case CV_FIRST_SRHSFUNC_ERR:
2068     sprintf(name,"CV_FIRST_SRHSFUNC_ERR");
2069     break;
2070   case CV_REPTD_SRHSFUNC_ERR:
2071     sprintf(name,"CV_REPTD_SRHSFUNC_ERR");
2072     break;
2073   case CV_UNREC_SRHSFUNC_ERR:
2074     sprintf(name,"CV_UNREC_SRHSFUNC_ERR");
2075     break;
2076   case CV_TOO_CLOSE:
2077     sprintf(name,"CV_TOO_CLOSE");
2078     break;
2079   case CV_NO_ADJ:
2080     sprintf(name,"CV_NO_ADJ");
2081     break;
2082   case CV_NO_FWD:
2083     sprintf(name,"CV_NO_FWD");
2084     break;
2085   case CV_NO_BCK:
2086     sprintf(name,"CV_NO_BCK");
2087     break;
2088   case CV_BAD_TB0:
2089     sprintf(name,"CV_BAD_TB0");
2090     break;
2091   case CV_REIFWD_FAIL:
2092     sprintf(name,"CV_REIFWD_FAIL");
2093     break;
2094   case CV_FWD_FAIL:
2095     sprintf(name,"CV_FWD_FAIL");
2096     break;
2097   case CV_GETY_BADT:
2098     sprintf(name,"CV_GETY_BADT");
2099     break;
2100   case CV_NLS_FAIL:
2101     sprintf(name,"CV_NLS_FAIL");
2102     break;
2103   default:
2104     sprintf(name,"NONE");
2105   }
2106 
2107   return(name);
2108 }
2109