1 /*
2  * -----------------------------------------------------------------
3  * $Revision$
4  * $Date$
5  * -----------------------------------------------------------------
6  * Programmer: Radu Serban @ LLNL
7  * -----------------------------------------------------------------
8  * SUNDIALS Copyright Start
9  * Copyright (c) 2002-2021, 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 adjoint module in the CVODES solver.
20  * -----------------------------------------------------------------
21  */
22 
23 /*
24  * =================================================================
25  * IMPORTED HEADER FILES
26  * =================================================================
27  */
28 
29 #include <stdio.h>
30 #include <stdlib.h>
31 
32 #include "cvodes_impl.h"
33 #include <sundials/sundials_types.h>
34 
35 /*
36  * =================================================================
37  * CVODEA PRIVATE CONSTANTS
38  * =================================================================
39  */
40 
41 #define ONE         RCONST(1.0)
42 
43 /*
44  * =================================================================
45  * EXPORTED FUNCTIONS IMPLEMENTATION
46  * =================================================================
47  */
48 
49 /*
50  * -----------------------------------------------------------------
51  * Optional input functions for ASA
52  * -----------------------------------------------------------------
53  */
54 
CVodeSetAdjNoSensi(void * cvode_mem)55 int CVodeSetAdjNoSensi(void *cvode_mem)
56 {
57   CVodeMem cv_mem;
58   CVadjMem ca_mem;
59 
60   /* Check if cvode_mem exists */
61   if (cvode_mem == NULL) {
62     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetAdjNoSensi", MSGCV_NO_MEM);
63     return(CV_MEM_NULL);
64   }
65   cv_mem = (CVodeMem) cvode_mem;
66 
67   /* Was ASA initialized? */
68   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
69     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetAdjNoSensi", MSGCV_NO_ADJ);
70     return(CV_NO_ADJ);
71   }
72   ca_mem = cv_mem->cv_adj_mem;
73 
74   ca_mem->ca_IMstoreSensi = SUNFALSE;
75 
76   return(CV_SUCCESS);
77 }
78 
79 /*
80  * -----------------------------------------------------------------
81  * Optional input functions for backward integration
82  * -----------------------------------------------------------------
83  */
84 
CVodeSetNonlinearSolverB(void * cvode_mem,int which,SUNNonlinearSolver NLS)85 int CVodeSetNonlinearSolverB(void *cvode_mem, int which, SUNNonlinearSolver NLS)
86 {
87   CVodeMem cv_mem;
88   CVadjMem ca_mem;
89   CVodeBMem cvB_mem;
90   void *cvodeB_mem;
91 
92   /* Check if cvode_mem exists */
93   if (cvode_mem == NULL) {
94     cvProcessError(NULL, CV_MEM_NULL, "CVODEA",
95                    "CVodeSetNonlinearSolverB", MSGCV_NO_MEM);
96     return(CV_MEM_NULL);
97   }
98   cv_mem = (CVodeMem) cvode_mem;
99 
100   /* Was ASA initialized? */
101   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
102     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA",
103                    "CVodeSetNonlinearSolverB", MSGCV_NO_ADJ);
104     return(CV_NO_ADJ);
105   }
106   ca_mem = cv_mem->cv_adj_mem;
107 
108   /* Check which */
109   if ( which >= ca_mem->ca_nbckpbs ) {
110     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA",
111                    "CVodeSetNonlinearSolverB", MSGCV_BAD_WHICH);
112     return(CV_ILL_INPUT);
113   }
114 
115   /* Find the CVodeBMem entry in the linked list corresponding to which */
116   cvB_mem = ca_mem->cvB_mem;
117   while (cvB_mem != NULL) {
118     if ( which == cvB_mem->cv_index ) break;
119     cvB_mem = cvB_mem->cv_next;
120   }
121 
122   cvodeB_mem = (void *) (cvB_mem->cv_mem);
123 
124   return(CVodeSetNonlinearSolver(cvodeB_mem, NLS));
125 }
126 
CVodeSetUserDataB(void * cvode_mem,int which,void * user_dataB)127 int CVodeSetUserDataB(void *cvode_mem, int which, void *user_dataB)
128 {
129   CVodeMem cv_mem;
130   CVadjMem ca_mem;
131   CVodeBMem cvB_mem;
132 
133   /* Check if cvode_mem exists */
134   if (cvode_mem == NULL) {
135     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetUserDataB", MSGCV_NO_MEM);
136     return(CV_MEM_NULL);
137   }
138   cv_mem = (CVodeMem) cvode_mem;
139 
140   /* Was ASA initialized? */
141   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
142     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetUserDataB", MSGCV_NO_ADJ);
143     return(CV_NO_ADJ);
144   }
145   ca_mem = cv_mem->cv_adj_mem;
146 
147   /* Check which */
148   if ( which >= ca_mem->ca_nbckpbs ) {
149     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetUserDataB", MSGCV_BAD_WHICH);
150     return(CV_ILL_INPUT);
151   }
152 
153   /* Find the CVodeBMem entry in the linked list corresponding to which */
154   cvB_mem = ca_mem->cvB_mem;
155   while (cvB_mem != NULL) {
156     if ( which == cvB_mem->cv_index ) break;
157     cvB_mem = cvB_mem->cv_next;
158   }
159 
160   cvB_mem->cv_user_data = user_dataB;
161 
162   return(CV_SUCCESS);
163 }
164 
CVodeSetMaxOrdB(void * cvode_mem,int which,int maxordB)165 int CVodeSetMaxOrdB(void *cvode_mem, int which, int maxordB)
166 {
167   CVodeMem cv_mem;
168   CVadjMem ca_mem;
169   CVodeBMem cvB_mem;
170   void *cvodeB_mem;
171   int flag;
172 
173 
174   /* Check if cvode_mem exists */
175   if (cvode_mem == NULL) {
176     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetMaxOrdB", MSGCV_NO_MEM);
177     return(CV_MEM_NULL);
178   }
179   cv_mem = (CVodeMem) cvode_mem;
180 
181   /* Was ASA initialized? */
182   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
183     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetMaxOrdB", MSGCV_NO_ADJ);
184     return(CV_NO_ADJ);
185   }
186   ca_mem = cv_mem->cv_adj_mem;
187 
188   /* Check which */
189   if ( which >= ca_mem->ca_nbckpbs ) {
190     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetMaxOrdB", MSGCV_BAD_WHICH);
191     return(CV_ILL_INPUT);
192   }
193 
194   /* Find the CVodeBMem entry in the linked list corresponding to which */
195   cvB_mem = ca_mem->cvB_mem;
196   while (cvB_mem != NULL) {
197     if ( which == cvB_mem->cv_index ) break;
198     cvB_mem = cvB_mem->cv_next;
199   }
200 
201   cvodeB_mem = (void *) (cvB_mem->cv_mem);
202 
203   flag = CVodeSetMaxOrd(cvodeB_mem, maxordB);
204 
205   return(flag);
206 }
207 
208 
CVodeSetMaxNumStepsB(void * cvode_mem,int which,long int mxstepsB)209 int CVodeSetMaxNumStepsB(void *cvode_mem, int which, long int mxstepsB)
210 {
211   CVodeMem cv_mem;
212   CVadjMem ca_mem;
213   CVodeBMem cvB_mem;
214   void *cvodeB_mem;
215   int flag;
216 
217   /* Check if cvode_mem exists */
218   if (cvode_mem == NULL) {
219     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetMaxNumStepsB", MSGCV_NO_MEM);
220     return(CV_MEM_NULL);
221   }
222   cv_mem = (CVodeMem) cvode_mem;
223 
224   /* Was ASA initialized? */
225   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
226     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetMaxNumStepsB", MSGCV_NO_ADJ);
227     return(CV_NO_ADJ);
228   }
229   ca_mem = cv_mem->cv_adj_mem;
230 
231   /* Check which */
232   if ( which >= ca_mem->ca_nbckpbs ) {
233     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetMaxNumStepsB", MSGCV_BAD_WHICH);
234     return(CV_ILL_INPUT);
235   }
236 
237   /* Find the CVodeBMem entry in the linked list corresponding to which */
238   cvB_mem = ca_mem->cvB_mem;
239   while (cvB_mem != NULL) {
240     if ( which == cvB_mem->cv_index ) break;
241     cvB_mem = cvB_mem->cv_next;
242   }
243 
244   cvodeB_mem = (void *) (cvB_mem->cv_mem);
245 
246   flag = CVodeSetMaxNumSteps(cvodeB_mem, mxstepsB);
247 
248   return(flag);
249 }
250 
CVodeSetStabLimDetB(void * cvode_mem,int which,booleantype stldetB)251 int CVodeSetStabLimDetB(void *cvode_mem, int which, booleantype stldetB)
252 {
253   CVodeMem cv_mem;
254   CVadjMem ca_mem;
255   CVodeBMem cvB_mem;
256   void *cvodeB_mem;
257   int flag;
258 
259   /* Check if cvode_mem exists */
260   if (cvode_mem == NULL) {
261     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetStabLimDetB", MSGCV_NO_MEM);
262     return(CV_MEM_NULL);
263   }
264   cv_mem = (CVodeMem) cvode_mem;
265 
266   /* Was ASA initialized? */
267   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
268     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetStabLimDetB", MSGCV_NO_ADJ);
269     return(CV_NO_ADJ);
270   }
271   ca_mem = cv_mem->cv_adj_mem;
272 
273   /* Check which */
274   if ( which >= ca_mem->ca_nbckpbs ) {
275     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetStabLimDetB", MSGCV_BAD_WHICH);
276     return(CV_ILL_INPUT);
277   }
278 
279   /* Find the CVodeBMem entry in the linked list corresponding to which */
280   cvB_mem = ca_mem->cvB_mem;
281   while (cvB_mem != NULL) {
282     if ( which == cvB_mem->cv_index ) break;
283     cvB_mem = cvB_mem->cv_next;
284   }
285 
286   cvodeB_mem = (void *) (cvB_mem->cv_mem);
287 
288   flag = CVodeSetStabLimDet(cvodeB_mem, stldetB);
289 
290   return(flag);
291 }
292 
CVodeSetInitStepB(void * cvode_mem,int which,realtype hinB)293 int CVodeSetInitStepB(void *cvode_mem, int which, realtype hinB)
294 {
295   CVodeMem cv_mem;
296   CVadjMem ca_mem;
297   CVodeBMem cvB_mem;
298   void *cvodeB_mem;
299   int flag;
300 
301   /* Check if cvode_mem exists */
302   if (cvode_mem == NULL) {
303     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetInitStepB", MSGCV_NO_MEM);
304     return(CV_MEM_NULL);
305   }
306   cv_mem = (CVodeMem) cvode_mem;
307 
308   /* Was ASA initialized? */
309   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
310     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetInitStepB", MSGCV_NO_ADJ);
311     return(CV_NO_ADJ);
312   }
313   ca_mem = cv_mem->cv_adj_mem;
314 
315   /* Check which */
316   if ( which >= ca_mem->ca_nbckpbs ) {
317     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetInitStepB", MSGCV_BAD_WHICH);
318     return(CV_ILL_INPUT);
319   }
320 
321   /* Find the CVodeBMem entry in the linked list corresponding to which */
322   cvB_mem = ca_mem->cvB_mem;
323   while (cvB_mem != NULL) {
324     if ( which == cvB_mem->cv_index ) break;
325     cvB_mem = cvB_mem->cv_next;
326   }
327 
328   cvodeB_mem = (void *) (cvB_mem->cv_mem);
329 
330   flag = CVodeSetInitStep(cvodeB_mem, hinB);
331 
332   return(flag);
333 }
334 
CVodeSetMinStepB(void * cvode_mem,int which,realtype hminB)335 int CVodeSetMinStepB(void *cvode_mem, int which, realtype hminB)
336 {
337   CVodeMem cv_mem;
338   CVadjMem ca_mem;
339   CVodeBMem cvB_mem;
340   void *cvodeB_mem;
341   int flag;
342 
343   /* Check if cvode_mem exists */
344   if (cvode_mem == NULL) {
345     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetMinStepB", MSGCV_NO_MEM);
346     return(CV_MEM_NULL);
347   }
348   cv_mem = (CVodeMem) cvode_mem;
349 
350   /* Was ASA initialized? */
351   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
352     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetMinStepB", MSGCV_NO_ADJ);
353     return(CV_NO_ADJ);
354   }
355   ca_mem = cv_mem->cv_adj_mem;
356 
357   /* Check which */
358   if ( which >= ca_mem->ca_nbckpbs ) {
359     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetMinStepB", MSGCV_BAD_WHICH);
360     return(CV_ILL_INPUT);
361   }
362 
363   /* Find the CVodeBMem entry in the linked list corresponding to which */
364   cvB_mem = ca_mem->cvB_mem;
365   while (cvB_mem != NULL) {
366     if ( which == cvB_mem->cv_index ) break;
367     cvB_mem = cvB_mem->cv_next;
368   }
369 
370   cvodeB_mem = (void *) (cvB_mem->cv_mem);
371 
372   flag = CVodeSetMinStep(cvodeB_mem, hminB);
373 
374   return(flag);
375 }
376 
CVodeSetMaxStepB(void * cvode_mem,int which,realtype hmaxB)377 int CVodeSetMaxStepB(void *cvode_mem, int which, realtype hmaxB)
378 {
379   CVodeMem cv_mem;
380   CVadjMem ca_mem;
381   CVodeBMem cvB_mem;
382   void *cvodeB_mem;
383   int flag;
384 
385   /* Check if cvode_mem exists */
386   if (cvode_mem == NULL) {
387     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetMaxStepB", MSGCV_NO_MEM);
388     return(CV_MEM_NULL);
389   }
390   cv_mem = (CVodeMem) cvode_mem;
391 
392   /* Was ASA initialized? */
393   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
394     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetMaxStepB", MSGCV_NO_ADJ);
395     return(CV_NO_ADJ);
396   }
397   ca_mem = cv_mem->cv_adj_mem;
398 
399   /* Check which */
400   if ( which >= ca_mem->ca_nbckpbs ) {
401     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetMaxStepB", MSGCV_BAD_WHICH);
402     return(CV_ILL_INPUT);
403   }
404 
405   /* Find the CVodeBMem entry in the linked list corresponding to which */
406   cvB_mem = ca_mem->cvB_mem;
407   while (cvB_mem != NULL) {
408     if ( which == cvB_mem->cv_index ) break;
409     cvB_mem = cvB_mem->cv_next;
410   }
411 
412   cvodeB_mem = (void *) (cvB_mem->cv_mem);
413 
414   flag = CVodeSetMaxStep(cvodeB_mem, hmaxB);
415 
416   return(flag);
417 }
418 
CVodeSetConstraintsB(void * cvode_mem,int which,N_Vector constraintsB)419 int CVodeSetConstraintsB(void *cvode_mem, int which, N_Vector constraintsB)
420 {
421   CVodeMem cv_mem;
422   CVadjMem ca_mem;
423   CVodeBMem cvB_mem;
424   void *cvodeB_mem;
425   int flag;
426 
427   /* Is cvode_mem valid? */
428   if (cvode_mem == NULL) {
429     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetConstraintsB", MSGCV_NO_MEM);
430     return(CV_MEM_NULL);
431   }
432   cv_mem = (CVodeMem) cvode_mem;
433 
434   /* Is ASA initialized? */
435   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
436     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetConstraintsB", MSGCV_NO_ADJ);
437   }
438   ca_mem = cv_mem->cv_adj_mem;
439 
440   /* Check the value of which */
441   if ( which >= ca_mem->ca_nbckpbs ) {
442     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetConstraintsB", MSGCV_BAD_WHICH);
443     return(CV_ILL_INPUT);
444   }
445 
446   /* Find the CVodeBMem entry in the linked list corresponding to 'which'. */
447   cvB_mem = ca_mem->cvB_mem;
448   while (cvB_mem != NULL) {
449     if ( which == cvB_mem->cv_index) break;
450     /* advance */
451     cvB_mem = cvB_mem->cv_next;
452   }
453   cvodeB_mem = (void *) cvB_mem->cv_mem;
454 
455   flag = CVodeSetConstraints(cvodeB_mem, constraintsB);
456   return(flag);
457 }
458 
459 /*
460  * CVodeSetQuad*B
461  *
462  * Wrappers for the backward phase around the corresponding
463  * CVODES quadrature optional input functions
464  */
465 
CVodeSetQuadErrConB(void * cvode_mem,int which,booleantype errconQB)466 int CVodeSetQuadErrConB(void *cvode_mem, int which, booleantype errconQB)
467 {
468   CVodeMem cv_mem;
469   CVadjMem ca_mem;
470   CVodeBMem cvB_mem;
471   void *cvodeB_mem;
472   int flag;
473 
474   /* Check if cvode_mem exists */
475   if (cvode_mem == NULL) {
476     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeSetQuadErrConB", MSGCV_NO_MEM);
477     return(CV_MEM_NULL);
478   }
479   cv_mem = (CVodeMem) cvode_mem;
480 
481   /* Was ASA initialized? */
482   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
483     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeSetQuadErrConB", MSGCV_NO_ADJ);
484     return(CV_NO_ADJ);
485   }
486   ca_mem = cv_mem->cv_adj_mem;
487 
488   /* Check which */
489   if ( which >= ca_mem->ca_nbckpbs ) {
490     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVodeSetQuadErrConB", MSGCV_BAD_WHICH);
491     return(CV_ILL_INPUT);
492   }
493 
494   /* Find the CVodeBMem entry in the linked list corresponding to which */
495   cvB_mem = ca_mem->cvB_mem;
496   while (cvB_mem != NULL) {
497     if ( which == cvB_mem->cv_index ) break;
498     cvB_mem = cvB_mem->cv_next;
499   }
500 
501   cvodeB_mem = (void *) (cvB_mem->cv_mem);
502 
503   flag = CVodeSetQuadErrCon(cvodeB_mem, errconQB);
504 
505   return(flag);
506 }
507 
508 /*
509  * -----------------------------------------------------------------
510  * Optional output functions for backward integration
511  * -----------------------------------------------------------------
512  */
513 
514 /*
515  * CVodeGetAdjCVodeBmem
516  *
517  * This function returns a (void *) pointer to the CVODES
518  * memory allocated for the backward problem. This pointer can
519  * then be used to call any of the CVodeGet* CVODES routines to
520  * extract optional output for the backward integration phase.
521  */
522 
CVodeGetAdjCVodeBmem(void * cvode_mem,int which)523 void *CVodeGetAdjCVodeBmem(void *cvode_mem, int which)
524 {
525   CVodeMem cv_mem;
526   CVadjMem ca_mem;
527   CVodeBMem cvB_mem;
528   void *cvodeB_mem;
529 
530   /* Check if cvode_mem exists */
531   if (cvode_mem == NULL) {
532     cvProcessError(NULL, 0, "CVODEA", "CVodeGetAdjCVodeBmem", MSGCV_NO_MEM);
533     return(NULL);
534   }
535   cv_mem = (CVodeMem) cvode_mem;
536 
537   /* Was ASA initialized? */
538   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
539     cvProcessError(cv_mem, 0, "CVODEA", "CVodeGetAdjCVodeBmem", MSGCV_NO_ADJ);
540     return(NULL);
541   }
542   ca_mem = cv_mem->cv_adj_mem;
543 
544   /* Check which */
545   if ( which >= ca_mem->ca_nbckpbs ) {
546     cvProcessError(cv_mem, 0, "CVODEA", "CVodeGetAdjCVodeBmem", MSGCV_BAD_WHICH);
547     return(NULL);
548   }
549 
550   /* Find the CVodeBMem entry in the linked list corresponding to which */
551   cvB_mem = ca_mem->cvB_mem;
552   while (cvB_mem != NULL) {
553     if ( which == cvB_mem->cv_index ) break;
554     cvB_mem = cvB_mem->cv_next;
555   }
556 
557   cvodeB_mem = (void *) (cvB_mem->cv_mem);
558 
559   return(cvodeB_mem);
560 }
561 
562 /*
563  * CVodeGetAdjCheckPointsInfo
564  *
565  * This routine loads an array of nckpnts structures of type CVadjCheckPointRec.
566  * The user must allocate space for ckpnt.
567  */
568 
CVodeGetAdjCheckPointsInfo(void * cvode_mem,CVadjCheckPointRec * ckpnt)569 int CVodeGetAdjCheckPointsInfo(void *cvode_mem, CVadjCheckPointRec *ckpnt)
570 {
571   CVodeMem cv_mem;
572   CVadjMem ca_mem;
573   CkpntMem ck_mem;
574   int i;
575 
576   /* Check if cvode_mem exists */
577   if (cvode_mem == NULL) {
578     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeGetAdjCheckPointsInfo", MSGCV_NO_MEM);
579     return(CV_MEM_NULL);
580   }
581   cv_mem = (CVodeMem) cvode_mem;
582 
583   /* Was ASA initialized? */
584   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
585     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeGetAdjCheckPointsInfo", MSGCV_NO_ADJ);
586     return(CV_NO_ADJ);
587   }
588   ca_mem = cv_mem->cv_adj_mem;
589 
590   ck_mem = ca_mem->ck_mem;
591 
592   i = 0;
593 
594   while (ck_mem != NULL) {
595 
596     ckpnt[i].my_addr = (void *) ck_mem;
597     ckpnt[i].next_addr = (void *) ck_mem->ck_next;
598     ckpnt[i].t0 = ck_mem->ck_t0;
599     ckpnt[i].t1 = ck_mem->ck_t1;
600     ckpnt[i].nstep = ck_mem->ck_nst;
601     ckpnt[i].order = ck_mem->ck_q;
602     ckpnt[i].step = ck_mem->ck_h;
603 
604     ck_mem = ck_mem->ck_next;
605     i++;
606 
607   }
608 
609   return(CV_SUCCESS);
610 
611 }
612 
613 
614 /*
615  * -----------------------------------------------------------------
616  * Undocumented Development User-Callable Functions
617  * -----------------------------------------------------------------
618  */
619 
620 
621 /*
622  * CVodeGetAdjDataPointHermite
623  *
624  * This routine returns the solution stored in the data structure
625  * at the 'which' data point. Cubic Hermite interpolation.
626  */
627 
CVodeGetAdjDataPointHermite(void * cvode_mem,int which,realtype * t,N_Vector y,N_Vector yd)628 int CVodeGetAdjDataPointHermite(void *cvode_mem, int which,
629                                 realtype *t, N_Vector y, N_Vector yd)
630 {
631   CVodeMem cv_mem;
632   CVadjMem ca_mem;
633   DtpntMem *dt_mem;
634   HermiteDataMem content;
635 
636   /* Check if cvode_mem exists */
637   if (cvode_mem == NULL) {
638     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeGetAdjDataPointHermite", MSGCV_NO_MEM);
639     return(CV_MEM_NULL);
640   }
641   cv_mem = (CVodeMem) cvode_mem;
642 
643   /* Was ASA initialized? */
644   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
645     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeGetAdjDataPointHermite", MSGCV_NO_ADJ);
646     return(CV_NO_ADJ);
647   }
648   ca_mem = cv_mem->cv_adj_mem;
649 
650   dt_mem = ca_mem->dt_mem;
651 
652   if (ca_mem->ca_IMtype != CV_HERMITE) {
653     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVadjGetDataPointHermite", MSGCV_WRONG_INTERP);
654     return(CV_ILL_INPUT);
655   }
656 
657   *t = dt_mem[which]->t;
658 
659   content = (HermiteDataMem) (dt_mem[which]->content);
660 
661   if (y != NULL)
662     N_VScale(ONE, content->y, y);
663 
664   if (yd != NULL)
665     N_VScale(ONE, content->yd, yd);
666 
667   return(CV_SUCCESS);
668 }
669 
670 /*
671  * CVodeGetAdjDataPointPolynomial
672  *
673  * This routine returns the solution stored in the data structure
674  * at the 'which' data point. Polynomial interpolation.
675  */
676 
CVodeGetAdjDataPointPolynomial(void * cvode_mem,int which,realtype * t,int * order,N_Vector y)677 int CVodeGetAdjDataPointPolynomial(void *cvode_mem, int which,
678                                    realtype *t, int *order, N_Vector y)
679 {
680   CVodeMem cv_mem;
681   CVadjMem ca_mem;
682   DtpntMem *dt_mem;
683   PolynomialDataMem content;
684 
685   /* Check if cvode_mem exists */
686   if (cvode_mem == NULL) {
687     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeGetAdjDataPointPolynomial", MSGCV_NO_MEM);
688     return(CV_MEM_NULL);
689   }
690   cv_mem = (CVodeMem) cvode_mem;
691 
692   /* Was ASA initialized? */
693   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
694     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeGetAdjDataPointPolynomial", MSGCV_NO_ADJ);
695     return(CV_NO_ADJ);
696   }
697   ca_mem = cv_mem->cv_adj_mem;
698 
699   dt_mem = ca_mem->dt_mem;
700 
701   if (ca_mem->ca_IMtype != CV_POLYNOMIAL) {
702     cvProcessError(cv_mem, CV_ILL_INPUT, "CVODEA", "CVadjGetDataPointPolynomial", MSGCV_WRONG_INTERP);
703     return(CV_ILL_INPUT);
704   }
705 
706   *t = dt_mem[which]->t;
707 
708   content = (PolynomialDataMem) (dt_mem[which]->content);
709 
710   if (y != NULL)
711     N_VScale(ONE, content->y, y);
712 
713   *order = content->order;
714 
715   return(CV_SUCCESS);
716 }
717 
718 
719 /*
720  * CVodeGetAdjCurrentCheckPoint
721  *
722  * Returns the address of the 'active' check point.
723  */
724 
CVodeGetAdjCurrentCheckPoint(void * cvode_mem,void ** addr)725 int CVodeGetAdjCurrentCheckPoint(void *cvode_mem, void **addr)
726 {
727   CVodeMem cv_mem;
728   CVadjMem ca_mem;
729 
730   /* Check if cvode_mem exists */
731   if (cvode_mem == NULL) {
732     cvProcessError(NULL, CV_MEM_NULL, "CVODEA", "CVodeGetAdjCurrentCheckPoint", MSGCV_NO_MEM);
733     return(CV_MEM_NULL);
734   }
735   cv_mem = (CVodeMem) cvode_mem;
736 
737   /* Was ASA initialized? */
738   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
739     cvProcessError(cv_mem, CV_NO_ADJ, "CVODEA", "CVodeGetAdjCurrentCheckPoint", MSGCV_NO_ADJ);
740     return(CV_NO_ADJ);
741   }
742   ca_mem = cv_mem->cv_adj_mem;
743 
744   *addr = (void *) ca_mem->ca_ckpntData;
745 
746   return(CV_SUCCESS);
747 }
748