1 /*
2  * -----------------------------------------------------------------
3  * Programmer(s): Daniel R. Reynolds @ SMU
4  *                Radu Serban and Aaron Collier @ LLNL
5  * -----------------------------------------------------------------
6  * SUNDIALS Copyright Start
7  * Copyright (c) 2002-2021, Lawrence Livermore National Security
8  * and Southern Methodist University.
9  * All rights reserved.
10  *
11  * See the top-level LICENSE and NOTICE files for details.
12  *
13  * SPDX-License-Identifier: BSD-3-Clause
14  * SUNDIALS Copyright End
15  * -----------------------------------------------------------------
16  * This file contains implementations of the banded difference
17  * quotient Jacobian-based preconditioner and solver routines for
18  * use with the CVSLS linear solver interface.
19  * -----------------------------------------------------------------
20  */
21 
22 #include <stdio.h>
23 #include <stdlib.h>
24 
25 #include "cvodes_impl.h"
26 #include "cvodes_bandpre_impl.h"
27 #include "cvodes_ls_impl.h"
28 #include <sundials/sundials_math.h>
29 
30 #define MIN_INC_MULT RCONST(1000.0)
31 #define ZERO         RCONST(0.0)
32 #define ONE          RCONST(1.0)
33 #define TWO          RCONST(2.0)
34 
35 /* Prototypes of cvBandPrecSetup and cvBandPrecSolve */
36 static int cvBandPrecSetup(realtype t, N_Vector y, N_Vector fy,
37                            booleantype jok, booleantype *jcurPtr,
38                            realtype gamma, void *bp_data);
39 static int cvBandPrecSolve(realtype t, N_Vector y, N_Vector fy,
40                            N_Vector r, N_Vector z,
41                            realtype gamma, realtype delta,
42                            int lr, void *bp_data);
43 
44 /* Prototype for cvBandPrecFree */
45 static int cvBandPrecFree(CVodeMem cv_mem);
46 
47 /* Prototype for difference quotient Jacobian calculation routine */
48 static int cvBandPrecDQJac(CVBandPrecData pdata, realtype t, N_Vector y,
49                            N_Vector fy, N_Vector ftemp, N_Vector ytemp);
50 
51 
52 /*================================================================
53   PART I - Forward Problems
54   ================================================================*/
55 
56 /*-----------------------------------------------------------------
57   Initialization, Free, and Get Functions
58   NOTE: The band linear solver assumes a serial/OpenMP/Pthreads
59         implementation of the NVECTOR package. Therefore,
60         CVBandPrecInit will first test for a compatible N_Vector
61         internal representation by checking that the function
62         N_VGetArrayPointer exists.
63   -----------------------------------------------------------------*/
CVBandPrecInit(void * cvode_mem,sunindextype N,sunindextype mu,sunindextype ml)64 int CVBandPrecInit(void *cvode_mem, sunindextype N,
65                    sunindextype mu, sunindextype ml)
66 {
67   CVodeMem cv_mem;
68   CVLsMem cvls_mem;
69   CVBandPrecData pdata;
70   sunindextype mup, mlp, storagemu;
71   int flag;
72 
73   if (cvode_mem == NULL) {
74     cvProcessError(NULL, CVLS_MEM_NULL, "CVSBANDPRE",
75                    "CVBandPrecInit", MSGBP_MEM_NULL);
76     return(CVLS_MEM_NULL);
77   }
78   cv_mem = (CVodeMem) cvode_mem;
79 
80   /* Test if the CVSLS linear solver interface has been attached */
81   if (cv_mem->cv_lmem == NULL) {
82     cvProcessError(cv_mem, CVLS_LMEM_NULL, "CVSBANDPRE",
83                    "CVBandPrecInit", MSGBP_LMEM_NULL);
84     return(CVLS_LMEM_NULL);
85   }
86   cvls_mem = (CVLsMem) cv_mem->cv_lmem;
87 
88   /* Test compatibility of NVECTOR package with the BAND preconditioner */
89   if(cv_mem->cv_tempv->ops->nvgetarraypointer == NULL) {
90     cvProcessError(cv_mem, CVLS_ILL_INPUT, "CVSBANDPRE",
91                    "CVBandPrecInit", MSGBP_BAD_NVECTOR);
92     return(CVLS_ILL_INPUT);
93   }
94 
95   /* Allocate data memory */
96   pdata = NULL;
97   pdata = (CVBandPrecData) malloc(sizeof *pdata);
98   if (pdata == NULL) {
99     cvProcessError(cv_mem, CVLS_MEM_FAIL, "CVSBANDPRE",
100                    "CVBandPrecInit", MSGBP_MEM_FAIL);
101     return(CVLS_MEM_FAIL);
102   }
103 
104   /* Load pointers and bandwidths into pdata block. */
105   pdata->cvode_mem = cvode_mem;
106   pdata->N = N;
107   pdata->mu = mup = SUNMIN(N-1, SUNMAX(0,mu));
108   pdata->ml = mlp = SUNMIN(N-1, SUNMAX(0,ml));
109 
110   /* Initialize nfeBP counter */
111   pdata->nfeBP = 0;
112 
113   /* Allocate memory for saved banded Jacobian approximation. */
114   pdata->savedJ = NULL;
115   pdata->savedJ = SUNBandMatrixStorage(N, mup, mlp, mup);
116   if (pdata->savedJ == NULL) {
117     free(pdata); pdata = NULL;
118     cvProcessError(cv_mem, CVLS_MEM_FAIL, "CVSBANDPRE",
119                    "CVBandPrecInit", MSGBP_MEM_FAIL);
120     return(CVLS_MEM_FAIL);
121   }
122 
123   /* Allocate memory for banded preconditioner. */
124   storagemu = SUNMIN(N-1, mup+mlp);
125   pdata->savedP = NULL;
126   pdata->savedP = SUNBandMatrixStorage(N, mup, mlp, storagemu);
127   if (pdata->savedP == NULL) {
128     SUNMatDestroy(pdata->savedJ);
129     free(pdata); pdata = NULL;
130     cvProcessError(cv_mem, CVLS_MEM_FAIL, "CVSBANDPRE",
131                    "CVBandPrecInit", MSGBP_MEM_FAIL);
132     return(CVLS_MEM_FAIL);
133   }
134 
135   /* Allocate memory for banded linear solver */
136   pdata->LS = NULL;
137   pdata->LS = SUNLinSol_Band(cv_mem->cv_tempv, pdata->savedP);
138   if (pdata->LS == NULL) {
139     SUNMatDestroy(pdata->savedP);
140     SUNMatDestroy(pdata->savedJ);
141     free(pdata); pdata = NULL;
142     cvProcessError(cv_mem, CVLS_MEM_FAIL, "CVSBANDPRE",
143                    "CVBandPrecInit", MSGBP_MEM_FAIL);
144     return(CVLS_MEM_FAIL);
145   }
146 
147   /* allocate memory for temporary N_Vectors */
148   pdata->tmp1 = NULL;
149   pdata->tmp1 = N_VClone(cv_mem->cv_tempv);
150   if (pdata->tmp1 == NULL) {
151     SUNLinSolFree(pdata->LS);
152     SUNMatDestroy(pdata->savedP);
153     SUNMatDestroy(pdata->savedJ);
154     free(pdata); pdata = NULL;
155     cvProcessError(cv_mem, CVLS_MEM_FAIL, "CVSBANDPRE",
156                     "CVBandPrecInit", MSGBP_MEM_FAIL);
157     return(CVLS_MEM_FAIL);
158   }
159   pdata->tmp2 = NULL;
160   pdata->tmp2 = N_VClone(cv_mem->cv_tempv);
161   if (pdata->tmp2 == NULL) {
162     SUNLinSolFree(pdata->LS);
163     SUNMatDestroy(pdata->savedP);
164     SUNMatDestroy(pdata->savedJ);
165     N_VDestroy(pdata->tmp1);
166     free(pdata); pdata = NULL;
167     cvProcessError(cv_mem, CVLS_MEM_FAIL, "CVSBANDPRE",
168                     "CVBandPrecInit", MSGBP_MEM_FAIL);
169     return(CVLS_MEM_FAIL);
170   }
171 
172   /* initialize band linear solver object */
173   flag = SUNLinSolInitialize(pdata->LS);
174   if (flag != SUNLS_SUCCESS) {
175     SUNLinSolFree(pdata->LS);
176     SUNMatDestroy(pdata->savedP);
177     SUNMatDestroy(pdata->savedJ);
178     N_VDestroy(pdata->tmp1);
179     N_VDestroy(pdata->tmp2);
180     free(pdata); pdata = NULL;
181     cvProcessError(cv_mem, CVLS_SUNLS_FAIL, "CVSBANDPRE",
182                     "CVBandPrecInit", MSGBP_SUNLS_FAIL);
183     return(CVLS_SUNLS_FAIL);
184   }
185 
186   /* make sure P_data is free from any previous allocations */
187   if (cvls_mem->pfree)
188     cvls_mem->pfree(cv_mem);
189 
190   /* Point to the new P_data field in the LS memory */
191   cvls_mem->P_data = pdata;
192 
193   /* Attach the pfree function */
194   cvls_mem->pfree = cvBandPrecFree;
195 
196   /* Attach preconditioner solve and setup functions */
197   flag = CVodeSetPreconditioner(cvode_mem, cvBandPrecSetup,
198                                 cvBandPrecSolve);
199   return(flag);
200 }
201 
202 
CVBandPrecGetWorkSpace(void * cvode_mem,long int * lenrwBP,long int * leniwBP)203 int CVBandPrecGetWorkSpace(void *cvode_mem, long int *lenrwBP,
204                            long int *leniwBP)
205 {
206   CVodeMem cv_mem;
207   CVLsMem cvls_mem;
208   CVBandPrecData pdata;
209   sunindextype lrw1, liw1;
210   long int lrw, liw;
211   int flag;
212 
213   if (cvode_mem == NULL) {
214     cvProcessError(NULL, CVLS_MEM_NULL, "CVSBANDPRE",
215                    "CVBandPrecGetWorkSpace", MSGBP_MEM_NULL);
216     return(CVLS_MEM_NULL);
217   }
218   cv_mem = (CVodeMem) cvode_mem;
219 
220   if (cv_mem->cv_lmem == NULL) {
221     cvProcessError(cv_mem, CVLS_LMEM_NULL, "CVSBANDPRE",
222                    "CVBandPrecGetWorkSpace", MSGBP_LMEM_NULL);
223     return(CVLS_LMEM_NULL);
224   }
225   cvls_mem = (CVLsMem) cv_mem->cv_lmem;
226 
227   if (cvls_mem->P_data == NULL) {
228     cvProcessError(cv_mem, CVLS_PMEM_NULL, "CVSBANDPRE",
229                    "CVBandPrecGetWorkSpace", MSGBP_PMEM_NULL);
230     return(CVLS_PMEM_NULL);
231   }
232   pdata = (CVBandPrecData) cvls_mem->P_data;
233 
234   /* sum space requirements for all objects in pdata */
235   *leniwBP = 4;
236   *lenrwBP = 0;
237   if (cv_mem->cv_tempv->ops->nvspace) {
238     N_VSpace(cv_mem->cv_tempv, &lrw1, &liw1);
239     *leniwBP += 2*liw1;
240     *lenrwBP += 2*lrw1;
241   }
242   if (pdata->savedJ->ops->space) {
243     flag = SUNMatSpace(pdata->savedJ, &lrw, &liw);
244     if (flag != 0) return(-1);
245     *leniwBP += liw;
246     *lenrwBP += lrw;
247   }
248   if (pdata->savedP->ops->space) {
249     flag = SUNMatSpace(pdata->savedP, &lrw, &liw);
250     if (flag != 0) return(-1);
251     *leniwBP += liw;
252     *lenrwBP += lrw;
253   }
254   if (pdata->LS->ops->space) {
255     flag = SUNLinSolSpace(pdata->LS, &lrw, &liw);
256     if (flag != 0) return(-1);
257     *leniwBP += liw;
258     *lenrwBP += lrw;
259   }
260 
261   return(CVLS_SUCCESS);
262 }
263 
264 
CVBandPrecGetNumRhsEvals(void * cvode_mem,long int * nfevalsBP)265 int CVBandPrecGetNumRhsEvals(void *cvode_mem, long int *nfevalsBP)
266 {
267   CVodeMem cv_mem;
268   CVLsMem cvls_mem;
269   CVBandPrecData pdata;
270 
271   if (cvode_mem == NULL) {
272     cvProcessError(NULL, CVLS_MEM_NULL, "CVSBANDPRE",
273                    "CVBandPrecGetNumRhsEvals", MSGBP_MEM_NULL);
274     return(CVLS_MEM_NULL);
275   }
276   cv_mem = (CVodeMem) cvode_mem;
277 
278   if (cv_mem->cv_lmem == NULL) {
279     cvProcessError(cv_mem, CVLS_LMEM_NULL, "CVSBANDPRE",
280                    "CVBandPrecGetNumRhsEvals", MSGBP_LMEM_NULL);
281     return(CVLS_LMEM_NULL);
282   }
283   cvls_mem = (CVLsMem) cv_mem->cv_lmem;
284 
285   if (cvls_mem->P_data == NULL) {
286     cvProcessError(cv_mem, CVLS_PMEM_NULL, "CVSBANDPRE",
287                    "CVBandPrecGetNumRhsEvals", MSGBP_PMEM_NULL);
288     return(CVLS_PMEM_NULL);
289   }
290   pdata = (CVBandPrecData) cvls_mem->P_data;
291 
292   *nfevalsBP = pdata->nfeBP;
293 
294   return(CVLS_SUCCESS);
295 }
296 
297 
298 /*-----------------------------------------------------------------
299   cvBandPrecSetup
300   -----------------------------------------------------------------
301   Together cvBandPrecSetup and cvBandPrecSolve use a banded
302   difference quotient Jacobian to create a preconditioner.
303   cvBandPrecSetup calculates a new J, if necessary, then
304   calculates P = I - gamma*J, and does an LU factorization of P.
305 
306   The parameters of cvBandPrecSetup are as follows:
307 
308   t       is the current value of the independent variable.
309 
310   y       is the current value of the dependent variable vector,
311           namely the predicted value of y(t).
312 
313   fy      is the vector f(t,y).
314 
315   jok     is an input flag indicating whether Jacobian-related
316           data needs to be recomputed, as follows:
317             jok == SUNFALSE means recompute Jacobian-related data
318                    from scratch.
319             jok == SUNTRUE means that Jacobian data from the
320                    previous PrecSetup call will be reused
321                    (with the current value of gamma).
322           A cvBandPrecSetup call with jok == SUNTRUE should only
323           occur after a call with jok == SUNFALSE.
324 
325   *jcurPtr is a pointer to an output integer flag which is
326            set by cvBandPrecSetup as follows:
327              *jcurPtr = SUNTRUE if Jacobian data was recomputed.
328              *jcurPtr = SUNFALSE if Jacobian data was not recomputed,
329                         but saved data was reused.
330 
331   gamma   is the scalar appearing in the Newton matrix.
332 
333   bp_data is a pointer to preconditoner data (set by cvBandPrecInit)
334 
335   The value to be returned by the cvBandPrecSetup function is
336     0  if successful, or
337     1  if the band factorization failed.
338   -----------------------------------------------------------------*/
cvBandPrecSetup(realtype t,N_Vector y,N_Vector fy,booleantype jok,booleantype * jcurPtr,realtype gamma,void * bp_data)339 static int cvBandPrecSetup(realtype t, N_Vector y, N_Vector fy,
340                            booleantype jok, booleantype *jcurPtr,
341                            realtype gamma, void *bp_data)
342 {
343   CVBandPrecData pdata;
344   CVodeMem cv_mem;
345   int retval;
346 
347   /* Assume matrix and lpivots have already been allocated. */
348   pdata = (CVBandPrecData) bp_data;
349   cv_mem = (CVodeMem) pdata->cvode_mem;
350 
351   if (jok) {
352 
353     /* If jok = SUNTRUE, use saved copy of J. */
354     *jcurPtr = SUNFALSE;
355     retval = SUNMatCopy(pdata->savedJ, pdata->savedP);
356     if (retval < 0) {
357       cvProcessError(cv_mem, -1, "CVBANDPRE",
358                      "cvBandPrecSetup", MSGBP_SUNMAT_FAIL);
359       return(-1);
360     }
361     if (retval > 0) {
362       return(1);
363     }
364 
365   } else {
366 
367     /* If jok = SUNFALSE, call CVBandPDQJac for new J value. */
368     *jcurPtr = SUNTRUE;
369     retval = SUNMatZero(pdata->savedJ);
370     if (retval < 0) {
371       cvProcessError(cv_mem, -1, "CVBANDPRE",
372                      "cvBandPrecSetup", MSGBP_SUNMAT_FAIL);
373       return(-1);
374     }
375     if (retval > 0) {
376       return(1);
377     }
378 
379     retval = cvBandPrecDQJac(pdata, t, y, fy,
380                              pdata->tmp1, pdata->tmp2);
381     if (retval < 0) {
382       cvProcessError(cv_mem, -1, "CVBANDPRE",
383                      "cvBandPrecSetup", MSGBP_RHSFUNC_FAILED);
384       return(-1);
385     }
386     if (retval > 0) {
387       return(1);
388     }
389 
390     retval = SUNMatCopy(pdata->savedJ, pdata->savedP);
391     if (retval < 0) {
392       cvProcessError(cv_mem, -1, "CVBANDPRE",
393                      "cvBandPrecSetup", MSGBP_SUNMAT_FAIL);
394       return(-1);
395     }
396     if (retval > 0) {
397       return(1);
398     }
399 
400   }
401 
402   /* Scale and add identity to get savedP = I - gamma*J. */
403   retval = SUNMatScaleAddI(-gamma, pdata->savedP);
404   if (retval) {
405     cvProcessError(cv_mem, -1, "CVBANDPRE",
406                    "cvBandPrecSetup", MSGBP_SUNMAT_FAIL);
407     return(-1);
408   }
409 
410   /* Do LU factorization of matrix and return error flag */
411   retval = SUNLinSolSetup_Band(pdata->LS, pdata->savedP);
412   return(retval);
413 }
414 
415 
416 /*-----------------------------------------------------------------
417   cvBandPrecSolve
418   -----------------------------------------------------------------
419   cvBandPrecSolve solves a linear system P z = r, where P is the
420   matrix computed by cvBandPrecond.
421 
422   The parameters of cvBandPrecSolve used here are as follows:
423 
424   r       is the right-hand side vector of the linear system.
425 
426   bp_data is a pointer to preconditoner data (set by CVBandPrecInit)
427 
428   z       is the output vector computed by cvBandPrecSolve.
429 
430   The value returned by the cvBandPrecSolve function is always 0,
431   indicating success.
432   -----------------------------------------------------------------*/
cvBandPrecSolve(realtype t,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,realtype gamma,realtype delta,int lr,void * bp_data)433 static int cvBandPrecSolve(realtype t, N_Vector y, N_Vector fy,
434                            N_Vector r, N_Vector z, realtype gamma,
435                            realtype delta, int lr, void *bp_data)
436 {
437   CVBandPrecData pdata;
438   int retval;
439 
440   /* Assume matrix and lpivots have already been allocated. */
441   pdata = (CVBandPrecData) bp_data;
442 
443   /* Call banded solver object to do the work */
444   retval = SUNLinSolSolve(pdata->LS, pdata->savedP, z, r, ZERO);
445   return(retval);
446 }
447 
448 
cvBandPrecFree(CVodeMem cv_mem)449 static int cvBandPrecFree(CVodeMem cv_mem)
450 {
451   CVLsMem cvls_mem;
452   CVBandPrecData pdata;
453 
454   if (cv_mem->cv_lmem == NULL) return(0);
455   cvls_mem = (CVLsMem) cv_mem->cv_lmem;
456 
457   if (cvls_mem->P_data == NULL) return(0);
458   pdata = (CVBandPrecData) cvls_mem->P_data;
459 
460   SUNLinSolFree(pdata->LS);
461   SUNMatDestroy(pdata->savedP);
462   SUNMatDestroy(pdata->savedJ);
463   N_VDestroy(pdata->tmp1);
464   N_VDestroy(pdata->tmp2);
465 
466   free(pdata);
467   pdata = NULL;
468 
469   return(0);
470 }
471 
472 
473 /*-----------------------------------------------------------------
474   cvBandPrecDQJac
475   -----------------------------------------------------------------
476   This routine generates a banded difference quotient approximation
477   to the Jacobian of f(t,y). It assumes that a band SUNMatrix is
478   stored column-wise, and that elements within each column are
479   contiguous. This makes it possible to get the address of a column
480   of J via the accessor function SUNBandMatrix_Column() and to
481   write a simple for loop to set each of the elements of a column
482   in succession.
483   -----------------------------------------------------------------*/
cvBandPrecDQJac(CVBandPrecData pdata,realtype t,N_Vector y,N_Vector fy,N_Vector ftemp,N_Vector ytemp)484 static int cvBandPrecDQJac(CVBandPrecData pdata, realtype t, N_Vector y,
485                            N_Vector fy, N_Vector ftemp, N_Vector ytemp)
486 {
487   CVodeMem cv_mem;
488   realtype fnorm, minInc, inc, inc_inv, yj, srur, conj;
489   sunindextype group, i, j, width, ngroups, i1, i2;
490   realtype *col_j, *ewt_data, *fy_data, *ftemp_data;
491   realtype *y_data, *ytemp_data, *cns_data;
492   int retval;
493 
494   /* initialize cns_data to avoid compiler warning */
495   cns_data = NULL;
496 
497   cv_mem = (CVodeMem) pdata->cvode_mem;
498 
499   /* Obtain pointers to the data for ewt, fy, ftemp, y, ytemp. */
500   ewt_data   = N_VGetArrayPointer(cv_mem->cv_ewt);
501   fy_data    = N_VGetArrayPointer(fy);
502   ftemp_data = N_VGetArrayPointer(ftemp);
503   y_data     = N_VGetArrayPointer(y);
504   ytemp_data = N_VGetArrayPointer(ytemp);
505   if (cv_mem->cv_constraintsSet)
506     cns_data  = N_VGetArrayPointer(cv_mem->cv_constraints);
507 
508   /* Load ytemp with y = predicted y vector. */
509   N_VScale(ONE, y, ytemp);
510 
511   /* Set minimum increment based on uround and norm of f. */
512   srur = SUNRsqrt(cv_mem->cv_uround);
513   fnorm = N_VWrmsNorm(fy, cv_mem->cv_ewt);
514   minInc = (fnorm != ZERO) ?
515     (MIN_INC_MULT * SUNRabs(cv_mem->cv_h) * cv_mem->cv_uround * pdata->N * fnorm) : ONE;
516 
517   /* Set bandwidth and number of column groups for band differencing. */
518   width = pdata->ml + pdata->mu + 1;
519   ngroups = SUNMIN(width, pdata->N);
520 
521   for (group = 1; group <= ngroups; group++) {
522 
523     /* Increment all y_j in group. */
524     for(j = group-1; j < pdata->N; j += width) {
525       inc = SUNMAX(srur*SUNRabs(y_data[j]), minInc/ewt_data[j]);
526       yj = y_data[j];
527 
528       /* Adjust sign(inc) again if yj has an inequality constraint. */
529       if (cv_mem->cv_constraintsSet) {
530         conj = cns_data[j];
531         if (SUNRabs(conj) == ONE)      {if ((yj+inc)*conj < ZERO)  inc = -inc;}
532         else if (SUNRabs(conj) == TWO) {if ((yj+inc)*conj <= ZERO) inc = -inc;}
533       }
534 
535       ytemp_data[j] += inc;
536     }
537 
538     /* Evaluate f with incremented y. */
539     retval = cv_mem->cv_f(t, ytemp, ftemp, cv_mem->cv_user_data);
540     pdata->nfeBP++;
541     if (retval != 0) return(retval);
542 
543     /* Restore ytemp, then form and load difference quotients. */
544     for (j = group-1; j < pdata->N; j += width) {
545       yj = y_data[j];
546       ytemp_data[j] = y_data[j];
547       col_j = SUNBandMatrix_Column(pdata->savedJ,j);
548       inc = SUNMAX(srur*SUNRabs(y_data[j]), minInc/ewt_data[j]);
549 
550       /* Adjust sign(inc) as before. */
551       if (cv_mem->cv_constraintsSet) {
552         conj = cns_data[j];
553         if (SUNRabs(conj) == ONE)      {if ((yj+inc)*conj < ZERO)  inc = -inc;}
554         else if (SUNRabs(conj) == TWO) {if ((yj+inc)*conj <= ZERO) inc = -inc;}
555       }
556 
557       inc_inv = ONE/inc;
558       i1 = SUNMAX(0, j-pdata->mu);
559       i2 = SUNMIN(j + pdata->ml, pdata->N - 1);
560       for (i=i1; i <= i2; i++)
561         SM_COLUMN_ELEMENT_B(col_j,i,j) = inc_inv * (ftemp_data[i] - fy_data[i]);
562     }
563   }
564 
565   return(0);
566 }
567 
568 
569 /*================================================================
570   PART II - Backward Problems
571   ================================================================*/
572 
573 /*---------------------------------------------------------------
574   User-Callable initialization function: wrapper for the backward
575   phase around the corresponding CVODES functions
576   ---------------------------------------------------------------*/
CVBandPrecInitB(void * cvode_mem,int which,sunindextype nB,sunindextype muB,sunindextype mlB)577 int CVBandPrecInitB(void *cvode_mem, int which, sunindextype nB,
578                     sunindextype muB, sunindextype mlB)
579 {
580   CVodeMem cv_mem;
581   CVadjMem ca_mem;
582   CVodeBMem cvB_mem;
583   void *cvodeB_mem;
584   int flag;
585 
586   /* Check if cvode_mem exists */
587   if (cvode_mem == NULL) {
588     cvProcessError(NULL, CVLS_MEM_NULL, "CVSBANDPRE",
589                    "CVBandPrecInitB", MSGBP_MEM_NULL);
590     return(CVLS_MEM_NULL);
591   }
592   cv_mem = (CVodeMem) cvode_mem;
593 
594   /* Was ASA initialized? */
595   if (cv_mem->cv_adjMallocDone == SUNFALSE) {
596     cvProcessError(cv_mem, CVLS_NO_ADJ, "CVSBANDPRE",
597                    "CVBandPrecInitB", MSGBP_NO_ADJ);
598     return(CVLS_NO_ADJ);
599   }
600   ca_mem = cv_mem->cv_adj_mem;
601 
602   /* Check which */
603   if ( which >= ca_mem->ca_nbckpbs ) {
604     cvProcessError(cv_mem, CVLS_ILL_INPUT, "CVSBANDPRE",
605                    "CVBandPrecInitB", MSGBP_BAD_WHICH);
606     return(CVLS_ILL_INPUT);
607   }
608 
609   /* Find the CVodeBMem entry in the linked list corresponding to which */
610   cvB_mem = ca_mem->cvB_mem;
611   while (cvB_mem != NULL) {
612     if ( which == cvB_mem->cv_index ) break;
613     /* advance */
614     cvB_mem = cvB_mem->cv_next;
615   }
616   /* cv_mem corresponding to 'which' problem. */
617   cvodeB_mem = (void *) (cvB_mem->cv_mem);
618 
619   /* Set pfree */
620   cvB_mem->cv_pfree = NULL;
621 
622   /* Initialize the band preconditioner for this backward problem. */
623   flag = CVBandPrecInit(cvodeB_mem, nB, muB, mlB);
624   return(flag);
625 }
626