Show cvodes.c syntax highlighted
/*
* -----------------------------------------------------------------
* $Revision: 1.3 $
* $Date: 2006/02/16 21:49:19 $
* -----------------------------------------------------------------
* Programmer(s): Alan C. Hindmarsh and Radu Serban @ LLNL
* -----------------------------------------------------------------
* Copyright (c) 2002, The Regents of the University of California.
* Produced at the Lawrence Livermore National Laboratory.
* All rights reserved.
* For details, see sundials/cvodes/LICENSE.
* -----------------------------------------------------------------
* This is the implementation file for the main CVODES integrator
* with sensitivity analysis capabilities.
* It is independent of the CVODES linear solver in use.
* -----------------------------------------------------------------
*/
/*=================================================================*/
/* Import Header Files */
/*=================================================================*/
#include <stdio.h>
#include <stdlib.h>
#include "Oscill8External.h" // emery added for oscill8
#include "cvodes_impl.h"
#include "sundialsmath.h"
#include "sundialstypes.h"
/*=================================================================*/
/* Macros */
/*=================================================================*/
/* Macro: loop */
#define loop for(;;)
/*=================================================================*/
/* CVODES Private Constants */
/*=================================================================*/
#define ZERO RCONST(0.0) /* real 0.0 */
#define TINY RCONST(1.0e-10) /* small number */
#define TENTH RCONST(0.1) /* real 0.1 */
#define FOURTH RCONST(0.25) /* real 0.25 */
#define HALF RCONST(0.5) /* real 0.5 */
#define ONE RCONST(1.0) /* real 1.0 */
#define TWO RCONST(2.0) /* real 2.0 */
#define THREE RCONST(3.0) /* real 3.0 */
#define FOUR RCONST(4.0) /* real 4.0 */
#define FIVE RCONST(5.0) /* real 5.0 */
#define TWELVE RCONST(12.0) /* real 12.0 */
#define HUN RCONST(100.0) /* real 100.0 */
/*=================================================================*/
/* CVODES Default Constants */
/*=================================================================*/
#define HMIN_DEFAULT ZERO /* hmin default value */
#define HMAX_INV_DEFAULT ZERO /* hmax_inv default value */
#define MXHNIL_DEFAULT 10 /* mxhnil default value */
#define MXSTEP_DEFAULT 500 /* mxstep default value */
/*=================================================================*/
/* CVODES Routine-Specific Constants */
/*=================================================================*/
/* CVodeGetDky and CVStep */
#define FUZZ_FACTOR RCONST(100.0)
/* CVHin */
#define HLB_FACTOR RCONST(100.0)
#define HUB_FACTOR RCONST(0.1)
#define H_BIAS HALF
#define MAX_ITERS 4
/* CVSet */
#define CORTES RCONST(0.1)
/* CVStep return values */
#define SUCCESS_STEP 0
#define REP_ERR_FAIL -1
#define REP_CONV_FAIL -2
#define SETUP_FAILED -3
#define SOLVE_FAILED -4
/* CVStep control constants */
#define PREDICT_AGAIN -5
#define DO_ERROR_TEST 1
/* CVStep */
#define THRESH RCONST(1.5)
#define ETAMX1 RCONST(10000.0)
#define ETAMX2 RCONST(10.0)
#define ETAMX3 RCONST(10.0)
#define ETAMXF RCONST(0.2)
#define ETAMIN RCONST(0.1)
#define ETACF RCONST(0.25)
#define ADDON RCONST(0.000001)
#define BIAS1 RCONST(6.0)
#define BIAS2 RCONST(6.0)
#define BIAS3 RCONST(10.0)
#define ONEPSM RCONST(1.000001)
#define SMALL_NST 10 /* nst > SMALL_NST => use ETAMX3 */
#define MXNCF 10 /* max no. of convergence failures during
one step try */
#define MXNEF 7 /* max no. of error test failures during
one step try */
#define MXNEF1 3 /* max no. of error test failures before
forcing a reduction of order */
#define SMALL_NEF 2 /* if an error failure occurs and
SMALL_NEF <= nef <= MXNEF1, then
reset eta = MIN(eta, ETAMXF) */
#define LONG_WAIT 10 /* number of steps to wait before
considering an order change when
q==1 and MXNEF1 error test failures
have occurred */
/* CVNls return values */
#define SOLVED 0
#define CONV_FAIL -1
#define SETUP_FAIL_UNREC -2
#define SOLVE_FAIL_UNREC -3
/* CVNls input flags */
#define FIRST_CALL 0
#define PREV_CONV_FAIL -1
#define PREV_ERR_FAIL -2
/* CVNls other constants */
#define NLS_MAXCOR 3 /* maximum no. of corrector iterations for the
nonlinear solver */
#define CRDOWN RCONST(0.3) /* constant used in the estimation of the
convergence rate (crate) of the
iterates for the nonlinear equation */
#define DGMAX RCONST(0.3) /* iter == CV_NEWTON, |gamma/gammap-1| > DGMAX
=> call lsetup */
#define RDIV TWO /* declare divergence if ratio del/delp > RDIV */
#define MSBP 20 /* max no. of steps between lsetup calls */
#define TRY_AGAIN 99 /* control constant for CVNlsNewton - should be
distinct from CVNls return values */
/* CVRcheck* return values */
#define INITROOT -1
#define CLOSERT -2
#define RTFOUND 1
/* CVSensRhs1DQ finite difference methods */
#define CENTERED1 0
#define CENTERED2 1
#define FORWARD1 2
#define FORWARD2 3
/*=================================================================*/
/* Private Helper Functions Prototypes */
/*=================================================================*/
static booleantype CVCheckNvector(N_Vector tmpl);
static int CVInitialSetup(CVodeMem cv_mem);
static booleantype CVAllocVectors(CVodeMem cv_mem, N_Vector tmpl);
static void CVFreeVectors(CVodeMem cv_mem);
static booleantype CVEwtSet(CVodeMem cv_mem, N_Vector ycur);
static booleantype CVEwtSetSS(CVodeMem cv_mem, N_Vector ycur);
static booleantype CVEwtSetSV(CVodeMem cv_mem, N_Vector ycur);
static booleantype CVHin(CVodeMem cv_mem, realtype tout);
static realtype CVUpperBoundH0(CVodeMem cv_mem, realtype tdist);
static realtype CVYddNorm(CVodeMem cv_mem, realtype hg);
static int CVStep(CVodeMem cv_mem);
static int CVsldet(CVodeMem cv_mem);
static void CVAdjustParams(CVodeMem cv_mem);
static void CVAdjustOrder(CVodeMem cv_mem, int deltaq);
static void CVAdjustAdams(CVodeMem cv_mem, int deltaq);
static void CVAdjustBDF(CVodeMem cv_mem, int deltaq);
static void CVIncreaseBDF(CVodeMem cv_mem);
static void CVDecreaseBDF(CVodeMem cv_mem);
static void CVRescale(CVodeMem cv_mem);
static void CVPredict(CVodeMem cv_mem);
static void CVSet(CVodeMem cv_mem);
static void CVSetAdams(CVodeMem cv_mem);
static realtype CVAdamsStart(CVodeMem cv_mem, realtype m[]);
static void CVAdamsFinish(CVodeMem cv_mem, realtype m[], realtype M[], realtype hsum);
static realtype CVAltSum(int iend, realtype a[], int k);
static void CVSetBDF(CVodeMem cv_mem);
static void CVSetTqBDF(CVodeMem cv_mem, realtype hsum, realtype alpha0,
realtype alpha0_hat, realtype xi_inv, realtype xistar_inv);
static int CVNls(CVodeMem cv_mem, int nflag);
static int CVNlsFunctional(CVodeMem cv_mem);
static int CVNlsNewton(CVodeMem cv_mem, int nflag);
static int CVNewtonIteration(CVodeMem cv_mem);
static int CVHandleNFlag(CVodeMem cv_mem, int *nflagPtr, realtype saved_t,
int *ncfPtr, long int *ncfnPtr);
static void CVRestore(CVodeMem cv_mem, realtype saved_t);
static booleantype CVDoErrorTest(CVodeMem cv_mem, int *nflagPtr, int *kflagPtr,
realtype saved_t, int *nefPtr, realtype *dsmPtr);
static void CVCompleteStep(CVodeMem cv_mem);
static void CVPrepareNextStep(CVodeMem cv_mem, realtype dsm);
static void CVSetEta(CVodeMem cv_mem);
static realtype CVComputeEtaqm1(CVodeMem cv_mem);
static realtype CVComputeEtaqp1(CVodeMem cv_mem);
static void CVChooseEta(CVodeMem cv_mem);
static void CVBDFStab(CVodeMem cv_mem);
static int CVHandleFailure(CVodeMem cv_mem,int kflag);
/*----------------*/
static int CVRcheck1(CVodeMem cv_mem);
static int CVRcheck2(CVodeMem cv_mem);
static int CVRcheck3(CVodeMem cv_mem);
static int CVRootfind(CVodeMem cv_mem);
/*----------------*/
static booleantype CVQuadAllocVectors(CVodeMem cv_mem, N_Vector tmpl);
static booleantype CVQuadEwtSet(CVodeMem cv_mem, N_Vector qcur);
static booleantype CVQuadEwtSetSS(CVodeMem cv_mem, N_Vector qcur);
static booleantype CVQuadEwtSetSV(CVodeMem cv_mem, N_Vector qcur);
static void CVQuadFreeVectors(CVodeMem cv_mem);
/*----------------*/
static booleantype CVQuadDoErrorTest(CVodeMem cv_mem, int *nflagPtr,
int *kflagPtr, realtype saved_t,
int *nefQPtr, realtype *dsmQPtr);
/*----------------*/
static realtype CVQuadUpdateNorm(CVodeMem cv_mem, realtype old_nrm,
N_Vector xQ, N_Vector wQ);
static realtype CVQuadUpdateDsm(CVodeMem cv_mem, realtype old_dsm,
realtype dsmQ);
/*----------------*/
static int CVSensTestTolerances(CVodeMem cv_mem);
static int CVSensSetTolerances(CVodeMem cv_mem);
static booleantype CVSensAllocAtol(CVodeMem cv_mem, void **atolSPtr);
static void CVSensFreeAtol(CVodeMem cv_mem, void *atolS);
static booleantype CVSensSetAtolSS(CVodeMem cv_mem, realtype *atolS);
static booleantype CVSensSetAtolSV(CVodeMem cv_mem, N_Vector *atolS);
/*----------------*/
static booleantype CVSensAllocVectors(CVodeMem cv_mem, N_Vector tmpl);
static void CVSensFreeVectors(CVodeMem cv_mem);
/*----------------*/
static booleantype CVSensEwtSet(CVodeMem cv_mem, N_Vector *yScur);
static booleantype CVSensEwtSetSS(CVodeMem cv_mem, N_Vector *yScur);
static booleantype CVSensEwtSetSV(CVodeMem cv_mem, N_Vector *yScur);
/*----------------*/
static int CVStgrNls(CVodeMem cv_mem);
static int CVStgrNlsFunctional(CVodeMem cv_mem);
static int CVStgrNlsNewton(CVodeMem cv_mem);
static int CVStgrNewtonIteration(CVodeMem cv_mem);
static int CVStgr1Nls(CVodeMem cv_mem, int is);
static int CVStgr1NlsFunctional(CVodeMem cv_mem, int is);
static int CVStgr1NlsNewton(CVodeMem cv_mem, int is);
static int CVStgr1NewtonIteration(CVodeMem cv_mem, int is);
static booleantype CVStgrDoErrorTest(CVodeMem cv_mem, int *nflagPtr, int *kflagPtr,
realtype saved_t, int *nefSPtr, realtype *dsmSPtr);
/*----------------*/
static realtype CVSensNorm(CVodeMem cv_mem, N_Vector *xS, N_Vector *wS);
static realtype CVSensUpdateNorm(CVodeMem cv_mem, realtype old_nrm,
N_Vector *xS, N_Vector *wS);
static realtype CVStgrUpdateDsm(CVodeMem cv_mem, realtype old_dsm,
realtype dsmS);
/*----------------*/
static void CVSensRhs(CVodeMem cv_mem, realtype time,
N_Vector ycur, N_Vector fcur,
N_Vector *yScur, N_Vector *fScur,
N_Vector temp1, N_Vector temp2);
static void CVSensRhs1(CVodeMem cv_mem, realtype time,
N_Vector ycur, N_Vector fcur,
int is, N_Vector yScur, N_Vector fScur,
N_Vector temp1, N_Vector temp2);
/*=================================================================*/
/* EXPORTED FUNCTIONS IMPLEMENTATION */
/*=================================================================*/
/*
* CVodeCreate
*
* CVodeCreate creates an internal memory block for a problem to
* be solved by CVODES.
* If successful, CVodeCreate returns a pointer to the problem memory.
* This pointer should be passed to CVodeMalloc.
* If an initialization error occurs, CVodeCreate prints an error
* message to standard err and returns NULL.
*/
void *CVodeCreate(int lmm, int iter)
{
int maxord;
CVodeMem cv_mem;
/* Test inputs */
if ((lmm != CV_ADAMS) && (lmm != CV_BDF)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_BAD_LMM);
return (NULL);
}
if ((iter != CV_FUNCTIONAL) && (iter != CV_NEWTON)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_BAD_ITER);
return (NULL);
}
cv_mem = (CVodeMem) malloc(sizeof(struct CVodeMemRec));
if (cv_mem == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_CVMEM_FAIL);
return (NULL);
}
maxord = (lmm == CV_ADAMS) ? ADAMS_Q_MAX : BDF_Q_MAX;
/* copy input parameters into cv_mem */
cv_mem->cv_lmm = lmm;
cv_mem->cv_iter = iter;
/* Set uround */
cv_mem->cv_uround = UNIT_ROUNDOFF;
/* Set default values for integrator optional inputs */
cv_mem->cv_f = NULL;
cv_mem->cv_f_data = NULL;
cv_mem->cv_errfp = stderr;
cv_mem->cv_qmax = maxord;
cv_mem->cv_mxstep = MXSTEP_DEFAULT;
cv_mem->cv_mxhnil = MXHNIL_DEFAULT;
cv_mem->cv_sldeton = FALSE;
cv_mem->cv_hin = ZERO;
cv_mem->cv_hmin = HMIN_DEFAULT;
cv_mem->cv_hmax_inv = HMAX_INV_DEFAULT;
cv_mem->cv_tstopset = FALSE;
cv_mem->cv_maxcor = NLS_MAXCOR;
cv_mem->cv_maxnef = MXNEF;
cv_mem->cv_maxncf = MXNCF;
cv_mem->cv_nlscoef = CORTES;
/* Set default values for quad. optional inputs */
cv_mem->cv_quad = FALSE;
cv_mem->cv_fQ = NULL;
cv_mem->cv_fQ_data = NULL;
cv_mem->cv_errconQ = FALSE;
cv_mem->cv_reltolQ = NULL;
cv_mem->cv_abstolQ = NULL;
/* Set default values for sensi. optional inputs */
cv_mem->cv_sensi = FALSE;
cv_mem->cv_fS_data = (void *)cv_mem;
cv_mem->cv_fS = CVSensRhsDQ;
cv_mem->cv_fS1 = CVSensRhs1DQ;
cv_mem->cv_fSDQ = TRUE;
cv_mem->cv_ifS = CV_ONESENS;
cv_mem->cv_rhomax = ZERO;
cv_mem->cv_pbar = NULL;
cv_mem->cv_plist = NULL;
cv_mem->cv_errconS = FALSE;
cv_mem->cv_maxcorS = NLS_MAXCOR;
cv_mem->cv_ncfS1 = NULL;
cv_mem->cv_ncfnS1 = NULL;
cv_mem->cv_nniS1 = NULL;
/* By default, CVODES sets sensi tolerances */
cv_mem->cv_setSensTol = TRUE;
cv_mem->cv_atolSallocated = FALSE;
cv_mem->cv_testSensTol = FALSE;
cv_mem->cv_reltolS = NULL;
cv_mem->cv_abstolS = NULL;
/* No mallocs have been done yet */
cv_mem->cv_MallocDone = FALSE;
cv_mem->cv_quadMallocDone = FALSE;
cv_mem->cv_sensMallocDone = FALSE;
/* Return pointer to CVODES memory block */
return((void *)cv_mem);
}
/*-----------------------------------------------------------------*/
#define iter (cv_mem->cv_iter)
#define lmm (cv_mem->cv_lmm)
#define errfp (cv_mem->cv_errfp)
/*-----------------------------------------------------------------*/
/*
* CVodeMalloc
*
* CVodeMalloc allocates and initializes memory for a problem. All
* problem inputs are checked for errors. If any error occurs during
* initialization, it is reported to the file whose file pointer is
* errfp and an error flag is returned. Otherwise, it returns CV_SUCCESS
*/
int CVodeMalloc(void *cvode_mem, CVRhsFn f, realtype t0, N_Vector y0,
int itol, realtype *reltol, void *abstol)
{
CVodeMem cv_mem;
booleantype nvectorOK, allocOK, neg_abstol;
long int lrw1, liw1;
int i,k;
/* Check cvode_mem */
if (cvode_mem==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_CVM_NO_MEM);
return(CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
/* Check for legal input parameters */
if (y0==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_Y0_NULL);
return(CV_ILL_INPUT);
}
if ((itol != CV_SS) && (itol != CV_SV)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ITOL);
return(CV_ILL_INPUT);
}
if (f == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_F_NULL);
return(CV_ILL_INPUT);
}
if (reltol == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_RELTOL_NULL);
return(CV_ILL_INPUT);
}
if (*reltol < ZERO) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_RELTOL);
return(CV_ILL_INPUT);
}
if (abstol == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_ABSTOL_NULL);
return(CV_ILL_INPUT);
}
/* Test if all required vector operations are implemented */
nvectorOK = CVCheckNvector(y0);
if(!nvectorOK) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_NVECTOR);
return(CV_ILL_INPUT);
}
/* Test absolute tolerances */
if (itol == CV_SS) {
neg_abstol = (*((realtype *)abstol) < ZERO);
} else {
neg_abstol = (N_VMin((N_Vector)abstol) < ZERO);
}
if (neg_abstol) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ABSTOL);
return(CV_ILL_INPUT);
}
/* Set space requirements for one N_Vector */
if (y0->ops->nvspace != NULL) {
N_VSpace(y0, &lrw1, &liw1);
} else {
lrw1 = 0;
liw1 = 0;
}
cv_mem->cv_lrw1 = lrw1;
cv_mem->cv_liw1 = liw1;
/* Allocate the vectors (using y0 as a template) */
allocOK = CVAllocVectors(cv_mem, y0);
if (!allocOK) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_MEM_FAIL);
return(CV_MEM_FAIL);
}
/* Copy tolerances into memory */
cv_mem->cv_itol = itol;
cv_mem->cv_reltol = reltol;
cv_mem->cv_abstol = abstol;
/* All error checking is complete at this point */
/* Copy the input parameters into CVODES state */
cv_mem->cv_f = f;
cv_mem->cv_tn = t0;
/* Set step parameters */
cv_mem->cv_q = 1;
cv_mem->cv_L = 2;
cv_mem->cv_qwait = cv_mem->cv_L;
cv_mem->cv_etamax = ETAMX1;
cv_mem->cv_qu = 0;
cv_mem->cv_hu = ZERO;
cv_mem->cv_tolsf = ONE;
/* Set the linear solver addresses to NULL.
(We check != NULL later, in CVode, if using CV_NEWTON.) */
cv_mem->cv_linit = NULL;
cv_mem->cv_lsetup = NULL;
cv_mem->cv_lsolve = NULL;
cv_mem->cv_lfree = NULL;
cv_mem->cv_lmem = NULL;
/* Set forceSetup to FALSE */
cv_mem->cv_forceSetup = FALSE;
/* Initialize zn[0] in the history array */
N_VScale(ONE, y0, cv_mem->cv_zn[0]);
/* Initialize all the counters */
cv_mem->cv_nst = 0;
cv_mem->cv_nfe = 0;
cv_mem->cv_ncfn = 0;
cv_mem->cv_netf = 0;
cv_mem->cv_nni = 0;
cv_mem->cv_nsetups = 0;
cv_mem->cv_nhnil = 0;
cv_mem->cv_nstlp = 0;
cv_mem->cv_nscon = 0;
cv_mem->cv_nge = 0;
/* Initialize root finding variables */
cv_mem->cv_glo = NULL;
cv_mem->cv_ghi = NULL;
cv_mem->cv_groot = NULL;
cv_mem->cv_iroots = NULL;
cv_mem->cv_gfun = NULL;
cv_mem->cv_g_data = NULL;
cv_mem->cv_nrtfn = 0;
/* Initialize Stablilty Limit Detection data */
/* NOTE: We do this even if stab lim det was not
turned on yet. This way, the user can turn it
on at any time */
cv_mem->cv_nor = 0;
for (i = 1; i <= 5; i++)
for (k = 1; k <= 3; k++)
cv_mem->cv_ssdat[i-1][k-1] = ZERO;
/* Problem has been successfully initialized */
cv_mem->cv_MallocDone = TRUE;
return(CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
/*
* CVodeReInit
*
* CVodeReInit re-initializes CVODES' memory for a problem, assuming
* it has already been allocated in a prior CVodeMalloc call.
* All problem specification inputs are checked for errors.
* If any error occurs during initialization, it is reported to the
* file whose file pointer is errfp.
* The return value is CV_SUCCESS = 0 if no errors occurred, or
* a negative value otherwise.
*/
int CVodeReInit(void *cvode_mem, CVRhsFn f, realtype t0, N_Vector y0,
int itol, realtype *reltol, void *abstol)
{
CVodeMem cv_mem;
booleantype neg_abstol;
int i,k;
/* Check cvode_mem */
if (cvode_mem==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_CVM_NO_MEM);
return(CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
/* Check if cvode_mem was allocated */
if (cv_mem->cv_MallocDone == FALSE) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_CVREI_NO_MALLOC);
return(CV_NO_MALLOC);
}
/* Check for legal input parameters */
if (y0 == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_Y0_NULL);
return (CV_ILL_INPUT);
}
if ((itol != CV_SS) && (itol != CV_SV)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ITOL);
return (CV_ILL_INPUT);
}
if (f == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_F_NULL);
return (CV_ILL_INPUT);
}
if (reltol == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_RELTOL_NULL);
return (CV_ILL_INPUT);
}
if (*reltol < ZERO) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_RELTOL);
return (CV_ILL_INPUT);
}
if (abstol == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_ABSTOL_NULL);
return (CV_ILL_INPUT);
}
if (itol == CV_SS) {
neg_abstol = (*((realtype *)abstol) < ZERO);
} else {
neg_abstol = (N_VMin((N_Vector)abstol) < ZERO);
}
if (neg_abstol) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ABSTOL);
return (CV_ILL_INPUT);
}
/* Copy tolerances into memory and set the ewt vector */
cv_mem->cv_itol = itol;
cv_mem->cv_reltol = reltol;
cv_mem->cv_abstol = abstol;
/* All error checking is complete at this point */
/* Copy the input parameters into CVODE state */
cv_mem->cv_f = f;
cv_mem->cv_tn = t0;
/* Set step parameters */
cv_mem->cv_q = 1;
cv_mem->cv_L = 2;
cv_mem->cv_qwait = cv_mem->cv_L;
cv_mem->cv_etamax = ETAMX1;
cv_mem->cv_qu = 0;
cv_mem->cv_hu = ZERO;
cv_mem->cv_tolsf = ONE;
/* Set forceSetup to FALSE */
cv_mem->cv_forceSetup = FALSE;
/* Initialize zn[0] in the history array */
N_VScale(ONE, y0, cv_mem->cv_zn[0]);
/* Initialize all the counters */
cv_mem->cv_nst = 0;
cv_mem->cv_nfe = 0;
cv_mem->cv_ncfn = 0;
cv_mem->cv_netf = 0;
cv_mem->cv_nni = 0;
cv_mem->cv_nsetups = 0;
cv_mem->cv_nhnil = 0;
cv_mem->cv_nstlp = 0;
cv_mem->cv_nscon = 0;
cv_mem->cv_nge = 0;
/* Initialize Stablilty Limit Detection data */
cv_mem->cv_nor = 0;
for (i = 1; i <= 5; i++)
for (k = 1; k <= 3; k++)
cv_mem->cv_ssdat[i-1][k-1] = ZERO;
/* Problem has been successfully re-initialized */
return(CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
#define gfun (cv_mem->cv_gfun)
#define glo (cv_mem->cv_glo)
#define ghi (cv_mem->cv_ghi)
#define groot (cv_mem->cv_groot)
#define iroots (cv_mem->cv_iroots)
/*-----------------------------------------------------------------*/
/*
* CVodeRootInit
*
* CVodeRootInit initializes a rootfinding problem to be solved
* during the integration of the ODE system. It loads the root
* function pointer and the number of root functions, and allocates
* workspace memory. The return value is CV_SUCCESS = 0 if no errors
* occurred, or a negative value otherwise.
*/
int CVodeRootInit(void *cvode_mem, CVRootFn g, int nrtfn)
{
CVodeMem cv_mem;
int nrt;
/* Check cvode_mem */
if (cvode_mem==NULL) {
EXTERNAL_LIBRARY_PRINTF(MSGCVS_ROOT_NO_MEM);
return(CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
nrt = (nrtfn < 0) ? 0 : nrtfn;
/* If rerunning CVodeRootInit() with a different number of root
functions (changing number of gfun components), then free
currently held memory resources */
if ((nrt != cv_mem->cv_nrtfn) && (cv_mem->cv_nrtfn > 0)) {
free(glo);
free(ghi);
free(groot);
free(iroots);
/* Linux version of free() routine doesn't set pointer to NULL */
glo = ghi = groot = NULL;
iroots = NULL;
}
/* If CVodeRootInit() was called with nrtfn == 0, then set cv_nrtfn to
zero and cv_gfun to NULL before returning */
if (nrt == 0) {
cv_mem->cv_nrtfn = nrt;
gfun = NULL;
return(CV_SUCCESS);
}
/* If rerunning CVodeRootInit() with the same number of root functions
(not changing number of gfun components), then check if the root
function argument has changed */
/* If g != NULL then return as currently reserved memory resources
will suffice */
if (nrt == cv_mem->cv_nrtfn) {
if (g != gfun) {
if (g == NULL) {
free(glo);
free(ghi);
free(groot);
free(iroots);
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_ROOT_FUNC_NULL);
return(CV_RTFUNC_NULL);
}
else {
gfun = g;
return(CV_SUCCESS);
}
}
else return(CV_SUCCESS);
}
/* Set variable values in CVode memory block */
cv_mem->cv_nrtfn = nrt;
if (g == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_ROOT_FUNC_NULL);
return(CV_RTFUNC_NULL);
}
else gfun = g;
/* Allocate necessary memory and return */
glo = (realtype *) malloc(nrt*sizeof(realtype));
if (glo == NULL) {
EXTERNAL_LIBRARY_PRINTF(MSGCVS_ROOT_MEM_FAIL);
return(CV_MEM_FAIL);
}
ghi = (realtype *) malloc(nrt*sizeof(realtype));
if (ghi == NULL) {
free(glo);
EXTERNAL_LIBRARY_PRINTF(MSGCVS_ROOT_MEM_FAIL);
return(CV_MEM_FAIL);
}
groot = (realtype *) malloc(nrt*sizeof(realtype));
if (groot == NULL) {
free(glo); free(ghi);
EXTERNAL_LIBRARY_PRINTF(MSGCVS_ROOT_MEM_FAIL);
return(CV_MEM_FAIL);
}
iroots = (int *) malloc(nrt*sizeof(int));
if (iroots == NULL) {
free(glo); free(ghi); free(groot);
EXTERNAL_LIBRARY_PRINTF(MSGCVS_ROOT_MEM_FAIL);
return(CV_MEM_FAIL);
}
return(CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
/*
* CVodeQuadMalloc
*
* CVodeQuadMalloc allocates and initializes quadrature related
* memory for a problem. All problem specification inputs are
* checked for errors. If any error occurs during initialization,
* it is reported to the file whose file pointer is errfp.
* The return value is CV_SUCCESS = 0 if no errors occurred, or
* a negative value otherwise.
*/
int CVodeQuadMalloc(void *cvode_mem, CVQuadRhsFn fQ, N_Vector yQ0)
{
CVodeMem cv_mem;
booleantype allocOK;
long int lrw1Q, liw1Q;
/* Check cvode_mem */
if (cvode_mem==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_QCVM_NO_MEM);
return(CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
/* Set space requirements for one N_Vector */
N_VSpace(yQ0, &lrw1Q, &liw1Q);
cv_mem->cv_lrw1Q = lrw1Q;
cv_mem->cv_liw1Q = liw1Q;
/* Allocate the vectors (using yQ0 as a template) */
allocOK = CVQuadAllocVectors(cv_mem, yQ0);
if (!allocOK) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_QCVM_MEM_FAIL);
return(CV_MEM_FAIL);
}
/* Initialize znQ[0] in the history array */
N_VScale(ONE, yQ0, cv_mem->cv_znQ[0]);
/* Copy the input parameters into CVODES state */
cv_mem->cv_fQ = fQ;
/* Initialize counters */
cv_mem->cv_nfQe = 0;
cv_mem->cv_netfQ = 0;
/* Quadrature integration turned ON */
cv_mem->cv_quad = TRUE;
cv_mem->cv_quadMallocDone = TRUE;
/* Quadrature initialization was successfull */
return(CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
/*
* CVodeQuadReInit
*
* CVodeQuadReInit re-initializes CVODES' quadrature related memory
* for a problem, assuming it has already been allocated in prior
* calls to CVodeMalloc and CvodeQuadMalloc.
* All problem specification inputs are checked for errors.
* If any error occurs during initialization, it is reported to the
* file whose file pointer is errfp.
* The return value is CV_SUCCESS = 0 if no errors occurred, or
* a negative value otherwise.
*/
int CVodeQuadReInit(void *cvode_mem, CVQuadRhsFn fQ, N_Vector yQ0)
{
CVodeMem cv_mem;
/* Check cvode_mem */
if (cvode_mem==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_QCVM_NO_MEM);
return(CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
/* Ckeck if quadrature was initialized? */
if (cv_mem->cv_quadMallocDone == FALSE) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_QREI_NO_QUAD);
return(CV_NO_QUAD);
}
/* Initialize znQ[0] in the history array */
N_VScale(ONE, yQ0, cv_mem->cv_znQ[0]);
/* Copy the input parameters into CVODE state */
cv_mem->cv_fQ = fQ;
/* Initialize counters */
cv_mem->cv_nfQe = 0;
cv_mem->cv_netfQ = 0;
/* Quadrature integration turned ON */
cv_mem->cv_quad = TRUE;
/* Quadrature re-initialization was successfull */
return(CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
#define ifS (cv_mem->cv_ifS)
#define fSDQ (cv_mem->cv_fSDQ)
#define stgr1alloc (cv_mem->cv_stgr1alloc)
#define nniS1 (cv_mem->cv_nniS1)
#define ncfnS1 (cv_mem->cv_ncfnS1)
#define ncfS1 (cv_mem->cv_ncfS1)
/*-----------------------------------------------------------------*/
/*
* CVodeSenMalloc
*
* CVodeSensMalloc allocates and initializes sensitivity related
* memory for a problem. All problem specification inputs are
* checked for errors. If any error occurs during initialization,
* it is reported to the file whose file pointer is errfp.
* The return value is CV_SUCCESS = 0 if no errors occurred, or
* a negative value otherwise.
*/
int CVodeSensMalloc(void *cvode_mem, int Ns, int ism,
realtype *p, int *plist, N_Vector *yS0)
{
CVodeMem cv_mem;
booleantype allocOK;
int is;
/* Check cvode_mem */
if (cvode_mem==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_SCVM_NO_MEM);
return(CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
/* Check if Ns is legal */
if (Ns<=0) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_NS);
return(CV_ILL_INPUT);
}
cv_mem->cv_Ns = Ns;
/* Check if ism is legal */
if ((ism!=CV_SIMULTANEOUS) && (ism!=CV_STAGGERED) && (ism!=CV_STAGGERED1)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ISM);
return(CV_ILL_INPUT);
}
cv_mem->cv_ism = ism;
/* Check if p is non-null */
if (p==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_P_NULL);
return(CV_ILL_INPUT);
}
cv_mem->cv_p = p;
cv_mem->cv_plist = plist;
/* Check if yS0 is non-null */
if ((cv_mem->cv_yS = yS0) == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_YS0_NULL);
return(CV_ILL_INPUT);
}
/* Allocate ncfS1, ncfnS1, and nniS1 if needed */
if (ism == CV_STAGGERED1) {
stgr1alloc = TRUE;
ncfS1 = (int*)malloc(Ns*sizeof(int));
ncfnS1 = (long int*)malloc(Ns*sizeof(long int));
nniS1 = (long int*)malloc(Ns*sizeof(long int));
if ( (ncfS1 == NULL) || (ncfnS1 == NULL) || (nniS1 == NULL) ) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SCVM_MEM_FAIL);
return(CV_MEM_FAIL);
}
} else {
stgr1alloc = FALSE;
}
/* Allocate the vectors (using yS0[0] as a template) */
allocOK = CVSensAllocVectors(cv_mem, yS0[0]);
if (!allocOK) {
if (stgr1alloc) {
free(ncfS1);
free(ncfnS1);
free(nniS1);
}
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SCVM_MEM_FAIL);
return(CV_MEM_FAIL);
}
/*----------------------------------------------
All error checking is complete at this point
-----------------------------------------------*/
/* Initialize znS[0] in the history array */
for (is=0; is<Ns; is++)
N_VScale(ONE, yS0[is], cv_mem->cv_znS[0][is]);
/* Initialize all sensitivity related counters */
cv_mem->cv_nfSe = 0;
cv_mem->cv_nfeS = 0;
cv_mem->cv_ncfnS = 0;
cv_mem->cv_netfS = 0;
cv_mem->cv_nniS = 0;
cv_mem->cv_nsetupsS = 0;
if (ism==CV_STAGGERED1)
for (is=0; is<Ns; is++) {
ncfnS1[is] = 0;
nniS1[is] = 0;
}
/* Sensitivities will be computed */
cv_mem->cv_sensi = TRUE;
cv_mem->cv_sensMallocDone = TRUE;
/* Sensitivity initialization was successfull */
return (CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
#define Ns (cv_mem->cv_Ns)
/*-----------------------------------------------------------------*/
/*
* CVodeSensReInit
*
* CVodeSensReInit re-initializes CVODES' sensitivity related memory
* for a problem, assuming it has already been allocated in prior
* calls to CVodeMalloc and CVodeSensMalloc.
* All problem specification inputs are checked for errors.
* The number of sensitivities Ns is assumed to be unchanged since
* the previous call to CVodeSensMalloc.
* If any error occurs during initialization, it is reported to the
* file whose file pointer is errfp.
* The return value is CV_SUCCESS = 0 if no errors occurred, or
* a negative value otherwise.
*/
int CVodeSensReInit(void *cvode_mem, int ism,
realtype *p, int *plist, N_Vector *yS0)
{
CVodeMem cv_mem;
int is;
/* Check cvode_mem */
if (cvode_mem==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_SCVM_NO_MEM);
return(CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
/* Was sensitivity initialized? */
if (cv_mem->cv_sensMallocDone == FALSE) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SREI_NO_SENSI);
return(CV_NO_SENS);
}
/* Check if ism is legal */
if ((ism!=CV_SIMULTANEOUS) && (ism!=CV_STAGGERED) && (ism!=CV_STAGGERED1)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ISM);
return(CV_ILL_INPUT);
}
cv_mem->cv_ism = ism;
/* Check if p is non-null */
if (p==NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_P_NULL);
return(CV_ILL_INPUT);
}
cv_mem->cv_p = p;
cv_mem->cv_plist = plist;
/* Check if yS0 is non-null */
if (yS0 == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_YS0_NULL);
return(CV_ILL_INPUT);
}
/* Allocate ncfS1, ncfnS1, and nniS1 if needed */
if ( (ism==CV_STAGGERED1) && (stgr1alloc==FALSE) ) {
stgr1alloc = TRUE;
ncfS1 = (int*)malloc(Ns*sizeof(int));
ncfnS1 = (long int*)malloc(Ns*sizeof(long int));
nniS1 = (long int*)malloc(Ns*sizeof(long int));
if ( (ncfS1==NULL) || (ncfnS1==NULL) || (nniS1==NULL) ) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SCVM_MEM_FAIL);
return(CV_MEM_FAIL);
}
}
/*----------------------------------------------
All error checking is complete at this point
-----------------------------------------------*/
/* Initialize znS[0] in the history array */
for (is=0; is<Ns; is++)
N_VScale(ONE, yS0[is], cv_mem->cv_znS[0][is]);
/* Initialize all sensitivity related counters */
cv_mem->cv_nfSe = 0;
cv_mem->cv_nfeS = 0;
cv_mem->cv_ncfnS = 0;
cv_mem->cv_netfS = 0;
cv_mem->cv_nniS = 0;
cv_mem->cv_nsetupsS = 0;
if (ism==CV_STAGGERED1)
for (is=0; is<Ns; is++) {
ncfnS1[is] = 0;
nniS1[is] = 0;
}
/* Problem has been successfully re-initialized */
cv_mem->cv_sensi = TRUE;
return (CV_SUCCESS);
}
/*=================================================================*/
/* Readibility Constants */
/*=================================================================*/
#define f (cv_mem->cv_f)
#define f_data (cv_mem->cv_f_data)
#define g_data (cv_mem->cv_g_data)
#define qmax (cv_mem->cv_qmax)
#define mxstep (cv_mem->cv_mxstep)
#define mxhnil (cv_mem->cv_mxhnil)
#define sldeton (cv_mem->cv_sldeton)
#define hin (cv_mem->cv_hin)
#define hmin (cv_mem->cv_hmin)
#define hmax_inv (cv_mem->cv_hmax_inv)
#define tstop (cv_mem->cv_tstop)
#define tstopset (cv_mem->cv_tstopset)
#define maxnef (cv_mem->cv_maxnef)
#define maxncf (cv_mem->cv_maxncf)
#define maxcor (cv_mem->cv_maxcor)
#define nlscoef (cv_mem->cv_nlscoef)
#define itol (cv_mem->cv_itol)
#define reltol (cv_mem->cv_reltol)
#define abstol (cv_mem->cv_abstol)
#define fQ (cv_mem->cv_fQ)
#define fQ_data (cv_mem->cv_fQ_data)
#define errconQ (cv_mem->cv_errconQ)
#define itolQ (cv_mem->cv_itolQ)
#define reltolQ (cv_mem->cv_reltolQ)
#define abstolQ (cv_mem->cv_abstolQ)
#define fS (cv_mem->cv_fS)
#define fS1 (cv_mem->cv_fS1)
#define fS_data (cv_mem->cv_fS_data)
#define rhomax (cv_mem->cv_rhomax)
#define pbar (cv_mem->cv_pbar)
#define errconS (cv_mem->cv_errconS)
#define maxcorS (cv_mem->cv_maxcorS)
#define itolS (cv_mem->cv_itolS)
#define reltolS (cv_mem->cv_reltolS)
#define abstolS (cv_mem->cv_abstolS)
#define ism (cv_mem->cv_ism)
#define p (cv_mem->cv_p)
#define plist (cv_mem->cv_plist)
#define uround (cv_mem->cv_uround)
#define zn (cv_mem->cv_zn)
#define ewt (cv_mem->cv_ewt)
#define y (cv_mem->cv_y)
#define acor (cv_mem->cv_acor)
#define tempv (cv_mem->cv_tempv)
#define ftemp (cv_mem->cv_ftemp)
#define q (cv_mem->cv_q)
#define qprime (cv_mem->cv_qprime)
#define next_q (cv_mem->cv_next_q)
#define qwait (cv_mem->cv_qwait)
#define L (cv_mem->cv_L)
#define h (cv_mem->cv_h)
#define hprime (cv_mem->cv_hprime)
#define next_h (cv_mem->cv_next_h)
#define eta (cv_mem->cv_eta)
#define etaqm1 (cv_mem->cv_etaqm1)
#define etaq (cv_mem->cv_etaq)
#define etaqp1 (cv_mem->cv_etaqp1)
#define nscon (cv_mem->cv_nscon)
#define hscale (cv_mem->cv_hscale)
#define tn (cv_mem->cv_tn)
#define tau (cv_mem->cv_tau)
#define tq (cv_mem->cv_tq)
#define l (cv_mem->cv_l)
#define rl1 (cv_mem->cv_rl1)
#define gamma (cv_mem->cv_gamma)
#define gammap (cv_mem->cv_gammap)
#define gamrat (cv_mem->cv_gamrat)
#define crate (cv_mem->cv_crate)
#define acnrm (cv_mem->cv_acnrm)
#define mnewt (cv_mem->cv_mnewt)
#define etamax (cv_mem->cv_etamax)
#define nst (cv_mem->cv_nst)
#define nfe (cv_mem->cv_nfe)
#define ncfn (cv_mem->cv_ncfn)
#define netf (cv_mem->cv_netf)
#define nni (cv_mem->cv_nni)
#define nsetups (cv_mem->cv_nsetups)
#define nhnil (cv_mem->cv_nhnil)
#define lrw1 (cv_mem->cv_lrw1)
#define liw1 (cv_mem->cv_liw1)
#define lrw (cv_mem->cv_lrw)
#define liw (cv_mem->cv_liw)
#define linit (cv_mem->cv_linit)
#define lsetup (cv_mem->cv_lsetup)
#define lsolve (cv_mem->cv_lsolve)
#define lfree (cv_mem->cv_lfree)
#define lmem (cv_mem->cv_lmem)
#define qu (cv_mem->cv_qu)
#define nstlp (cv_mem->cv_nstlp)
#define h0u (cv_mem->cv_h0u)
#define hu (cv_mem->cv_hu)
#define saved_tq5 (cv_mem->cv_saved_tq5)
#define jcur (cv_mem->cv_jcur)
#define tolsf (cv_mem->cv_tolsf)
#define setupNonNull (cv_mem->cv_setupNonNull)
#define forceSetup (cv_mem->cv_forceSetup)
#define nor (cv_mem->cv_nor)
#define ssdat (cv_mem->cv_ssdat)
#define nrtfn (cv_mem->cv_nrtfn)
#define tlo (cv_mem->cv_tlo)
#define thi (cv_mem->cv_thi)
#define tretlast (cv_mem->cv_tretlast)
#define toutc (cv_mem->cv_toutc)
#define troot (cv_mem->cv_troot)
#define ttol (cv_mem->cv_ttol)
#define taskc (cv_mem->cv_taskc)
#define irfnd (cv_mem->cv_irfnd)
#define nge (cv_mem->cv_nge)
#define sensi (cv_mem->cv_sensi)
#define znS (cv_mem->cv_znS)
#define ewtS (cv_mem->cv_ewtS)
#define acorS (cv_mem->cv_acorS)
#define yS (cv_mem->cv_yS)
#define tempvS (cv_mem->cv_tempvS)
#define ftempS (cv_mem->cv_ftempS)
#define crateS (cv_mem->cv_crateS)
#define acnrmS (cv_mem->cv_acnrmS)
#define nfSe (cv_mem->cv_nfSe)
#define nfeS (cv_mem->cv_nfeS)
#define nniS (cv_mem->cv_nniS)
#define ncfnS (cv_mem->cv_ncfnS)
#define netfS (cv_mem->cv_netfS)
#define nsetupsS (cv_mem->cv_nsetupsS)
#define stgr1alloc (cv_mem->cv_stgr1alloc)
#define setSensTol (cv_mem->cv_setSensTol)
#define testSensTol (cv_mem->cv_testSensTol)
#define atolSallocated (cv_mem->cv_atolSallocated)
#define quad (cv_mem->cv_quad)
#define znQ (cv_mem->cv_znQ)
#define ewtQ (cv_mem->cv_ewtQ)
#define acorQ (cv_mem->cv_acorQ)
#define yQ (cv_mem->cv_yQ)
#define tempvQ (cv_mem->cv_tempvQ)
#define acnrmQ (cv_mem->cv_acnrmQ)
#define nfQe (cv_mem->cv_nfQe)
#define netfQ (cv_mem->cv_netfQ)
#define lrw1Q (cv_mem->cv_lrw1Q)
#define liw1Q (cv_mem->cv_liw1Q)
/*-----------------------------------------------------------------*/
/*
* CVode
*
* This routine is the main driver of the CVODE package.
*
* It integrates over a time interval defined by the user, by calling
* CVStep to do internal time steps.
*
* The first time that CVode is called for a successfully initialized
* problem, it computes a tentative initial step size h.
*
* CVode supports four modes, specified by itask: CV_NORMAL, CV_ONE_STEP,
* CV_NORMAL_TSTOP, and CV_ONE_STEP_TSTOP.
* In the CV_NORMAL mode, the solver steps until it reaches or passes tout
* and then interpolates to obtain y(tout).
* In the CV_ONE_STEP mode, it takes one internal step and returns.
* CV_NORMAL_TSTOP and CV_ONE_STEP_TSTOP are similar to CV_NORMAL and CV_ONE_STEP,
* respectively, but the integration never proceeds past tstop (which
* must have been defined through a call to CVodeSetStopTime).
*/
int CVode(void *cvode_mem, realtype tout, N_Vector yout,
realtype *tret, int itask)
{
CVodeMem cv_mem;
N_Vector wrk1, wrk2;
long int nstloc;
int kflag, istate, ier, task, irfndp;
booleantype istop, hOK, ewtsetOK, ewtSsetOK, ewtQsetOK;
int is;
realtype troundoff, rh, nrm;
/* Check if cvode_mem exists */
if (cvode_mem == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_CVODE_NO_MEM);
return (CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
/* Check if cvode_mem was allocated */
if (cv_mem->cv_MallocDone == FALSE) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_CVODE_NO_MALLOC);
return(CV_NO_MALLOC);
}
/* Check for yout != NULL */
if ((y = yout) == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_YOUT_NULL);
return (CV_ILL_INPUT);
}
/* Check for tret != NULL */
if (tret == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_TRET_NULL);
return (CV_ILL_INPUT);
}
//@@edc BUG!!! this is set below!!!
//tretlast = *tret = tn;
/* Check for valid itask */
if ((itask != CV_NORMAL) &&
(itask != CV_ONE_STEP) &&
(itask != CV_NORMAL_TSTOP) &&
(itask != CV_ONE_STEP_TSTOP) ) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ITASK);
return (CV_ILL_INPUT);
}
/* Split itask into task and istop */
if ((itask == CV_NORMAL_TSTOP) || (itask == CV_ONE_STEP_TSTOP)) {
if ( tstopset == FALSE ) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_NO_TSTOP);
return(CV_ILL_INPUT);
}
istop = TRUE;
} else {
istop = FALSE;
}
if ((itask == CV_NORMAL) || (itask == CV_NORMAL_TSTOP)) {
task = CV_NORMAL; toutc = tout;
} else {
task = CV_ONE_STEP;
}
taskc = task;
/* If doing FSA, deal with the sensitivity tolerances */
if (sensi) {
if (testSensTol) {
ier = CVSensTestTolerances(cv_mem);
if (ier != CV_SUCCESS) return (ier);
testSensTol = FALSE;
}
if (setSensTol) {
ier = CVSensSetTolerances(cv_mem);
if (ier != CV_SUCCESS) return (ier);
setSensTol = FALSE;
}
}
/* Begin first call block */
if (nst == 0) {
/* Check inputs for corectness */
ier = CVInitialSetup(cv_mem);
if (ier!= CV_SUCCESS) return (ier);
/*
Call f at (t0,y0), set zn[1] = y'(t0),
set initial h (from H0 or CVHin), and scale zn[1] by h.
Also check for zeros of root function g at and near t0.
If computing sensitivities, call fS at (t0,y0,yS0), set
znS[1][is] = yS'(t0), is=1,...,Ns, and scale znS[1][is] by h.
If computing any quadratures, call fQ at (t0,znQ[0]), set
znQ[1] = fQ, and scale znQ[1] by h.
*/
f(tn, zn[0], zn[1], f_data);
nfe++;
if (sensi) {
wrk1 = tempv;
wrk2 = ftemp;
CVSensRhs(cv_mem, tn, zn[0], zn[1], znS[0], znS[1], tempv, ftemp);
}
if (quad) {
fQ(tn, zn[0], znQ[1], fQ_data);
nfQe++;
}
h = hin;
if ( (h != ZERO) && ((tout-tn)*h < ZERO) ) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_H0);
return (CV_ILL_INPUT);
}
if (h == ZERO) {
hOK = CVHin(cv_mem, tout);
if (!hOK) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_TOO_CLOSE);
return (CV_ILL_INPUT);
}
}
rh = ABS(h)*hmax_inv;
if (rh > ONE) h /= rh;
if (ABS(h) < hmin) h *= hmin/ABS(h);
/* Check for approach to tstop */
if (istop) {
if ( (tstop - tn)*h < ZERO ) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_TSTOP, tn);
return(CV_ILL_INPUT);
}
if ( (tn + h - tstop)*h > ZERO )
h = tstop - tn;
}
hscale = h;
h0u = h;
hprime = h;
N_VScale(h, zn[1], zn[1]);
if (sensi)
for (is=0; is<Ns; is++)
N_VScale(h, znS[1][is], znS[1][is]);
if (quad)
N_VScale(h, znQ[1], znQ[1]);
if (nrtfn > 0) {
ier = CVRcheck1(cv_mem);
if (ier != CV_SUCCESS) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_INIT_ROOT);
return(CV_ILL_INPUT);
}
}
} /* end first call block */
/* At following steps, perform stop tests */
if (nst > 0) {
/* First check for a root in the last step taken, other than the
last root found, if any. If task = CV_ONE_STEP and y(tn) was not
returned because of an intervening root, return y(tn) now. */
if (nrtfn > 0) {
irfndp = irfnd;
ier = CVRcheck2(cv_mem);
if (ier == CLOSERT) {
tretlast = *tret = tlo;
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_CLOSE_ROOTS, tlo);
return(CV_ILL_INPUT);
}
if (ier == RTFOUND) {
tretlast = *tret = tlo;
return(CV_ROOT_RETURN);
}
if (tn != tretlast) { /* Check remaining interval for roots */
ier = CVRcheck3(cv_mem);
if (ier == CV_SUCCESS) { /* no root found */
irfnd = 0;
if (irfndp == 1 && task == CV_ONE_STEP) {
tretlast = *tret = tn;
N_VScale(ONE, zn[0], yout);
return(CV_SUCCESS);
}
}
if (ier == RTFOUND) { /* a new root was found */
irfnd = 1;
tretlast = *tret = tlo;
return(CV_ROOT_RETURN);
}
}
} /* end of root stop check */
/* Test for tn past tstop */
if ( istop && ((tstop - tn)*h < ZERO) ) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_TSTOP, tn);
return (CV_ILL_INPUT);
}
/* In CV_NORMAL mode, test if tout was reached */
if ( (task == CV_NORMAL) && ((tn-tout)*h >= ZERO) ) {
tretlast = *tret = tout;
ier = CVodeGetDky(cv_mem, tout, 0, yout);
if (ier != CV_SUCCESS) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_TOUT, tout);
return (CV_ILL_INPUT);
}
return (CV_SUCCESS);
}
/* In CV_ONE_STEP mode, test if tn was returned */
if (task == CV_ONE_STEP && tretlast != tn) {
tretlast = *tret = tn;
N_VScale(ONE, zn[0], yout);
return(CV_SUCCESS);
}
/* Test for tn at tstop or near tstop */
if ( istop ) {
troundoff = FUZZ_FACTOR*uround*(ABS(tn) + ABS(h));
if ( ABS(tn - tstop) <= troundoff) {
ier = CVodeGetDky(cv_mem, tstop, 0, yout);
if (ier != CV_SUCCESS) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_TSTOP, tn);
return (CV_ILL_INPUT);
}
tretlast = *tret = tstop;
tn = tstop;
return (CV_TSTOP_RETURN);
}
if ( (tn + hprime - tstop)*h > ZERO ) {
hprime = tstop - tn;
eta = hprime/h;
}
} /* end of istop tests block */
} /* end stopping tests block at nst>0 */
/* Start looping for internal steps */
nstloc = 0;
loop {
next_h = h;
next_q = q;
/* Reset and check ewt */
if (nst > 0) {
ewtsetOK = CVEwtSet(cv_mem, zn[0]);
if (sensi)
ewtSsetOK = CVSensEwtSet(cv_mem, znS[0]);
else
ewtSsetOK = TRUE;
if (quad && errconQ)
ewtQsetOK = CVQuadEwtSet(cv_mem, znQ[0]);
else
ewtQsetOK = TRUE;
if ( (!ewtsetOK) || (!ewtSsetOK) || (!ewtQsetOK) ) {
if(!ewtsetOK) EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_EWT_NOW_BAD, tn);
if(!ewtSsetOK) EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_EWTS_NOW_BAD, tn);
if(!ewtQsetOK) EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_EWTQ_NOW_BAD, tn);
istate = CV_ILL_INPUT;
tretlast = *tret = tn;
N_VScale(ONE, zn[0], yout);
break;
}
}
/* Check for too many steps */
if (nstloc >= mxstep) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_MAX_STEPS, tn);
istate = CV_TOO_MUCH_WORK;
tretlast = *tret = tn;
N_VScale(ONE, zn[0], yout);
break;
}
/* Check for too much accuracy requested */
nrm = N_VWrmsNorm(zn[0], ewt);
if (quad && errconQ) {
nrm = CVQuadUpdateNorm(cv_mem, nrm, znQ[0], ewtQ);
}
if (sensi && errconS) {
nrm = CVSensUpdateNorm(cv_mem, nrm, znS[0], ewtS);
}
tolsf = uround * nrm;
if (tolsf > ONE) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_TOO_MUCH_ACC, tn);
istate = CV_TOO_MUCH_ACC;
tretlast = *tret = tn;
N_VScale(ONE, zn[0], yout);
tolsf *= TWO;
break;
} else {
tolsf = ONE;
}
/* Check for h below roundoff level in tn */
if (tn + hprime == tn) {
nhnil++;
if (nhnil <= mxhnil) EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_HNIL, tn, hprime);
if (nhnil == mxhnil)
{
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_HNIL_DONE);
// edc added this to avoid infinite loop... for some reason, this error
// wasn't properly causing a break out of the loop... CHECK THIS!!
istate = CVHandleFailure(cv_mem, REP_ERR_FAIL);
break;
}
}
/* Call CVStep to take a step */
kflag = CVStep(cv_mem);
/* Process failed step cases, and exit loop */
if (kflag != SUCCESS_STEP) {
istate = CVHandleFailure(cv_mem, kflag);
tretlast = *tret = tn;
N_VScale(ONE, zn[0], yout);
break;
}
nstloc++;
/* Check for root in last step taken. */
if (nrtfn > 0) {
ier = CVRcheck3(cv_mem);
if (ier == RTFOUND) { /* a new root was found */
irfnd = 1;
tretlast = *tret = tlo;
return(CV_ROOT_RETURN);
}
}
/* Check if tn is at tstop or near tstop */
if ( istop ) {
troundoff = FUZZ_FACTOR*uround*(ABS(tn) + ABS(h));
if ( ABS(tn - tstop) <= troundoff) {
(void) CVodeGetDky(cv_mem, tstop, 0, yout);
tretlast = *tret = tstop;
tn = tstop;
istate = CV_TSTOP_RETURN;
break;
}
if ( (tn + hprime - tstop)*h > ZERO ) {
hprime = tstop - tn;
eta = hprime/h;
}
}
/* Check if in one-step mode, and if so copy y and exit loop */
if (task == CV_ONE_STEP) {
istate = CV_SUCCESS;
tretlast = *tret = tn;
N_VScale(ONE, zn[0], yout);
next_q = qprime;
next_h = hprime;
break;
}
/* Check if tout reached, and if so interpolate and exit loop */
if ((tn-tout)*h >= ZERO) {
istate = CV_SUCCESS;
tretlast = *tret = tout;
(void) CVodeGetDky(cv_mem, tout, 0, yout);
next_q = qprime;
next_h = hprime;
break;
}
} /* end looping for internal steps */
/* Load optional output */
if (sensi && (ism==CV_STAGGERED1)) {
nniS = 0;
ncfnS = 0;
for (is=0; is<Ns; is++) {
nniS += nniS1[is];
ncfnS += ncfnS1[is];
}
}
return (istate);
}
/*-----------------------------------------------------------------*/
/*
* CVodeGetDky
*
* This routine computes the k-th derivative of the interpolating
* polynomial at the time t and stores the result in the vector dky.
* The formula is:
* q
* dky = SUM c(j,k) * (t - tn)^(j-k) * h^(-j) * zn[j] ,
* j=k
* where c(j,k) = j*(j-1)*...*(j-k+1), q is the current order, and
* zn[j] is the j-th column of the Nordsieck history array.
*
* This function is called by CVode with k = 0 and t = tout, but
* may also be called directly by the user.
*/
int CVodeGetDky(void *cvode_mem, realtype t, int k, N_Vector dky)
{
realtype s, c, r;
realtype tfuzz, tp, tn1;
int i, j;
CVodeMem cv_mem;
/* Check all inputs for legality */
if (cvode_mem == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_DKY_NO_MEM);
return (CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
if (dky == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_DKY);
return (CV_BAD_DKY);
}
if ((k < 0) || (k > q)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_K);
return (CV_BAD_K);
}
/* Allow for some slack */
tfuzz = FUZZ_FACTOR * uround * (ABS(tn) + ABS(hu));
if (hu < ZERO) tfuzz = -tfuzz;
tp = tn - hu - tfuzz;
tn1 = tn + tfuzz;
if ((t-tp)*(t-tn1) > ZERO) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_T, t, tn-hu, tn);
return (CV_BAD_T);
}
/* Sum the differentiated interpolating polynomial */
s = (t - tn) / h;
for (j=q; j >= k; j--) {
c = ONE;
for (i=j; i >= j-k+1; i--) c *= i;
if (j == q) {
N_VScale(c, zn[q], dky);
} else {
N_VLinearSum(c, zn[j], s, dky, dky);
}
}
if (k == 0) return (CV_SUCCESS);
r = RPowerI(h,-k);
N_VScale(r, dky, dky);
return (CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
/*
* CVodeGetQuad
*
* This routine extracts quadrature solution into yQout.
* This is just a wrapper that calls CvodeGEtQuadDky with k=0
*/
int CVodeGetQuad(void *cvode_mem, realtype t, N_Vector yQout)
{
return (CVodeGetQuadDky(cvode_mem,t,0,yQout));
}
/*-----------------------------------------------------------------*/
/*
* CVodeGetQuadDky
*
* CVodeQuadDky computes the kth derivative of the yQ function at
* time t, where tn-hu <= t <= tn, tn denotes the current
* internal time reached, and hu is the last internal step size
* successfully used by the solver. The user may request
* k=0, 1, ..., qu, where qu is the current order.
* The derivative vector is returned in dky. This vector
* must be allocated by the caller. It is only legal to call this
* function after a successful return from CVode with quadrature
* computation enabled.
*/
int CVodeGetQuadDky(void *cvode_mem, realtype t, int k, N_Vector dkyQ)
{
realtype s, c, r;
realtype tfuzz, tp, tn1;
int i, j;
CVodeMem cv_mem;
/* Check all inputs for legality */
if (cvode_mem == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_QDKY_NO_MEM);
return (CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
if(quad != TRUE) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_QDKY_NO_QUAD);
return (CV_NO_QUAD);
}
if (dkyQ == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_QBAD_DKY);
return (CV_BAD_DKY);
}
if ((k < 0) || (k > q)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_QBAD_K);
return (CV_BAD_K);
}
/* Allow for some slack */
tfuzz = FUZZ_FACTOR * uround * (ABS(tn) + ABS(hu));
if (hu < ZERO) tfuzz = -tfuzz;
tp = tn - hu - tfuzz;
tn1 = tn + tfuzz;
if ((t-tp)*(t-tn1) > ZERO) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_QBAD_T, t, tn-hu, tn);
// Oscill8 temporary fix: returning here causes the root finder
// to use invalid values of Y; CVODE doesn't check return value
// of this function for root finding and this isn't good!
//return (CV_BAD_T);
EXTERNAL_LIBRARY_FPRINTF_STDERR("CVode (Oscill8 note): skipping unchecked return from CVodeGetQuadDky routine");
}
/* Sum the differentiated interpolating polynomial */
s = (t - tn) / h;
for (j=q; j >= k; j--) {
c = ONE;
for (i=j; i >= j-k+1; i--) c *= i;
if (j == q) {
N_VScale(c, znQ[q], dkyQ);
} else {
N_VLinearSum(c, znQ[j], s, dkyQ, dkyQ);
}
}
if (k == 0) return (CV_SUCCESS);
r = RPowerI(h,-k);
N_VScale(r, dkyQ, dkyQ);
return (CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
/*
* CVodeGetSens
*
* This routine extracts sensitivity solution into ySout.
* This is just a wrapper that calls CvodeSensDky with k=0
*/
int CVodeGetSens(void *cvode_mem, realtype t, N_Vector *ySout)
{
return (CVodeGetSensDky(cvode_mem,t,0,ySout));
}
/*-----------------------------------------------------------------*/
/*
* CVodeGetSens1
*
* This routine extracts the is-th sensitivity solution into ySout.
* This is just a wrapper that calls CvodeSensDky1 with k=0
*/
int CVodeGetSens1(void *cvode_mem, realtype t, int is, N_Vector ySout)
{
return (CVodeGetSensDky1(cvode_mem,t,0,is,ySout));
}
/*-----------------------------------------------------------------*/
/*
* CVodeGetSensDky
*
* If the user calls directly CVodeSensDky then s must be allocated
* prior to this call. When CVodeSensDky is called by
* CVodeGetSens, only ier=CV_SUCCESS, ier=CV_NO_SENS, or
* ier=CV_BAD_T are possible.
*/
int CVodeGetSensDky(void *cvode_mem, realtype t, int k, N_Vector *dkyS)
{
int ier=CV_SUCCESS;
int is;
CVodeMem cv_mem;
if (cvode_mem == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_SDKY_NO_MEM);
return (CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
if (dkyS == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SBAD_DKYA);
return (CV_BAD_DKY);
}
for (is=0; is<Ns; is++) {
ier = CVodeGetSensDky1(cvode_mem,t,k,is+1,dkyS[is]);
if (ier!=CV_SUCCESS) break;
}
return (ier);
}
/*-----------------------------------------------------------------*/
/*
* CVodeGetSensDky1
*
* CVodeSensDky1 computes the kth derivative of the yS[is] function at
* time t, where tn-hu <= t <= tn, tn denotes the current
* internal time reached, and hu is the last internal step size
* successfully used by the solver. The user may request
* is=1, 2, ..., Ns and k=0, 1, ..., qu, where qu is the current
* order. The derivative vector is returned in dky. This vector
* must be allocated by the caller. It is only legal to call this
* function after a successful return from CVode with sensitivity
* computation enabled.
*/
int CVodeGetSensDky1(void *cvode_mem, realtype t, int k, int is,
N_Vector dkyS)
{
realtype s, c, r;
realtype tfuzz, tp, tn1;
int i, j;
CVodeMem cv_mem;
/* Check all inputs for legality */
if (cvode_mem == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR( MSGCVS_SDKY_NO_MEM);
return (CV_MEM_NULL);
}
cv_mem = (CVodeMem) cvode_mem;
if(sensi != TRUE) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SDKY_NO_SENSI);
return (CV_NO_SENS);
}
if (dkyS == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SBAD_DKY);
return (CV_BAD_DKY);
}
if ((k < 0) || (k > q)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SBAD_K);
return (CV_BAD_K);
}
if ((is < 1) || (is > Ns)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SBAD_IS);
return (CV_BAD_IS);
}
is--;
/* Allow for some slack */
tfuzz = FUZZ_FACTOR * uround * (ABS(tn) + ABS(hu));
if (hu < ZERO) tfuzz = -tfuzz;
tp = tn - hu - tfuzz;
tn1 = tn + tfuzz;
if ((t-tp)*(t-tn1) > ZERO) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_SBAD_T);
return (CV_BAD_T);
}
/* Sum the differentiated interpolating polynomial */
s = (t - tn) / h;
for (j=q; j >= k; j--) {
c = ONE;
for (i=j; i >= j-k+1; i--) c *= i;
if (j == q) {
N_VScale(c, znS[q][is], dkyS);
} else {
N_VLinearSum(c, znS[j][is], s, dkyS, dkyS);
}
}
if (k == 0) return (CV_SUCCESS);
r = RPowerI(h,-k);
N_VScale(r, dkyS, dkyS);
return (CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
/*
* CVodeFree
*
* This routine frees the problem memory allocated by CVodeMalloc.
* Such memory includes all the vectors allocated by CVAllocVectors,
* and the memory lmem for the linear solver (deallocated by a call
* to lfree), as well as (if Ns!=0) all memory allocated for
* sensitivity computations by CVodeSensMalloc.
*/
void CVodeFree(void *cvode_mem)
{
CVodeMem cv_mem;
cv_mem = (CVodeMem) cvode_mem;
if (cvode_mem == NULL) return;
CVFreeVectors(cv_mem);
CVodeQuadFree(cv_mem);
CVodeSensFree(cv_mem);
if (iter == CV_NEWTON && lfree != NULL) lfree(cv_mem);
if (nrtfn > 0) {
free(glo);
free(ghi);
free(groot);
free(iroots);
}
free(cv_mem);
}
/*-----------------------------------------------------------------*/
/*
* CVodeQuadFree
*
* CVodeQuadFree frees the problem memory in cvode_mem allocated
* for quadrature integration. Its only argument is the pointer
* cvode_mem returned by CVodeCreate.
*/
void CVodeQuadFree(void *cvode_mem)
{
CVodeMem cv_mem;
if (cvode_mem == NULL) return;
cv_mem = (CVodeMem) cvode_mem;
if(quad) {
CVQuadFreeVectors(cv_mem);
quad = FALSE;
}
}
/*-----------------------------------------------------------------*/
/*
* CVodeSensFree
*
* CVodeSensFree frees the problem memory in cvode_mem allocated
* for sensitivity analysis. Its only argument is the pointer
* cvode_mem returned by CVodeCreate.
*/
void CVodeSensFree(void *cvode_mem)
{
CVodeMem cv_mem;
if (cvode_mem == NULL) return;
cv_mem = (CVodeMem) cvode_mem;
if(sensi) {
if (atolSallocated)
CVSensFreeAtol(cv_mem, abstolS);
if (stgr1alloc) {
free(ncfS1);
free(ncfnS1);
free(nniS1);
}
CVSensFreeVectors(cv_mem);
sensi = FALSE;
}
}
/*=================================================================*/
/* PRIVATE FUNCTIONS IMPLEMENTATION */
/*=================================================================*/
/*
* CVCheckNvector
* This routine checks if all required vector operations are present.
* If any of them is missing it returns FALSE.
*/
static booleantype CVCheckNvector(N_Vector tmpl)
{
if((tmpl->ops->nvclone == NULL) ||
(tmpl->ops->nvdestroy == NULL) ||
(tmpl->ops->nvlinearsum == NULL) ||
(tmpl->ops->nvconst == NULL) ||
(tmpl->ops->nvprod == NULL) ||
(tmpl->ops->nvdiv == NULL) ||
(tmpl->ops->nvscale == NULL) ||
(tmpl->ops->nvabs == NULL) ||
(tmpl->ops->nvinv == NULL) ||
(tmpl->ops->nvaddconst == NULL) ||
(tmpl->ops->nvmaxnorm == NULL) ||
(tmpl->ops->nvwrmsnorm == NULL) ||
(tmpl->ops->nvmin == NULL))
return(FALSE);
else
return(TRUE);
}
/*-----------------------------------------------------------------*/
/*
* CVAllocVectors
*
* This routine allocates the CVODE vectors ewt, acor, tempv, ftemp, and
* zn[0], ..., zn[maxord]. The length of the vectors is the input
* parameter neq and the maximum order (needed to allocate zn) is the
* input parameter maxord. If all memory allocations are successful,
* CVAllocVectors returns TRUE. Otherwise all allocated memory is freed
* and CVAllocVectors returns FALSE.
* This routine also sets the optional outputs lrw and liw, which are
* (respectively) the lengths of the real and integer work spaces
* allocated here.
*/
static booleantype CVAllocVectors(CVodeMem cv_mem, N_Vector tmpl)
{
int i, j;
/* Allocate ewt, acor, tempv, ftemp */
ewt = N_VClone(tmpl);
if (ewt == NULL) return (FALSE);
acor = N_VClone(tmpl);
if (acor == NULL) {
N_VDestroy(ewt);
return (FALSE);
}
tempv = N_VClone(tmpl);
if (tempv == NULL) {
N_VDestroy(ewt);
N_VDestroy(acor);
return (FALSE);
}
ftemp = N_VClone(tmpl);
if (ftemp == NULL) {
N_VDestroy(tempv);
N_VDestroy(ewt);
N_VDestroy(acor);
return (FALSE);
}
/* Allocate zn[0] ... zn[maxord] */
for (j=0; j <= qmax; j++) {
zn[j] = N_VClone(tmpl);
if (zn[j] == NULL) {
N_VDestroy(ewt);
N_VDestroy(acor);
N_VDestroy(tempv);
N_VDestroy(ftemp);
for (i=0; i < j; i++) N_VDestroy(zn[i]);
return (FALSE);
}
}
/* Set solver workspace lengths */
lrw = (qmax + 5)*lrw1;
liw = (qmax + 5)*liw1;
return (TRUE);
}
/*-----------------------------------------------------------------*/
/*
* CVFreeVectors
*
* This routine frees the CVODE vectors allocated in CVAllocVectors.
*/
static void CVFreeVectors(CVodeMem cv_mem)
{
int j;
N_VDestroy(ewt);
N_VDestroy(acor);
N_VDestroy(tempv);
N_VDestroy(ftemp);
for (j=0; j <= qmax; j++) N_VDestroy(zn[j]);
}
/*-----------------------------------------------------------------*/
/*
* CVInitialSetup
*
* This routine performs input consistency checks at the first step.
* If needed, it also checks the linear solver module and calls the
* linear solver initialization routine.
*/
static int CVInitialSetup(CVodeMem cv_mem)
{
int ier;
booleantype ewtsetOK;
/* Solver initial setup */
ewtsetOK = CVEwtSet(cv_mem, zn[0]);
if (!ewtsetOK) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_EWT);
return(CV_ILL_INPUT);
}
/* Quadrature initial setup */
if (quad && errconQ) {
if ( (reltolQ == NULL) || (abstolQ == NULL) ) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_NO_QUADTOL);
return(CV_ILL_INPUT);
}
/* Load ewtQ */
ewtsetOK = CVQuadEwtSet(cv_mem, znQ[0]);
if (!ewtsetOK) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_EWTQ);
return (CV_ILL_INPUT);
}
}
if (!quad) errconQ = FALSE;
/* Forward sensitivity initial setup */
if (sensi) {
/* Check if ism and ifS agree */
if ((ism==CV_STAGGERED1) && (ifS==CV_ALLSENS)) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ISM_IFS);
return (CV_ILL_INPUT);
}
/* Load ewtS */
ewtsetOK = CVSensEwtSet(cv_mem, znS[0]);
if (!ewtsetOK) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_EWTS);
return (CV_ILL_INPUT);
}
}
if (!sensi) errconS = FALSE;
/* Check if lsolve function exists (if needed)
and call linit function (if it exists) */
if (iter == CV_NEWTON) {
if (lsolve == NULL) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_LSOLVE_NULL);
return (CV_ILL_INPUT);
}
if (linit != NULL) {
ier = linit(cv_mem);
if (ier != 0) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_LINIT_FAIL);
return (CV_ILL_INPUT);
}
}
}
return(CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
/*
* CVEwtSet
*
* This routine is responsible for setting the error weight vector ewt,
* according to tol_type, as follows:
*
* (1) ewt[i] = 1 / (*rtol * ABS(ycur[i]) + *atol), i=0,...,neq-1
* if tol_type = CV_SS
* (2) ewt[i] = 1 / (*rtol * ABS(ycur[i]) + atol[i]), i=0,...,neq-1
* if tol_type = CV_SV
*
* CVEwtSet returns TRUE if ewt is successfully set as above to a
* positive vector and FALSE otherwise. In the latter case, ewt is
* considered undefined after the FALSE return from CVEwtSet.
*
* All the real work is done in the routines CVEwtSetSS, CVEwtSetSV.
*/
static booleantype CVEwtSet(CVodeMem cv_mem, N_Vector ycur)
{
booleantype flag=TRUE;
switch (itol) {
case CV_SS:
flag = CVEwtSetSS(cv_mem, ycur);
break;
case CV_SV:
flag = CVEwtSetSV(cv_mem, ycur);
break;
}
return(flag);
}
/*-----------------------------------------------------------------*/
/*
* CVEwtSetSS
*
* This routine sets ewt as decribed above in the case tol_type = CV_SS.
* It tests for non-positive components before inverting. CVEwtSetSS
* returns TRUE if ewt is successfully set to a positive vector
* and FALSE otherwise. In the latter case, ewt is considered
* undefined after the FALSE return from CVEwtSetSS.
*/
static booleantype CVEwtSetSS(CVodeMem cv_mem, N_Vector ycur)
{
realtype rtoli, atoli;
rtoli = *reltol;
atoli = *((realtype *)abstol);
N_VAbs(ycur, tempv);
N_VScale(rtoli, tempv, tempv);
N_VAddConst(tempv, atoli, tempv);
if (N_VMin(tempv) <= ZERO) return (FALSE);
N_VInv(tempv, ewt);
return (TRUE);
}
/*-----------------------------------------------------------------*/
/*
* CVEwtSetSV
*
* This routine sets ewt as decribed above in the case tol_type = CV_SV.
* It tests for non-positive components before inverting. CVEwtSetSV
* returns TRUE if ewt is successfully set to a positive vector
* and FALSE otherwise. In the latter case, ewt is considered
* undefined after the FALSE return from CVEwtSetSV.
*/
static booleantype CVEwtSetSV(CVodeMem cv_mem, N_Vector ycur)
{
realtype rtoli;
rtoli = *reltol;
N_VAbs(ycur, tempv);
N_VLinearSum(rtoli, tempv, ONE, (N_Vector)abstol, tempv);
if (N_VMin(tempv) <= ZERO) return (FALSE);
N_VInv(tempv, ewt);
return (TRUE);
}
/*
* -----------------------------------------------------------------
* PRIVATE FUNCTIONS FOR QUADRATURES
* -----------------------------------------------------------------
*/
/*
* CVodeQuadAllocVectors
*
* NOTE: Space for ewtQ is allocated even when errconQ=FALSE,
* although in this case, ewtQ is never used. The reason for this
* decision is to allow the user to re-initialize the quadrature
* computation with errconQ=TRUE, after an initialization with
* errconQ=FALSE, without new memory allocation within
* CVodeQuadReInit.
*/
static booleantype CVQuadAllocVectors(CVodeMem cv_mem, N_Vector tmpl)
{
int i, j;
/* Allocate ewtQ */
ewtQ = N_VClone(tmpl);
if (ewtQ == NULL) {
return (FALSE);
}
/* Allocate acorQ */
acorQ = N_VClone(tmpl);
if (acorQ == NULL) {
N_VDestroy(ewtQ);
return (FALSE);
}
/* Allocate yQ */
yQ = N_VClone(tmpl);
if (yQ == NULL) {
N_VDestroy(ewtQ);
N_VDestroy(acorQ);
return (FALSE);
}
/* Allocate tempvQ */
tempvQ = N_VClone(tmpl);
if (tempvQ == NULL) {
N_VDestroy(ewtQ);
N_VDestroy(acorQ);
N_VDestroy(yQ);
return (FALSE);
}
/* Allocate zQn[0] ... zQn[maxord] */
for (j=0; j <= qmax; j++) {
znQ[j] = N_VClone(tmpl);
if (znQ[j] == NULL) {
N_VDestroy(ewtQ);
N_VDestroy(acorQ);
N_VDestroy(yQ);
N_VDestroy(tempvQ);
for (i=0; i < j; i++) N_VDestroy(znQ[i]);
return (FALSE);
}
}
/* Update solver workspace lengths */
lrw += (qmax + 5)*lrw1Q;
liw += (qmax + 5)*liw1Q;
return(TRUE);
}
/*-----------------------------------------------------------------*/
static booleantype CVQuadEwtSet(CVodeMem cv_mem, N_Vector qcur)
{
booleantype flag=TRUE;
switch (itolQ) {
case CV_SS:
flag = CVQuadEwtSetSS(cv_mem, qcur);
break;
case CV_SV:
flag = CVQuadEwtSetSV(cv_mem, qcur);
break;
}
return(flag);
}
/*-----------------------------------------------------------------*/
static booleantype CVQuadEwtSetSS(CVodeMem cv_mem, N_Vector qcur)
{
realtype rtoli, atoli;
rtoli = *reltolQ;
atoli = *((realtype *)abstolQ);
N_VAbs(qcur, tempvQ);
N_VScale(rtoli, tempvQ, tempvQ);
N_VAddConst(tempvQ, atoli, tempvQ);
if (N_VMin(tempvQ) <= ZERO) return (FALSE);
N_VInv(tempvQ, ewtQ);
return (TRUE);
}
/*-----------------------------------------------------------------*/
static booleantype CVQuadEwtSetSV(CVodeMem cv_mem, N_Vector qcur)
{
realtype rtoli;
rtoli = *reltolQ;
N_VAbs(qcur, tempvQ);
N_VLinearSum(rtoli, tempvQ, ONE, (N_Vector)abstolQ, tempvQ);
if (N_VMin(tempvQ) <= ZERO) return (FALSE);
N_VInv(tempvQ, ewtQ);
return (TRUE);
}
/*-----------------------------------------------------------------*/
static void CVQuadFreeVectors(CVodeMem cv_mem)
{
int j;
N_VDestroy(ewtQ);
N_VDestroy(acorQ);
N_VDestroy(yQ);
N_VDestroy(tempvQ);
for (j=0; j<=qmax; j++) N_VDestroy(znQ[j]);
}
/*
* -----------------------------------------------------------------
* PRIVATE FUNCTIONS FOR SENSITIVITY ANALYSIS
* -----------------------------------------------------------------
*/
static int CVSensTestTolerances(CVodeMem cv_mem)
{
int is;
booleantype neg_abstol;
realtype *atolSS;
N_Vector *atolSV;
neg_abstol = TRUE;
if (*reltolS<ZERO) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_RELTOLS);
return(CV_ILL_INPUT);
}
if (itolS == CV_SS) {
atolSS = (realtype *) abstolS;
for (is=0; is<Ns; is++)
if (atolSS[is] < ZERO) {neg_abstol = TRUE; break;}
} else {
atolSV = (N_Vector *) abstolS;
for (is=0; is<Ns; is++)
if (N_VMin(atolSV[is]) < ZERO) {neg_abstol = TRUE; break;}
}
if (neg_abstol) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_ABSTOLS);
return(CV_ILL_INPUT);
}
return(CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
static int CVSensSetTolerances(CVodeMem cv_mem)
{
booleantype allocOK, setOK;
setOK = TRUE;
itolS = itol;
reltolS = reltol;
if (!atolSallocated) {
allocOK = CVSensAllocAtol(cv_mem, &abstolS);
if (!allocOK) {
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_ATOLS_MEM_FAIL);
return(CV_ILL_INPUT);
}
atolSallocated = TRUE;
}
switch(itolS) {
case CV_SS:
setOK = CVSensSetAtolSS(cv_mem, (realtype *) abstolS);
break;
case CV_SV:
setOK = CVSensSetAtolSV(cv_mem, (N_Vector *) abstolS);
break;
}
if (!setOK) {
CVSensFreeAtol(cv_mem, abstolS);
EXTERNAL_LIBRARY_FPRINTF_STDERR(MSGCVS_BAD_PBAR);
return(CV_ILL_INPUT);
}
return(CV_SUCCESS);
}
/*-----------------------------------------------------------------*/
static booleantype CVSensAllocAtol(CVodeMem cv_mem, void **atolSPtr)
{
switch (itolS) {
case CV_SS:
*atolSPtr = (void *)malloc(Ns*sizeof(realtype));
break;
case CV_SV:
*atolSPtr = (void *)N_VCloneVectorArray(Ns, tempv);
lrw += Ns*lrw1;
liw += Ns*liw1;
break;
}
if (*atolSPtr==NULL) return (FALSE);
else return (TRUE);
}
/*-----------------------------------------------------------------*/
static void CVSensFreeAtol(CVodeMem cv_mem, void *atolS)
{
switch (itolS) {
case CV_SS:
free((realtype*)atolS);
break;
case CV_SV:
N_VDestroyVectorArray((N_Vector *)atolS, Ns);
break;
}
}
/*-----------------------------------------------------------------*/
static booleantype CVSensSetAtolSS(CVodeMem cv_mem, realtype *atolS)
{
int is, which;
realtype pb, rpb;
for (is=0; is<Ns; is++) {
if (plist!=NULL) which = abs(plist[is]) - 1;
else which = is;
if (pbar == NULL) pb = ONE;
else pb = ABS(pbar[which]);
if (pb == ZERO) return (FALSE);
rpb = ONE/pb;
atolS[is] = *((realtype *)abstol) * rpb;
}
return (TRUE);
}
/*-----------------------------------------------------------------*/
static booleantype CVSensSetAtolSV(CVodeMem cv_mem, N_Vector *atolS)
{
int is, which;
realtype pb, rpb;
for (is=0; is<Ns; is++) {
if (plist!=NULL) which = abs(plist[is]) - 1;
else which = is;
if (pbar == NULL) pb = ONE;
else pb = ABS(pbar[which]);
if (pb == ZERO) return (FALSE);
rpb = ONE/pb;
N_VScale(rpb, (N_Vector)abstol, atolS[is]);
}
return (TRUE);
}
/*-----------------------------------------------------------------*/
/*
* CVSensAllocVectors
*
* Create (through duplication) N_Vectors used for sensitivity analysis,
* using the N_Vector 'tmpl' as a template.
*/
static booleantype CVSensAllocVectors(CVodeMem cv_mem, N_Vector tmpl)
{
int i, j;
/* Allocate ewtS */
ewtS = N_VCloneVectorArray(Ns, tmpl);
if (ewtS == NULL) {
return (FALSE);
}
/* Allocate acorS */
acorS = N_VCloneVectorArray(Ns, tmpl);
if (acorS == NULL) {
N_VDestroyVectorArray(ewtS, Ns);
return (FALSE);
}
/* Allocate tempvS */
tempvS = N_VCloneVectorArray(Ns, tmpl);
if (tempvS == NULL) {
N_VDestroyVectorArray(ewtS, Ns);
N_VDestroyVectorArray(acorS, Ns);
return (FALSE);
}
/* Allocate ftempS */
ftempS = N_VCloneVectorArray(Ns, tmpl);
if (ftempS == NULL) {
N_VDestroyVectorArray(ewtS, Ns);
N_VDestroyVectorArray(acorS, Ns);
N_VDestroyVectorArray(tempvS, Ns);
return (FALSE);
}
/* Allocate znS */
for (j=0; j<=qmax; j++) {
znS[j] = N_VCloneVectorArray(Ns, tmpl);
if (znS[j] == NULL) {
N_VDestroyVectorArray(ewtS, Ns);
N_VDestroyVectorArray(acorS, Ns);
N_VDestroyVectorArray(tempvS, Ns);
N_VDestroyVectorArray(ftempS, Ns);
for (i=0; i<j; i++) N_VDestroyVectorArray(znS[i], Ns);
return (FALSE);
}
}
/* Update solver workspace lengths */
lrw += (qmax + 5)*Ns*lrw1;
liw += (qmax + 5)*Ns*liw1;
return (TRUE);
}
/*-----------------------------------------------------------------*/
static void CVSensFreeVectors(CVodeMem cv_mem)
{
int j;
N_VDestroyVectorArray(ewtS, Ns);
N_VDestroyVectorArray(acorS, Ns);
N_VDestroyVectorArray(tempvS, Ns);
N_VDestroyVectorArray(ftempS, Ns);
for (j=0; j<=qmax; j++) N_VDestroyVectorArray(znS[j], Ns);
}
/*-----------------------------------------------------------------*/
static booleantype CVSensEwtSet(CVodeMem cv_mem, N_Vector *yScur)
{
booleantype flag=TRUE;
switch (itolS) {
case CV_SS:
flag = CVSensEwtSetSS(cv_mem, yScur);
break;
case CV_SV:
flag = CVSensEwtSetSV(cv_mem, yScur);
break;
}
return(flag);
}
/*-----------------------------------------------------------------*/
static booleantype CVSensEwtSetSS(CVodeMem cv_mem, N_Vector *yScur)
{
int is;
realtype rtoli, atoli;
for (is=0; is<Ns; is++) {
rtoli = *reltolS;
atoli = ((realtype *)abstolS)[is];
N_VAbs(yScur[is], tempv);
N_VScale(rtoli, tempv, tempv);
N_VAddConst(tempv, atoli, tempv);
if (N_VMin(tempv) <= ZERO) return (FALSE);
N_VInv(tempv, ewtS[is]);
}
return (TRUE);
}
/*-----------------------------------------------------------------*/
static booleantype CVSensEwtSetSV(CVodeMem cv_mem, N_Vector *yScur)
{
int is;
realtype rtoli;
for (is=0; is<Ns; is++) {
rtoli = *reltolS;
N_VAbs(yScur[is], tempv);
N_VLinearSum(rtoli, tempv, ONE, ((N_Vector *)abstolS)[is], tempv);
if (N_VMin(tempv) <= ZERO) return (FALSE);
N_VInv(tempv, ewtS[is]);
}
return (TRUE);
}
/*
* -----------------------------------------------------------------
* PRIVATE FUNCTIONS FOR CVODE
* -----------------------------------------------------------------
*/
/*
* CVHin
*
* This routine computes a tentative initial step size h0.
* If tout is too close to tn (= t0), then CVHin returns FALSE and
* h remains uninitialized. Otherwise, CVHin sets h to the chosen
* value h0 and returns TRUE.
* The algorithm used seeks to find h0 as a solution of
* (WRMS norm of (h0^2 ydd / 2)) = 1,
* where ydd = estimated second derivative of y.
*/
static booleantype CVHin(CVodeMem cv_mem, realtype tout)
{
int sign, count;
realtype tdiff, tdist, tround, hlb, hub;
realtype hg, hgs, hnew, hrat, h0, yddnrm;
/* Test for tout too close to tn */
if ((tdiff = tout-tn) == ZERO) return (FALSE);
sign = (tdiff > ZERO) ? 1 : -1;
tdist = ABS(tdiff);
tround = uround * MAX(ABS(tn), ABS(tout));
if (tdist < TWO*tround) return (FALSE);
/* Set lower and upper bounds on h0, and take geometric mean
Exit with this value if the bounds cross each other */
hlb = HLB_FACTOR * tround;
hub = CVUpperBoundH0(cv_mem, tdist);
hg = RSqrt(hlb*hub);
if (hub < hlb) {
if (sign == -1) hg = -hg;
h = hg;
return (TRUE);
}
/* Loop up to MAX_ITERS times to find h0.
Stop if new and previous values differ by a factor < 2.
Stop if hnew/hg > 2 after one iteration, as this probably means
that the ydd value is bad because of cancellation error. */
count = 0;
loop {
hgs = hg*sign;
yddnrm = CVYddNorm(cv_mem, hgs);
hnew = (yddnrm*hub*hub > TWO) ? RSqrt(TWO/yddnrm) : RSqrt(hg*hub);
count++;
if (count >= MAX_ITERS) break;
hrat = hnew/hg;
if ((hrat > HALF) && (hrat < TWO)) break;
if ((count >= 2) && (hrat > TWO)) {
hnew = hg;
break;
}
hg = hnew;
}
/* Apply bounds, bias factor, and attach sign */
h0 = H_BIAS*hnew;
if (h0 < hlb) h0 = hlb;
if (h0 > hub) h0 = hub;
if (sign == -1) h0 = -h0;
h = h0;
return (TRUE);
}
/*-----------------------------------------------------------------*/
/*
* CVUpperBoundH0
*
* This routine sets an upper bound on abs(h0) based on
* tdist = tn - t0 and the values of y[i]/y'[i].
*/
static realtype CVUpperBoundH0(CVodeMem cv_mem, realtype tdist)
{
booleantype vectorAtol, vectorAtolQ, vectorAtolS;
realtype atoli, hub_inv, hub;
N_Vector temp1, temp2;
realtype hubQ_inv;
N_Vector tempQ1, tempQ2;
realtype *atolSS=NULL, hubS_inv;
N_Vector *atolSV=NULL;
int is;
vectorAtol = (itol == CV_SV);
vectorAtolQ = (itolQ == CV_SV);
vectorAtolS = (itolS == CV_SV);
temp1 = tempv;
temp2 = acor;
N_VAbs(zn[0], temp1);
N_VAbs(zn[1], temp2);
if (vectorAtol) {
N_VLinearSum(HUB_FACTOR, temp1, ONE, (N_Vector)abstol, temp1);
} else {
atoli = *((realtype *) abstol);
N_VScale(HUB_FACTOR, temp1, temp1);
N_VAddConst(temp1, atoli, temp1);
}
N_VDiv(temp2, temp1, temp1);
hub_inv = N_VMaxNorm(temp1);
if (quad && errconQ) {
tempQ1 = tempvQ;
tempQ2 = acorQ;
N_VAbs(znQ[0], tempQ1);
N_VAbs(znQ[1], tempQ2);
if (vectorAtolQ) {
N_VLinearSum(HUB_FACTOR, tempQ1, ONE, (N_Vector)abstolQ, tempQ1);
} else {
atoli = *((realtype *) abstolQ);
N_VScale(HUB_FACTOR, tempQ1, tempQ1);
N_VAddConst(tempQ1, atoli, tempQ1);
}
N_VDiv(tempQ2, tempQ1, tempQ1);
hubQ_inv = N_VMaxNorm(tempQ1);
if (hubQ_inv > hub_inv) hub_inv = hubQ_inv;
}
if (sensi && errconS) {
if (vectorAtolS) atolSV = (N_Vector *)abstolS;
else atolSS = (realtype *)abstolS;
for (is=0; is<Ns; is++) {
N_VAbs(znS[0][is], temp1);
N_VAbs(znS[1][is], temp2);
if (vectorAtolS) {
N_VLinearSum(HUB_FACTOR, temp1, ONE, atolSV[is], temp1);
} else {
N_VScale(HUB_FACTOR, temp1, temp1);
N_VAddConst(temp1, atolSS[is], temp1);
}
N_VDiv(temp2, temp1, temp1);
hubS_inv = N_VMaxNorm(temp1);
if (hubS_inv > hub_inv) hub_inv = hubS_inv;
}
}
hub = HUB_FACTOR*tdist;
if (hub*hub_inv > ONE) hub = ONE/hub_inv;
return (hub);
}
/*-----------------------------------------------------------------*/
/*
* CVYddNorm
*
* This routine computes an estimate of the second derivative of y
* using a difference quotient, and returns its WRMS norm.
*/
static realtype CVYddNorm(CVodeMem cv_mem, realtype hg)
{
realtype yddnrm;
int is;
N_Vector wrk1, wrk2;
/* y <- h*y'(t) + y(t) */
N_VLinearSum(hg, zn[1], ONE, zn[0], y);
if (sensi && errconS)
for (is=0; is<Ns; is++)
N_VLinearSum(hg, znS[1][is], ONE, znS[0][is], yS[is]);
/* tempv <- f(t+h, h*y'(t)+y(t)) */
f(tn+hg, y, tempv, f_data);
nfe++;
if (quad && errconQ) {
fQ(tn+hg, y, tempvQ, fQ_data);
nfQe++;
}
if (sensi && errconS) {
wrk1 = ftemp;
wrk2 = acor;
CVSensRhs(cv_mem, tn+hg, y, tempv, yS, tempvS, wrk1, wrk2);
}
/* tempv <- f(t+h, h*y'(t)+y(t)) - y'(t) */
/* tempv <- ydd */
N_VLinearSum(ONE, tempv, -ONE, zn[1], tempv);
N_VScale(ONE/hg, tempv, tempv);
if (quad && errconQ) {
N_VLinearSum(ONE, tempvQ, -ONE, znQ[1], tempvQ);
N_VScale(ONE/hg, tempvQ, tempvQ);
}
if (sensi && errconS)
for (is=0; is<Ns; is++) {
N_VLinearSum(ONE, tempvS[is], -ONE, znS[1][is], tempvS[is]);
N_VScale(ONE/hg, tempvS[is], tempvS[is]);
}
/* Estimate ||y''|| */
yddnrm = N_VWrmsNorm(tempv, ewt);
if (quad && errconQ) {
yddnrm = CVQuadUpdateNorm(cv_mem, yddnrm, tempvQ, ewtQ);
}
if (sensi && errconS) {
yddnrm = CVSensUpdateNorm(cv_mem, yddnrm, tempvS, ewtS);
}
return (yddnrm);
}
/*-----------------------------------------------------------------*/
/*
* CVStep
*
* This routine performs one internal cvode step, from tn to tn + h.
* It calls other routines to do all the work.
*
* The main operations done here are as follows:
* - preliminary adjustments if a new step size was chosen;
* - prediction of the Nordsieck history array zn at tn + h;
* - setting of multistep method coefficients and test quantities;
* - solution of the nonlinear system;
* - testing the local error;
* - updating zn and other state data if successful;
* - resetting stepsize and order for the next step.
* - if SLDET is on, check for stability, reduce order if necessary.
* On a failure in the nonlinear system solution or error test, the
* step may be reattempted, depending on the nature of the failure.
*/
static int CVStep(CVodeMem cv_mem)
{
realtype saved_t, dsm, dsmS, dsmQ;
int ncf, nef, kflag, nflag, ncfS, nefS, nefQ;
booleantype do_sensi_stg, do_sensi_stg1, passed;
int is;
saved_t = tn;
ncf = nef = 0;
nflag = FIRST_CALL;
if (quad) nefQ = 0;
/* Are we computing sensitivities with a staggered approach? */
do_sensi_stg = (sensi && (ism==CV_STAGGERED));
do_sensi_stg1 = (sensi && (ism==CV_STAGGERED1));
if (do_sensi_stg) {
ncfS = 0; /* local conv. failure counter for sensitivities */
nefS = 0; /* local err.test failure counter for sensi. is */
}
if (do_sensi_stg1) {
for (is=0; is<Ns; is++)
ncfS1[is] = 0;
nefS = 0;
}
if ((nst > 0) && (hprime != h)) CVAdjustParams(cv_mem);
/* Looping point for attempts to take a step */
loop {
CVPredict(cv_mem);
CVSet(cv_mem);
nflag = CVNls(cv_mem, nflag);
kflag = CVHandleNFlag(cv_mem, &nflag, saved_t, &ncf, &ncfn);
/* Go back in loop if we need to predict again */
if (kflag == PREDICT_AGAIN) continue;
/* Return if nonlinear solve failed and recovery not possible. */
if (kflag != DO_ERROR_TEST) return (kflag);
passed = CVDoErrorTest(cv_mem, &nflag, &kflag, saved_t, &nef, &dsm);
/* Return if error test failed and recovery not possible. */
if ((!passed) && (kflag == REP_ERR_FAIL)) return (kflag);
/* Retry step if error test failed, nflag == PREV_ERR_FAIL */
if (!passed) continue;
/* passed = TRUE, kflag = DO_ERROR_TEST, nflag = SOLVED */
/* Correct the quadrature variables */
if (quad) {
/* Save quadrature correction in acorQ */
fQ(tn, y, acorQ, fQ_data);
N_VLinearSum(h, acorQ, -ONE, znQ[1], acorQ);
N_VScale(rl1, acorQ, acorQ);
/* Apply correction to quadrature variables */
N_VLinearSum(ONE, znQ[0], ONE, acorQ, yQ);
/* Error test on quadratures */
if (errconQ) {
acnrmQ = N_VWrmsNorm(acorQ, ewtQ);
passed = CVQuadDoErrorTest(cv_mem, &nflag, &kflag, saved_t, &nefQ, &dsmQ);
if ((!passed) && (kflag == REP_ERR_FAIL)) return (kflag);
if (!passed) continue;
/* update 'dsm' with 'dsmQ' (to be used in CVPrepareNextStep) */
dsm = CVQuadUpdateDsm(cv_mem, dsm, dsmQ);
}
}
/* CV_STAGGERED approach for sensitivities */
if (do_sensi_stg) {
/* Reset counters for states */
ncf = nef = 0;
/* Evaluate f at converged y */
f(tn, y, ftemp, f_data);
nfe++;
/* Nonlinear solve for sensitivities (all-at-once) */
nflag = CVStgrNls(cv_mem);
kflag = CVHandleNFlag(cv_mem, &nflag, saved_t, &ncfS, &ncfnS);
if (kflag == PREDICT_AGAIN) continue;
if (kflag != DO_ERROR_TEST) return (kflag);
/* Error test on sensitivities */
if (errconS) {
passed = CVStgrDoErrorTest(cv_mem,&nflag,&kflag,saved_t,&nefS,&dsmS);
if ((!passed) && (kflag == REP_ERR_FAIL)) return (kflag);
if (!passed) continue;
/* update 'dsm' with 'dsmS' (to be used in CVPrepareNextStep) */
dsm = CVStgrUpdateDsm(cv_mem, dsm, dsmS);
}
}
/* CV_STAGGERED1 approach for sensitivities */
if (do_sensi_stg1) {
/* Reset counters for states */
ncf = nef = 0;
/* Evaluate f at converged y */
f(tn, y, ftemp, f_data);
nfe++;
/* Nonlinear solve for sensitivities (one-by-one) */
for (is=0; is<Ns; is++) {
nflag = CVStgr1Nls(cv_mem, is);
kflag = CVHandleNFlag(cv_mem, &nflag, saved_t, &ncfS1[is], &ncfnS1[is]);
if (kflag != DO_ERROR_TEST) break;
}
if (kflag == PREDICT_AGAIN) continue;
if (kflag != DO_ERROR_TEST) return (kflag);
/* Error test on sensitivities */
if (errconS) {
acnrmS = CVSensNorm(cv_mem, acorS, ewtS);
passed = CVStgrDoErrorTest(cv_mem,&nflag,&kflag,saved_t,&nefS,&dsmS);
if ((!passed) && (kflag == REP_ERR_FAIL)) return (kflag);
if (!passed) continue;
/* update 'dsm' with 'dsmS' (to be used in CVPrepareNextStep) */
dsm = CVStgrUpdateDsm(cv_mem, dsm, dsmS);
}
}
/* Everything went fine; exit loop */
break;
}
/* Nonlinear system solve and error test were both successful.
Update data, and consider change of step and/or order. */
CVCompleteStep(cv_mem);
CVPrepareNextStep(cv_mem, dsm);
/* If Stablilty Limit Detection is turned on, call stability limit
detection routine for possible order reduction. */
if (sldeton) CVBDFStab(cv_mem);
etamax = (nst <= SMALL_NST) ? ETAMX2 : ETAMX3;
/* Finally, we rescale the acor array to be the
estimated local error vector. */
N_VScale(ONE/tq[2], acor, acor);
if (quad)
N_VScale(ONE/tq[2], acorQ, acorQ);
if (sensi)
for (is=0; is<Ns; is++)
N_VScale(ONE/tq[2], acorS[is], acorS[is]);
return (SUCCESS_STEP);
}
/*-----------------------------------------------------------------*/
/*
* CVAdjustParams
*
* This routine is called when a change in step size was decided upon,
* and it handles the required adjustments to the history array zn.
* If there is to be a change in order, we call CVAdjustOrder and reset
* q, L = q+1, and qwait. Then in any case, we call CVRescale, which
* resets h and rescales the Nordsieck array.
*/
static void CVAdjustParams(CVodeMem cv_mem)
{
if (qprime != q) {
CVAdjustOrder(cv_mem, qprime-q);
q = qprime;
L = q+1;
qwait = L;
}
CVRescale(cv_mem);
}
/*-----------------------------------------------------------------*/
/*
* CVAdjustOrder
*
* This routine is a high level routine which handles an order
* change by an amount deltaq (= +1 or -1). If a decrease in order
* is requested and q==2, then the routine returns immediately.
* Otherwise CVAdjustAdams or CVAdjustBDF is called to handle the
* order change (depending on the value of lmm).
*/
static void CVAdjustOrder(CVodeMem cv_mem, int deltaq)
{
if ((q==2) && (deltaq != 1)) return;
switch(lmm){
case CV_ADAMS:
CVAdjustAdams(cv_mem, deltaq);
break;
case CV_BDF:
CVAdjustBDF(cv_mem, deltaq);
break;
}
}
/*-----------------------------------------------------------------*/
/*
* CVAdjustAdams
*
* This routine adjusts the history array on a change of order q by
* deltaq, in the case that lmm == CV_ADAMS.
*/
static void CVAdjustAdams(CVodeMem cv_mem, int deltaq)
{
int i, j;
int is;
realtype xi, hsum;
/* On an order increase, set new column of zn to zero and return */
if (deltaq==1) {
N_VConst(ZERO, zn[L]);
if (quad)
N_VConst(ZERO, znQ[L]);
if (sensi)
for (is=0; is<Ns; is++)
N_VConst(ZERO, znS[L][is]);
return;
}
/*
On an order decrease, each zn[j] is adjusted by a multiple of zn[q].
The coeffs. in the adjustment are the coeffs. of the polynomial:
x
q * INT { u * ( u + xi_1 ) * ... * ( u + xi_{q-2} ) } du
0
where xi_j = [t_n - t_(n-j)]/h => xi_0 = 0
*/
for (i=0; i <= qmax; i++) l[i] = ZERO;
l[1] = ONE;
hsum = ZERO;
for (j=1; j <= q-2; j++) {
hsum += tau[j];
xi = hsum / hscale;
for (i=j+1; i >= 1; i--) l[i] = l[i]*xi + l[i-1];
}
for (j=1; j <= q-2; j++) l[j+1] = q * (l[j] / (j+1));
for (j=2; j < q; j++)
N_VLinearSum(-l[j], zn[q], ONE, zn[j], zn[j]);
if (quad)
for (j=2; j < q; j++)
N_VLinearSum(-l[j], znQ[q], ONE, znQ[j], znQ[j]);
if (sensi)
for (is=0; is<Ns; is++)
for (j=2; j < q; j++)
N_VLinearSum(-l[j], znS[q][is], ONE, znS[j][is], znS[j][is]);
}
/*-----------------------------------------------------------------*/
/*
* CVAdjustBDF
*
* This is a high level routine which handles adjustments to the
* history array on a change of order by deltaq in the case that
* lmm == CV_BDF. CVAdjustBDF calls CVIncreaseBDF if deltaq = +1 and
* CVDecreaseBDF if deltaq = -1 to do the actual work.
*/
static void CVAdjustBDF(CVodeMem cv_mem, int deltaq)
{
switch(deltaq) {
case 1:
CVIncreaseBDF(cv_mem);
return;
case -1:
CVDecreaseBDF(cv_mem);
return;
}
}
/*-----------------------------------------------------------------*/
/*
* CVIncreaseBDF
*
* This routine adjusts the history array on an increase in the
* order q in the case that lmm == CV_BDF.
* A new column zn[q+1] is set equal to a multiple of the saved
* vector (= acor) in zn[qmax]. Then each zn[j] is adjusted by
* a multiple of zn[q+1]. The coefficients in the adjustment are the
* coefficients of the polynomial x*x*(x+xi_1)*...*(x+xi_j),
* where xi_j = [t_n - t_(n-j)]/h.
*/
static void CVIncreaseBDF(CVodeMem cv_mem)
{
realtype alpha0, alpha1, prod, xi, xiold, hsum, A1;
int i, j;
int is;
for (i=0; i <= qmax; i++) l[i] = ZERO;
l[2] = alpha1 = prod = xiold = ONE;
alpha0 = -ONE;
hsum = hscale;
if (q > 1) {
for (j=1; j < q; j++) {
hsum += tau[j+1];
xi = hsum / hscale;
prod *= xi;
alpha0 -= ONE / (j+1);
alpha1 += ONE / xi;
for (i=j+2; i >= 2; i--) l[i] = l[i]*xiold + l[i-1];
xiold = xi;
}
}
A1 = (-alpha0 - alpha1) / prod;
/*
zn[qmax] contains the value Delta_n = y_n - y_n(0)
This value was stored there at the previous successful
step (in CVCompleteStep)
A1 contains dbar = (1/xi* - 1/xi_q)/prod(xi_j)
*/
N_VScale(A1, zn[qmax], zn[L]);
for (j=2; j <= q; j++)
N_VLinearSum(l[j], zn[L], ONE, zn[j], zn[j]);
if (quad) {
N_VScale(A1, znQ[qmax], znQ[L]);
for (j=2; j <= q; j++)
N_VLinearSum(l[j], znQ[L], ONE, znQ[j], znQ[j]);
}
if (sensi) {
for (is=0; is<Ns; is++) {
N_VScale(A1, znS[qmax][is], znS[L][is]);
for (j=2; j <= q; j++)
N_VLinearSum(l[j], znS[L][is], ONE, znS[j][is], znS[j][is]);
}
}
}
/*-----------------------------------------------------------------*/
/*
* CVDecreaseBDF
*
* This routine adjusts the history array on a decrease in the
* order q in the case that lmm == CV_BDF.
* Each zn[j] is adjusted by a multiple of zn[q]. The coefficients
* in the adjustment are the coefficients of the polynomial
* x*x*(x+xi_1)*...*(x+xi_j), where xi_j = [t_n - t_(n-j)]/h.
*/
static void CVDecreaseBDF(CVodeMem cv_mem)
{
realtype hsum, xi;
int i, j;
int is;
for (i=0; i <= qmax; i++) l[i] = ZERO;
l[2] = ONE;
hsum = ZERO;
for (j=1; j <= q-2; j++) {
hsum += tau[j];
xi = hsum /hscale;
for (i=j+2; i >= 2; i--) l[i] = l[i]*xi + l[i-1];
}
for (j=2; j < q; j++)
N_VLinearSum(-l[j], zn[q], ONE, zn[j], zn[j]);
if (quad) {
for (j=2; j < q; j++)
N_VLinearSum(-l[j], znQ[q], ONE, znQ[j], znQ[j]);
}
if (sensi) {
for (is=0; is<Ns; is++)
for (j=2; j < q; j++)
N_VLinearSum(-l[j], znS[q][is], ONE, znS[j][is], znS[j][is]);
}
}
/*-----------------------------------------------------------------*/
/*
* CVRescale
*
* This routine rescales the Nordsieck array by multiplying the
* jth column zn[j] by eta^j, j = 1, ..., q. Then the value of
* h is rescaled by eta, and hscale is reset to h.
*/
static void CVRescale(CVodeMem cv_mem)
{
int j;
int is;
realtype factor;
factor = eta;
for (j=1; j <= q; j++) {
N_VScale(factor, zn[j], zn[j]);
if (quad)
N_VScale(factor, znQ[j], znQ[j]);
if (sensi)
for (is=0; is<Ns; is++)
N_VScale(factor, znS[j][is], znS[j][is]);
factor *= eta;
}
h = hscale * eta;
hscale = h;
nscon = 0;
}
/*-----------------------------------------------------------------*/
/*
* CVPredict
*
* This routine advances tn by the tentative step size h, and computes
* the predicted array z_n(0), which is overwritten on zn. The
* prediction of zn is done by repeated additions.
*/
static void CVPredict(CVodeMem cv_mem)
{
int j, k;
int is;
tn += h;
for (k = 1; k <= q; k++)
for (j = q; j >= k; j--)
N_VLinearSum(ONE, zn[j-1], ONE, zn[j], zn[j-1]);
if (quad) {
for (k = 1; k <= q; k++)
for (j = q; j >= k; j--)
N_VLinearSum(ONE, znQ[j-1], ONE, znQ[j], znQ[j-1]);
}
if (sensi) {
for (is=0; is<Ns; is++) {
for (k = 1; k <= q; k++)
for (j = q; j >= k; j--)
N_VLinearSum(ONE, znS[j-1][is], ONE, znS[j][is], znS[j-1][is]);
}
}
}
/*-----------------------------------------------------------------*/
/*
* CVSet
*
* This routine is a high level routine which calls CVSetAdams or
* CVSetBDF to set the polynomial l, the test quantity array tq,
* and the related variables rl1, gamma, and gamrat.
*/
static void CVSet(CVodeMem cv_mem)
{
switch(lmm) {
case CV_ADAMS:
CVSetAdams(cv_mem);
break;
case CV_BDF:
CVSetBDF(cv_mem);
break;
}
rl1 = ONE / l[1];
gamma = h * rl1;
if (nst == 0) gammap = gamma;
gamrat = (nst > 0) ? gamma / gammap : ONE; /* protect x / x != 1.0 */
}
/*-----------------------------------------------------------------*/
/*
* CVSetAdams
*
* This routine handles the computation of l and tq for the
* case lmm == CV_ADAMS.
*
* The components of the array l are the coefficients of a
* polynomial Lambda(x) = l_0 + l_1 x + ... + l_q x^q, given by
* q-1
* (d/dx) Lambda(x) = c * PRODUCT (1 + x / xi_i) , where
* i=1
* Lambda(-1) = 0, Lambda(0) = 1, and c is a normalization factor.
* Here xi_i = [t_n - t_(n-i)] / h.
*
* The array tq is set to test quantities used in the convergence
* test, the error test, and the selection of h at a new order.
*/
static void CVSetAdams(CVodeMem cv_mem)
{
realtype m[L_MAX], M[3], hsum;
if (q == 1) {
l[0] = l[1] = tq[1] = tq[5] = ONE;
tq[2] = TWO;
tq[3] = TWELVE;
tq[4] = nlscoef * tq[2]; /* = 0.1 * tq[2] */
return;
}
hsum = CVAdamsStart(cv_mem, m);
M[0] = CVAltSum(q-1, m, 1);
M[1] = CVAltSum(q-1, m, 2);
CVAdamsFinish(cv_mem, m, M, hsum);
}
/*-----------------------------------------------------------------*/
/*
* CVAdamsStart
*
* This routine generates in m[] the coefficients of the product
* polynomial needed for the Adams l and tq coefficients for q > 1.
*/
static realtype CVAdamsStart(CVodeMem cv_mem, realtype m[])
{
realtype hsum, xi_inv, sum;
int i, j;
hsum = h;
m[0] = ONE;
for (i=1; i <= q; i++) m[i] = ZERO;
for (j=1; j < q; j++) {
if ((j==q-1) && (qwait == 1)) {
sum = CVAltSum(q-2, m, 2);
tq[1] = m[q-2] / (q * sum);
}
xi_inv = h / hsum;
for (i=j; i >= 1; i--) m[i] += m[i-1] * xi_inv;
hsum += tau[j];
/* The m[i] are coefficients of product(1 to j) (1 + x/xi_i) */
}
return (hsum);
}
/*-----------------------------------------------------------------*/
/*
* CVAdamsFinish
*
* This routine completes the calculation of the Adams l and tq.
*/
static void CVAdamsFinish(CVodeMem cv_mem, realtype m[], realtype M[], realtype hsum)
{
int i;
realtype M0_inv, xi, xi_inv;
M0_inv = ONE / M[0];
l[0] = ONE;
for (i=1; i <= q; i++) l[i] = M0_inv * (m[i-1] / i);
xi = hsum / h;
xi_inv = ONE / xi;
tq[2] = xi * M[0] / M[1];
tq[5] = xi / l[q];
if (qwait == 1) {
for (i=q; i >= 1; i--) m[i] += m[i-1] * xi_inv;
M[2] = CVAltSum(q, m, 2);
tq[3] = L * M[0] / M[2];
}
tq[4] = nlscoef * tq[2];
}
/*-----------------------------------------------------------------*/
/*
* CVAltSum
*
* CVAltSum returns the value of the alternating sum
* sum (i= 0 ... iend) [ (-1)^i * (a[i] / (i + k)) ].
* If iend < 0 then CVAltSum returns 0.
* This operation is needed to compute the integral, from -1 to 0,
* of a polynomial x^(k-1) M(x) given the coefficients of M(x).
*/
static realtype CVAltSum(int iend, realtype a[], int k)
{
int i, sign;
realtype sum;
if (iend < 0) return (ZERO);
sum = ZERO;
sign = 1;
for (i=0; i <= iend; i++) {
sum += sign * (a[i] / (i+k));
sign = -sign;
}
return (sum);
}
/*-----------------------------------------------------------------*/
/*
* CVSetBDF
*
* This routine computes the coefficients l and tq in the case
* lmm == CV_BDF. CVSetBDF calls CVSetTqBDF to set the test
* quantity array tq.
*
* The components of the array l are the coefficients of a
* polynomial Lambda(x) = l_0 + l_1 x + ... + l_q x^q, given by
* q-1
* Lambda(x) = (1 + x / xi*_q) * PRODUCT (1 + x / xi_i) , where
* i=1
* xi_i = [t_n - t_(n-i)] / h.
*
* The array tq is set to test quantities used in the convergence
* test, the error test, and the selection of h at a new order.
*/
static void CVSetBDF(CVodeMem cv_mem)
{
realtype alpha0, alpha0_hat, xi_inv, xistar_inv, hsum;
int i,j;
l[0] = l[1] = xi_inv = xistar_inv = ONE;
for (i=2; i <= q; i++) l[i] = ZERO;
alpha0 = alpha0_hat = -ONE;
hsum = h;
if (q > 1) {
for (j=2; j < q; j++) {
hsum += tau[j-1];
xi_inv = h / hsum;
alpha0 -= ONE / j;
for (i=j; i >= 1; i--) l[i] += l[i-1]*xi_inv;
/* The l[i] are coefficients of product(1 to j) (1 + x/xi_i) */
}
/* j = q */
alpha0 -= ONE / q;
xistar_inv = -l[1] - alpha0;
hsum += tau[q-1];
xi_inv = h / hsum;
alpha0_hat = -l[1] - xi_inv;
for (i=q; i >= 1; i--) l[i] += l[i-1]*xistar_inv;
}
CVSetTqBDF(cv_mem, hsum, alpha0, alpha0_hat, xi_inv, xistar_inv);
}
/*-----------------------------------------------------------------*/
/*
* CVSetTqBDF
*
* This routine sets the test quantity array tq in the case
* lmm == CV_BDF.
*/
static void CVSetTqBDF(CVodeMem cv_mem, realtype hsum, realtype alpha0,
realtype alpha0_hat, realtype xi_inv, realtype xistar_inv)
{
realtype A1, A2, A3, A4, A5, A6;
realtype C, CPrime, CPrimePrime;
A1 = ONE - alpha0_hat + alpha0;
A2 = ONE + q * A1;
tq[2] = ABS(alpha0 * (A2 / A1));
tq[5] = ABS((A2) / (l[q] * xi_inv/xistar_inv));
if (qwait == 1) {
C = xistar_inv / l[q];
A3 = alpha0 + ONE / q;
A4 = alpha0_hat + xi_inv;
CPrime = A3 / (ONE - A4 + A3);
tq[1] = ABS(CPrime / C);
hsum += tau[q];
xi_inv = h / hsum;
A5 = alpha0 - (ONE / (q+1));
A6 = alpha0_hat - xi_inv;
CPrimePrime = A2 / (ONE - A6 + A5);
tq[3] = ABS(CPrimePrime * xi_inv * (q+2) * A5);
}
tq[4] = nlscoef * tq[2];
}
/*-----------------------------------------------------------------*/
/*
* CVnls
*
* This routine attempts to solve the nonlinear system associated
* with a single implicit step of the linear multistep method.
* Depending on iter, it calls CVNlsFunctional or CVNlsNewton
* to do the work.
*/
static int CVNls(CVodeMem cv_mem, int nflag)
{
int flag=SOLVED;
switch(iter) {
case CV_FUNCTIONAL:
flag = CVNlsFunctional(cv_mem);
break;
case CV_NEWTON:
flag = CVNlsNewton(cv_mem, nflag);
break;
}
return(flag);
}
/*-----------------------------------------------------------------*/
/*
* CVNlsFunctional
*
* This routine attempts to solve the nonlinear system using
* functional iteration (no matrices involved).
*
* This routine also handles the functional iteration of the
* combined system (states + sensitivities) when sensitivities are
* computed using the CV_SIMULTANEOUS approach.
*/
static int CVNlsFunctional(CVodeMem cv_mem)
{
int m;
realtype del, delS, Del, Delp, dcon;
int is;
booleantype do_sensi_sim;
N_Vector wrk1, wrk2;
/* Are we computing sensitivities with the CV_SIMULTANEOUS approach? */
do_sensi_sim = (sensi && (ism==CV_SIMULTANEOUS));
/* Initialize counter and evaluate f at predicted y */
crate = ONE;
m = 0;
f(tn, zn[0], tempv, f_data);
nfe++;
if (do_sensi_sim) {
wrk1 = ftemp;
wrk2 = ftempS[0];
CVSensRhs(cv_mem, tn, zn[0], tempv, znS[0], tempvS, wrk1, wrk2);
}
/* Initialize correction to zero */
N_VConst(ZERO, acor);
if (do_sensi_sim) {
for (is=0; is<Ns; is++)
N_VConst(ZERO,acorS[is]);
}
/* Loop until convergence; accumulate corrections in acor */
loop {
nni++;
/* Correct y directly from the last f value */
N_VLinearSum(h, tempv, -ONE, zn[1], tempv);
N_VScale(rl1, tempv, tempv);
N_VLinearSum(ONE, zn[0], ONE, tempv, y);
if (do_sensi_sim)
for (is=0; is<Ns; is++) {
N_VLinearSum(h, tempvS[is], -ONE, znS[1][is], tempvS[is]);
N_VScale(rl1, tempvS[is], tempvS[is]);
N_VLinearSum(ONE, znS[0][is], ONE, tempvS[is], yS[is]);
}
/* Get WRMS norm of current correction to use in convergence test */
N_VLinearSum(ONE, tempv, -ONE, acor, acor);
if (do_sensi_sim)
for (is=0; is<Ns; is++)
N_VLinearSum(ONE, tempvS[is], -ONE, acorS[is], acorS[is]);
del = N_VWrmsNorm(acor, ewt);
if (do_sensi_sim)
delS = CVSensUpdateNorm(cv_mem, del, acorS, ewtS);
N_VScale(ONE, tempv, acor);
if (do_sensi_sim)
for (is=0; is<Ns; is++)
N_VScale(ONE, tempvS[is], acorS[is]);
/* Test for convergence. If m > 0, an estimate of the convergence
rate constant is stored in crate, and used in the test.
Recall that, even when errconS=FALSE, all variables are used in the
convergence test. Hence, we use Del (and not del). However, acnrm
is used in the error test and thus it has different forms
depending on errconS (and this explains why we have to carry around
del and delS)
*/
Del = (do_sensi_sim) ? delS : del;
if (m > 0) crate = MAX(CRDOWN * crate, Del / Delp);
dcon = Del * MIN(ONE, crate) / tq[4];
if (dcon <= ONE) {
if (m == 0)
if (do_sensi_sim && errconS) acnrm = delS;
else acnrm = del;
else {
acnrm = N_VWrmsNorm(acor, ewt);
if (do_sensi_sim && errconS)
acnrm = CVSensUpdateNorm(cv_mem, acnrm, acorS, ewtS);
}
return (SOLVED); /* Convergence achieved */
}
/* Stop at maxcor iterations or if iter. seems to be diverging */
m++;
if ((m==maxcor) || ((m >= 2) && (Del > RDIV * Delp)))
return (CONV_FAIL);
/* Save norm of correction, evaluate f, and loop again */
Delp = Del;
f(tn, y, tempv, f_data);
nfe++;
if (do_sensi_sim) {
wrk1 = ftemp;
wrk2 = ftempS[0];
CVSensRhs(cv_mem, tn, y, tempv, yS, tempvS, wrk1, wrk2);
}
} /* end loop */
}
/*-----------------------------------------------------------------*/
/*
* CVNlsNewton
*
* This routine handles the Newton iteration. It calls lsetup if
* indicated, calls CVNewtonIteration to perform the iteration, and
* retries a failed attempt at Newton iteration if that is indicated.
* See return values at top of this file.
*
* This routine also handles the Newton iteration of the combined
* system when sensitivities are computed using the CV_SIMULTANEOUS
* approach. Since in that case we use a quasi-Newton on the
* combined system (by approximating the Jacobian matrix by its
* block diagonal) and thus only solve linear systems with
* multiple right hand sides (all sharing the same coefficient
* matrix - whatever iteration matrix we decide on) we set-up
* the linear solver to handle N equations at a time.
*/
static int CVNlsNewton(CVodeMem cv_mem, int nflag)
{
N_Vector vtemp1, vtemp2, vtemp3, wrk1, wrk2;
int convfail, ier;
booleantype callSetup, do_sensi_sim;
int is;
/* Are we computing sensitivities with the CV_SIMULTANEOUS approach? */
do_sensi_sim = (sensi && (ism==CV_SIMULTANEOUS));
vtemp1 = acor; /* rename acor as vtemp1 for readability */
vtemp2 = y; /* rename y as vtemp2 for readability */
vtemp3 = tempv; /* rename tempv as vtemp3 for readability */
/* Set flag convfail, input to lsetup for its evaluation decision */
convfail = ((nflag == FIRST_CALL) || (nflag == PREV_ERR_FAIL)) ?
CV_NO_FAILURES : CV_FAIL_OTHER;
/* Decide whether or not to call setup routine (if one exists) */
if (setupNonNull) {
callSetup = (nflag == PREV_CONV_FAIL) || (nflag == PREV_ERR_FAIL) ||
(nst == 0) || (nst >= nstlp + MSBP) || (ABS(gamrat-ONE) > DGMAX);
/* Decide whether to force a call to setup */
if (forceSetup) {
callSetup = TRUE;
convfail = CV_FAIL_OTHER;
}
} else {
crate = ONE;
crateS = ONE; /* if NO lsetup all conv. rates are set to ONE */
callSetup = FALSE;
}
/* Looping point for the solution of the nonlinear system.
Evaluate f at the predicted y, call lsetup if indicated, and
call CVNewtonIteration for the Newton iteration itself. */
loop {
f(tn, zn[0], ftemp, f_data);
nfe++;
if (do_sensi_sim) {
wrk1 = tempv;
wrk2 = tempvS[0];
CVSensRhs(cv_mem, tn, zn[0], ftemp, znS[0], ftempS, wrk1, wrk2);
}
if (callSetup) {
ier = lsetup(cv_mem, convfail, zn[0], ftemp, &jcur,
vtemp1, vtemp2, vtemp3);
nsetups++;
callSetup = FALSE;
forceSetup = FALSE;
gamrat = ONE;
gammap = gamma;
crate = ONE;
crateS = ONE; /* after lsetup all conv. rates are reset to ONE */
nstlp = nst;
/* Return if lsetup failed */
if (ier < 0) return (SETUP_FAIL_UNREC);
if (ier > 0) return (CONV_FAIL);
}
/* Set acor to zero and load prediction into y vector */
N_VConst(ZERO, acor);
N_VScale(ONE, zn[0], y);
if (do_sensi_sim)
for (is=0; is<Ns; is++) {
N_VConst(ZERO, acorS[is]);
N_VScale(ONE, znS[0][is], yS[is]);
}
/* Do the Newton iteration */
ier = CVNewtonIteration(cv_mem);
/* If there is a convergence failure and the Jacobian-related
data appears not to be current, loop again with a call to lsetup
in which convfail=CV_FAIL_BAD_J. Otherwise return. */
if (ier != TRY_AGAIN) return (ier);
callSetup = TRUE;
convfail = CV_FAIL_BAD_J;
}
}
/*-----------------------------------------------------------------*/
/*
* CVNewtonIteration
*
* This routine performs the Newton iteration. If the iteration succeeds,
* it returns the value SOLVED. If not, it may signal the CVNlsNewton
* routine to call lsetup again and reattempt the iteration, by
* returning the value TRY_AGAIN. (In this case, CVNlsNewton must set
* convfail to CV_FAIL_BAD_J before calling setup again).
* Otherwise, this routine returns one of the appropriate values
* SOLVE_FAIL_UNREC or CONV_FAIL back to CVNlsNewton.
*
* If sensitivities are computed using the CV_SIMULTANEOUS approach, this
* routine performs a quasi-Newton on the combined nonlinear system.
* The iteration matrix of the combined system is block diagonal with
* each block being the iteration matrix of the original system. Thus
* we solve linear systems with the same matrix but different right
* hand sides.
*/
static int CVNewtonIteration(CVodeMem cv_mem)
{
int m, ret;
realtype del, delS, Del, Delp, dcon;
N_Vector b, *bS=NULL, wrk1, wrk2;
booleantype do_sensi_sim;
int is;
/* Are we computing sensitivities with the CV_SIMULTANEOUS approach? */
do_sensi_sim = (sensi && (ism==CV_SIMULTANEOUS));
mnewt = m = 0;
/* At this point, ftemp <- f(t_n, y_n(0))
ftempS <- fS(t_n, y_n(0), s_n(0))
acor <- 0
acorS <- 0
y <- y_n(0)
yS <- yS_n(0) */
/* Looping point for Newton iteration */
loop {
/* Evaluate the residual of the nonlinear system */
N_VLinearSum(rl1, zn[1], ONE, acor, tempv);
N_VLinearSum(gamma, ftemp, -ONE, tempv, tempv);
/* Call the lsolve function */
b = tempv;
ret = lsolve(cv_mem, b, ewt, y, ftemp);
nni++;
if (ret < 0) return (SOLVE_FAIL_UNREC);
/* If lsolve had a recoverable failure and Jacobian data is
not current, signal to try the solution again */
if (ret > 0) {
if ((!jcur) && (setupNonNull)) return (TRY_AGAIN);
return (CONV_FAIL);
}
/* Solve the sensitivity linear systems and do the same
tests on the return value of lsolve. */
if (do_sensi_sim) {
for (is=0; is<Ns; is++) {
N_VLinearSum(rl1, znS[1][is], ONE, acorS[is], tempvS[is]);
N_VLinearSum(gamma, ftempS[is], -ONE, tempvS[is], tempvS[is]);
}
bS = tempvS;
for (is=0; is<Ns; is++) {
ret = lsolve(cv_mem, bS[is], ewtS[is], y, ftemp);
if (ret < 0) return (SOLVE_FAIL_UNREC);
if (ret > 0) {
if ((!jcur) && (setupNonNull)) return (TRY_AGAIN);
return (CONV_FAIL);
}
}
}
/* Get WRMS norm of correction; add correction to acor and y */
del = N_VWrmsNorm(b, ewt);
N_VLinearSum(ONE, acor, ONE, b, acor);
N_VLinearSum(ONE, zn[0], ONE, acor, y);
if (do_sensi_sim) {
delS = CVSensUpdateNorm(cv_mem, del, bS, ewtS);
for (is=0; is<Ns; is++) {
N_VLinearSum(ONE, acorS[is], ONE, bS[is], acorS[is]);
N_VLinearSum(ONE, znS[0][is], ONE, acorS[is], yS[is]);
}
}
/* Test for convergence. If m > 0, an estimate of the convergence
rate constant is stored in crate, and used in the test. */
Del = (do_sensi_sim) ? delS : del;
if (m > 0) crate = MAX(CRDOWN * crate, Del/Delp);
dcon = Del * MIN(ONE, crate) / tq[4];
if (dcon <= ONE) {
if (m == 0)
if (do_sensi_sim && errconS) acnrm = delS;
else acnrm = del;
else {
acnrm = N_VWrmsNorm(acor, ewt);
if (do_sensi_sim && errconS)
acnrm = CVSensUpdateNorm(cv_mem, acnrm, acorS, ewtS);
}
jcur = FALSE;
return (SOLVED); /* Convergence achieved */
}
mnewt = ++m;
/* Stop at maxcor iterations or if iter. seems to be diverging.
If still not converged and Jacobian data is not current,
signal to try the solution again */
if ((m == maxcor) || ((m >= 2) && (Del > RDIV * Delp))) {
if ((!jcur) && (setupNonNull)) return (TRY_AGAIN);
return (CONV_FAIL);
}
/* Save norm of correction, evaluate f, and loop again */
Delp = Del;
f(tn, y, ftemp, f_data);
nfe++;
if (do_sensi_sim) {
wrk1 = tempv;
wrk2 = tempvS[0];
CVSensRhs(cv_mem, tn, y, ftemp, yS, ftempS, wrk1, wrk2);
}
} /* end loop */
}
/*-----------------------------------------------------------------*/
/*
* CVHandleFlag
*
* This routine takes action on the return value nflag = *nflagPtr
* returned by CVNls, as follows:
*
* If CVNls succeeded in solving the nonlinear system, then
* CVHandleNFlag returns the constant DO_ERROR_TEST, which tells CVStep
* to perform the error test.
*
* If the nonlinear system was not solved successfully, then ncfn and
* ncf = *ncfPtr are incremented and Nordsieck array zn is restored.
*
* If the solution of the nonlinear system failed due to an
* unrecoverable failure by setup, we return the value SETUP_FAILED.
*
* If it failed due to an unrecoverable failure in solve, then we return
* the value SOLVE_FAILED.
*
* Otherwise, a recoverable failure occurred when solving the
* nonlinear system (CVNls returned nflag == CONV_FAIL).
* In this case, we return the value REP_CONV_FAIL if ncf is now
* equal to maxncf or |h| = hmin.
* If not, we set *nflagPtr = PREV_CONV_FAIL and return the value
* PREDICT_AGAIN, telling CVStep to reattempt the step.
*/
static int CVHandleNFlag(CVodeMem cv_mem, int *nflagPtr, realtype saved_t,
int *ncfPtr, long int *ncfnPtr)
{
int nflag;
nflag = *nflagPtr;
if (nflag == SOLVED) return (DO_ERROR_TEST);
/* The nonlinear soln. failed; increment ncfn and restore zn */
(*ncfnPtr)++;
CVRestore(cv_mem, saved_t);
/* Return if lsetup or lsolve failed unrecoverably */
if (nflag == SETUP_FAIL_UNREC) return (SETUP_FAILED);
if (nflag == SOLVE_FAIL_UNREC) return (SOLVE_FAILED);
/* At this point, nflag == CONV_FAIL; increment ncf */
(*ncfPtr)++;
etamax = ONE;
/* If we had maxncf failures or |h| = hmin, return REP_CONV_FAIL */
if ((ABS(h) <= hmin*ONEPSM) || (*ncfPtr == maxncf))
return (REP_CONV_FAIL);
/* Reduce step size; return to reattempt the step */
eta = MAX(ETACF, hmin / ABS(h));
*nflagPtr = PREV_CONV_FAIL;
CVRescale(cv_mem);
return (PREDICT_AGAIN);
}
/*-----------------------------------------------------------------*/
/*
* CVRestore
*
* This routine restores the value of tn to saved_t and undoes the
* prediction. After execution of CVRestore, the Nordsieck array zn has
* the same values as before the call to CVPredict.
*/
static void CVRestore(CVodeMem cv_mem, realtype saved_t)
{
int j, k;
int is;
tn = saved_t;
for (k = 1; k <= q; k++)
for (j = q; j >= k; j--)
N_VLinearSum(ONE, zn[j-1], -ONE, zn[j], zn[j-1]);
if (quad) {
for (k = 1; k <= q; k++)
for (j = q; j >= k; j--)
N_VLinearSum(ONE, znQ[j-1], -ONE, znQ[j], znQ[j-1]);
}
if (sensi) {
for (is=0; is<Ns; is++) {
for (k = 1; k <= q; k++)
for (j = q; j >= k; j--)
N_VLinearSum(ONE, znS[j-1][is], -ONE, znS[j][is], znS[j-1][is]);
}
}
}
/*-----------------------------------------------------------------*/
/*
* CVDoErrorTest
*
* This routine performs the local error test.
* The weighted local error norm dsm is loaded into *dsmPtr, and
* the test dsm ?<= 1 is made.
*
* If the test passes, CVDoErrorTest returns TRUE.
*
* If the test fails, we undo the step just taken (call CVRestore),
* set *nflagPtr to PREV_ERR_FAIL, and return FALSE.
*
* If maxnef error test failures have occurred or if ABS(h) = hmin,
* we set *kflagPtr = REP_ERR_FAIL. (Otherwise *kflagPtr has the
* value last returned by CVHandleNFlag.)
*
* If more than MXNEF1 error test failures have occurred, an order
* reduction is forced. If already at order 1 restart by reloading
* zn from scratch. Note that if sensitivities are computed, znS is
* also reloaded, no matter what 'ism' or 'errconS' are. Same for
* quadratures.
*/
static booleantype CVDoErrorTest(CVodeMem cv_mem, int *nflagPtr,
int *kflagPtr, realtype saved_t,
int *nefPtr, realtype *dsmPtr)
{
realtype dsm;
int is;
N_Vector wrk1, wrk2;
dsm = acnrm / tq[2];
/* If est. local error norm dsm passes test, return TRUE */
*dsmPtr = dsm;
if (dsm <= ONE) return (TRUE);
/* Test failed; increment counters, set nflag, and restore zn array */
(*nefPtr)++;
netf++;
*nflagPtr = PREV_ERR_FAIL;
CVRestore(cv_mem, saved_t);
/* At maxnef failures or |h| = hmin, return with kflag = REP_ERR_FAIL */
if ((ABS(h) <= hmin*ONEPSM) || (*nefPtr == maxnef)) {
*kflagPtr = REP_ERR_FAIL;
return (FALSE);
}
/* Set etamax = 1 to prevent step size increase at end of this step */
etamax = ONE;
/* Set h ratio eta from dsm, rescale, and return for retry of step */
if (*nefPtr <= MXNEF1) {
eta = ONE / (RPowerR(BIAS2*dsm,ONE/L) + ADDON);
eta = MAX(ETAMIN, MAX(eta, hmin / ABS(h)));
if (*nefPtr >= SMALL_NEF) eta = MIN(eta, ETAMXF);
CVRescale(cv_mem);
return (FALSE);
}
/* After MXNEF1 failures, force an order reduction and retry step */
if (q > 1) {
eta = MAX(ETAMIN, hmin / ABS(h));
CVAdjustOrder(cv_mem,-1);
L = q;
q--;
qwait = L;
CVRescale(cv_mem);
return (FALSE);
}
/* If already at order 1, restart: reload zn from scratch */
eta = MAX(ETAMIN, hmin / ABS(h));
h *= eta;
hscale = h;
qwait = LONG_WAIT;
nscon = 0;
f(tn, zn[0], tempv, f_data);
nfe++;
N_VScale(h, tempv, zn[1]);
if (quad) {
fQ(tn, zn[0], tempvQ, fQ_data);
nfQe++;
N_VScale(h, tempvQ, znQ[1]);
}
if (sensi) {
wrk1 = ftemp;
wrk2 = ftempS[0];
CVSensRhs(cv_mem, tn, zn[0], tempv, znS[0], tempvS, wrk1, wrk2);
for (is=0; is<Ns; is++)
N_VScale(h, tempvS[is], znS[1][is]);
}
return (FALSE);
}
/*-----------------------------------------------------------------*/
static booleantype CVQuadDoErrorTest(CVodeMem cv_mem, int *nflagPtr,
int *kflagPtr, realtype saved_t,
int *nefQPtr, realtype *dsmQPtr)
{
realtype dsmQ;
int is;
N_Vector wrk1, wrk2;
dsmQ = acnrmQ / tq[2];
/* If dsmQ passes the test, return TRUE */
*dsmQPtr = dsmQ;
if (dsmQ <= ONE) return (TRUE);
/* Test failed; increment counters, set nflag, and restore zn array */
(*nefQPtr)++;
netfQ++;
*nflagPtr = PREV_ERR_FAIL;
CVRestore(cv_mem, saved_t);
/* At maxnef failures or |h| = hmin, return with kflag = REP_ERR_FAIL */
if ((ABS(h) <= hmin*ONEPSM) || (*nefQPtr == maxnef)) {
*kflagPtr = REP_ERR_FAIL;
return (FALSE);
}
/* Set etamax = 1 to prevent step size increase at end of this step */
etamax = ONE;
/* Set h ratio eta from dsmQ, rescale, and return for retry of step */
if (*nefQPtr <= MXNEF1) {
eta = ONE / (RPowerR(BIAS2*dsmQ,ONE/L) + ADDON);
eta = MAX(ETAMIN, MAX(eta, hmin / ABS(h)));
if (*nefQPtr >= SMALL_NEF) eta = MIN(eta, ETAMXF);
CVRescale(cv_mem);
return (FALSE);
}
/* After MXNEF1 failures, force an order reduction and retry step */
if (q > 1) {
eta = MAX(ETAMIN, hmin / ABS(h));
CVAdjustOrder(cv_mem,-1);
L = q;
q--;
qwait = L;
CVRescale(cv_mem);
return (FALSE);
}
/* If already at order 1, restart: reload zn and znQ from scratch */
eta = MAX(ETAMIN, hmin / ABS(h));
h *= eta;
hscale = h;
qwait = LONG_WAIT;
nscon = 0;
f(tn, zn[0], tempv, f_data);
nfe++;
N_VScale(h, tempv, zn[1]);