xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
17cb94ee6SHong Zhang /*
2db268bfaSHong Zhang     Provides a PETSc interface to SUNDIALS/CVODE solver.
37cb94ee6SHong Zhang     The interface to PVODE (old version of CVODE) was originally contributed
47cb94ee6SHong Zhang     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.
528adb3f7SHong Zhang 
628adb3f7SHong Zhang     Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c
77cb94ee6SHong Zhang */
81c7d2463SJed Brown #include <../src/ts/impls/implicit/sundials/sundials.h> /*I "petscts.h" I*/
97cb94ee6SHong Zhang 
107cb94ee6SHong Zhang /*
117cb94ee6SHong Zhang       TSPrecond_Sundials - function that we provide to SUNDIALS to
127cb94ee6SHong Zhang                         evaluate the preconditioner.
137cb94ee6SHong Zhang */
14*9371c9d4SSatish Balay PetscErrorCode TSPrecond_Sundials(realtype tn, N_Vector y, N_Vector fy, booleantype jok, booleantype *jcurPtr, realtype _gamma, void *P_data, N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3) {
157cb94ee6SHong Zhang   TS           ts    = (TS)P_data;
167cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
17dcbc6d53SSean Farley   PC           pc;
180679a0aeSJed Brown   Mat          J, P;
190679a0aeSJed Brown   Vec          yy = cvode->w1, yydot = cvode->ydot;
200679a0aeSJed Brown   PetscReal    gm = (PetscReal)_gamma;
217cb94ee6SHong Zhang   PetscScalar *y_data;
227cb94ee6SHong Zhang 
237cb94ee6SHong Zhang   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &P, NULL, NULL));
257cb94ee6SHong Zhang   y_data = (PetscScalar *)N_VGetArrayPointer(y);
269566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(yy, y_data));
279566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(yydot)); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
280679a0aeSJed Brown   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
299566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ts->ptime, yy, yydot, 1 / gm, J, P, PETSC_FALSE));
309566063dSJacob Faibussowitsch   PetscCall(VecResetArray(yy));
319566063dSJacob Faibussowitsch   PetscCall(MatScale(P, gm)); /* turn into I-gm*Jrest, J is not used by Sundials  */
327cb94ee6SHong Zhang   *jcurPtr = TRUE;
339566063dSJacob Faibussowitsch   PetscCall(TSSundialsGetPC(ts, &pc));
349566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc, J, P));
357cb94ee6SHong Zhang   PetscFunctionReturn(0);
367cb94ee6SHong Zhang }
377cb94ee6SHong Zhang 
387cb94ee6SHong Zhang /*
397cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
407cb94ee6SHong Zhang */
41*9371c9d4SSatish Balay PetscErrorCode TSPSolve_Sundials(realtype tn, N_Vector y, N_Vector fy, N_Vector r, N_Vector z, realtype _gamma, realtype delta, int lr, void *P_data, N_Vector vtemp) {
427cb94ee6SHong Zhang   TS           ts    = (TS)P_data;
437cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
44dcbc6d53SSean Farley   PC           pc;
457cb94ee6SHong Zhang   Vec          rr = cvode->w1, zz = cvode->w2;
467cb94ee6SHong Zhang   PetscScalar *r_data, *z_data;
477cb94ee6SHong Zhang 
487cb94ee6SHong Zhang   PetscFunctionBegin;
497cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
507cb94ee6SHong Zhang   r_data = (PetscScalar *)N_VGetArrayPointer(r);
517cb94ee6SHong Zhang   z_data = (PetscScalar *)N_VGetArrayPointer(z);
529566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(rr, r_data));
539566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(zz, z_data));
544ac7836bSHong Zhang 
557cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
569566063dSJacob Faibussowitsch   PetscCall(TSSundialsGetPC(ts, &pc));
579566063dSJacob Faibussowitsch   PetscCall(PCApply(pc, rr, zz));
589566063dSJacob Faibussowitsch   PetscCall(VecResetArray(rr));
599566063dSJacob Faibussowitsch   PetscCall(VecResetArray(zz));
607cb94ee6SHong Zhang   PetscFunctionReturn(0);
617cb94ee6SHong Zhang }
627cb94ee6SHong Zhang 
637cb94ee6SHong Zhang /*
647cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
657cb94ee6SHong Zhang */
66*9371c9d4SSatish Balay int TSFunction_Sundials(realtype t, N_Vector y, N_Vector ydot, void *ctx) {
677cb94ee6SHong Zhang   TS           ts = (TS)ctx;
685fcef5e4SJed Brown   DM           dm;
69942e3340SBarry Smith   DMTS         tsdm;
705fcef5e4SJed Brown   TSIFunction  ifunction;
71ce94432eSBarry Smith   MPI_Comm     comm;
727cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
730679a0aeSJed Brown   Vec          yy = cvode->w1, yyd = cvode->w2, yydot = cvode->ydot;
747cb94ee6SHong Zhang   PetscScalar *y_data, *ydot_data;
757cb94ee6SHong Zhang 
767cb94ee6SHong Zhang   PetscFunctionBegin;
779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
785ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
797cb94ee6SHong Zhang   y_data    = (PetscScalar *)N_VGetArrayPointer(y);
807cb94ee6SHong Zhang   ydot_data = (PetscScalar *)N_VGetArrayPointer(ydot);
819566063dSJacob Faibussowitsch   PetscCallAbort(comm, VecPlaceArray(yy, y_data));
829566063dSJacob Faibussowitsch   PetscCallAbort(comm, VecPlaceArray(yyd, ydot_data));
834ac7836bSHong Zhang 
845fcef5e4SJed Brown   /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
859566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
869566063dSJacob Faibussowitsch   PetscCall(DMGetDMTS(dm, &tsdm));
879566063dSJacob Faibussowitsch   PetscCall(DMTSGetIFunction(dm, &ifunction, NULL));
885fcef5e4SJed Brown   if (!ifunction) {
899566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(ts, t, yy, yyd));
90bc0cc02bSJed Brown   } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */
919566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(yydot));
929566063dSJacob Faibussowitsch     PetscCallAbort(comm, TSComputeIFunction(ts, t, yy, yydot, yyd, PETSC_FALSE));
939566063dSJacob Faibussowitsch     PetscCall(VecScale(yyd, -1.));
94bc0cc02bSJed Brown   }
959566063dSJacob Faibussowitsch   PetscCallAbort(comm, VecResetArray(yy));
969566063dSJacob Faibussowitsch   PetscCallAbort(comm, VecResetArray(yyd));
974ac7836bSHong Zhang   PetscFunctionReturn(0);
987cb94ee6SHong Zhang }
997cb94ee6SHong Zhang 
1007cb94ee6SHong Zhang /*
101b4eba00bSSean Farley        TSStep_Sundials - Calls Sundials to integrate the ODE.
1027cb94ee6SHong Zhang */
103*9371c9d4SSatish Balay PetscErrorCode TSStep_Sundials(TS ts) {
1047cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
105b4eba00bSSean Farley   PetscInt     flag;
106be5899b3SLisandro Dalcin   long int     nits, lits, nsteps;
1077cb94ee6SHong Zhang   realtype     t, tout;
1087cb94ee6SHong Zhang   PetscScalar *y_data;
1097cb94ee6SHong Zhang   void        *mem;
1107cb94ee6SHong Zhang 
1117cb94ee6SHong Zhang   PetscFunctionBegin;
11216016685SHong Zhang   mem  = cvode->mem;
1137cb94ee6SHong Zhang   tout = ts->max_time;
1149566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ts->vec_sol, &y_data));
1157cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data, cvode->y);
1169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ts->vec_sol, NULL));
117186e87acSLisandro Dalcin 
1189687d888SLisandro Dalcin   /* We would like to TSPreStage() and TSPostStage()
1199687d888SLisandro Dalcin    * before each stage solve but CVode does not appear to support this. */
120*9371c9d4SSatish Balay   if (cvode->monitorstep) flag = CVode(mem, tout, cvode->y, &t, CV_ONE_STEP);
121*9371c9d4SSatish Balay   else flag = CVode(mem, tout, cvode->y, &t, CV_NORMAL);
1229f94935aSHong Zhang 
1239f94935aSHong Zhang   if (flag) { /* display error message */
1249f94935aSHong Zhang     switch (flag) {
125*9371c9d4SSatish Balay     case CV_ILL_INPUT: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ILL_INPUT"); break;
126*9371c9d4SSatish Balay     case CV_TOO_CLOSE: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_CLOSE"); break;
1273c7fefeeSJed Brown     case CV_TOO_MUCH_WORK: {
1289f94935aSHong Zhang       PetscReal tcur;
1299566063dSJacob Faibussowitsch       PetscCall(CVodeGetNumSteps(mem, &nsteps));
1309566063dSJacob Faibussowitsch       PetscCall(CVodeGetCurrentTime(mem, &tcur));
13163a3b9bcSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_MUCH_WORK. At t=%g, nsteps %ld exceeds maxstep %" PetscInt_FMT ". Increase '-ts_max_steps <>' or modify TSSetMaxSteps()", (double)tcur, nsteps, ts->max_steps);
1323c7fefeeSJed Brown     } break;
133*9371c9d4SSatish Balay     case CV_TOO_MUCH_ACC: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_MUCH_ACC"); break;
134*9371c9d4SSatish Balay     case CV_ERR_FAILURE: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ERR_FAILURE"); break;
135*9371c9d4SSatish Balay     case CV_CONV_FAILURE: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_CONV_FAILURE"); break;
136*9371c9d4SSatish Balay     case CV_LINIT_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LINIT_FAIL"); break;
137*9371c9d4SSatish Balay     case CV_LSETUP_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSETUP_FAIL"); break;
138*9371c9d4SSatish Balay     case CV_LSOLVE_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSOLVE_FAIL"); break;
139*9371c9d4SSatish Balay     case CV_RHSFUNC_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RHSFUNC_FAIL"); break;
140*9371c9d4SSatish Balay     case CV_FIRST_RHSFUNC_ERR: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_FIRST_RHSFUNC_ERR"); break;
141*9371c9d4SSatish Balay     case CV_REPTD_RHSFUNC_ERR: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_REPTD_RHSFUNC_ERR"); break;
142*9371c9d4SSatish Balay     case CV_UNREC_RHSFUNC_ERR: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_UNREC_RHSFUNC_ERR"); break;
143*9371c9d4SSatish Balay     case CV_RTFUNC_FAIL: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RTFUNC_FAIL"); break;
144*9371c9d4SSatish Balay     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, flag %d", flag);
1459f94935aSHong Zhang     }
1469f94935aSHong Zhang   }
1479f94935aSHong Zhang 
148be5899b3SLisandro Dalcin   /* log inner nonlinear and linear iterations */
1499566063dSJacob Faibussowitsch   PetscCall(CVodeGetNumNonlinSolvIters(mem, &nits));
1509566063dSJacob Faibussowitsch   PetscCall(CVSpilsGetNumLinIters(mem, &lits));
151*9371c9d4SSatish Balay   ts->snes_its += nits;
152*9371c9d4SSatish Balay   ts->ksp_its = lits;
153be5899b3SLisandro Dalcin 
1547cb94ee6SHong Zhang   /* copy the solution from cvode->y to cvode->update and sol */
1559566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(cvode->w1, y_data));
1569566063dSJacob Faibussowitsch   PetscCall(VecCopy(cvode->w1, cvode->update));
1579566063dSJacob Faibussowitsch   PetscCall(VecResetArray(cvode->w1));
1589566063dSJacob Faibussowitsch   PetscCall(VecCopy(cvode->update, ts->vec_sol));
159186e87acSLisandro Dalcin 
160186e87acSLisandro Dalcin   ts->time_step = t - ts->ptime;
161186e87acSLisandro Dalcin   ts->ptime     = t;
162186e87acSLisandro Dalcin 
1639566063dSJacob Faibussowitsch   PetscCall(CVodeGetNumSteps(mem, &nsteps));
1642808aa04SLisandro Dalcin   if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */
165b4eba00bSSean Farley   PetscFunctionReturn(0);
166b4eba00bSSean Farley }
167b4eba00bSSean Farley 
168*9371c9d4SSatish Balay static PetscErrorCode TSInterpolate_Sundials(TS ts, PetscReal t, Vec X) {
169b4eba00bSSean Farley   TS_Sundials *cvode = (TS_Sundials *)ts->data;
170b4eba00bSSean Farley   N_Vector     y;
171b4eba00bSSean Farley   PetscScalar *x_data;
172b4eba00bSSean Farley   PetscInt     glosize, locsize;
173b4eba00bSSean Farley 
174b4eba00bSSean Farley   PetscFunctionBegin;
175b4eba00bSSean Farley   /* get the vector size */
1769566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &glosize));
1779566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &locsize));
1789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x_data));
179b4eba00bSSean Farley 
1804d78c9e0SClaas Abert   /* Initialize N_Vec y with x_data */
181918687b8SPatrick Sanan   if (cvode->use_dense) {
182918687b8SPatrick Sanan     PetscMPIInt size;
183918687b8SPatrick Sanan 
1849566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1853c633725SBarry Smith     PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case");
186918687b8SPatrick Sanan     y = N_VMake_Serial(locsize, (realtype *)x_data);
187918687b8SPatrick Sanan   } else {
1884d78c9e0SClaas Abert     y = N_VMake_Parallel(cvode->comm_sundials, locsize, glosize, (realtype *)x_data);
189918687b8SPatrick Sanan   }
190918687b8SPatrick Sanan 
1913c633725SBarry Smith   PetscCheck(y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Interpolated y is not allocated");
192b4eba00bSSean Farley 
1939566063dSJacob Faibussowitsch   PetscCall(CVodeGetDky(cvode->mem, t, 0, y));
1949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x_data));
1957cb94ee6SHong Zhang   PetscFunctionReturn(0);
1967cb94ee6SHong Zhang }
1977cb94ee6SHong Zhang 
198*9371c9d4SSatish Balay PetscErrorCode TSReset_Sundials(TS ts) {
199277b19d0SLisandro Dalcin   TS_Sundials *cvode = (TS_Sundials *)ts->data;
200277b19d0SLisandro Dalcin 
201277b19d0SLisandro Dalcin   PetscFunctionBegin;
2029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cvode->update));
2039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cvode->ydot));
2049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cvode->w1));
2059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cvode->w2));
206bbd56ea5SKarl Rupp   if (cvode->mem) CVodeFree(&cvode->mem);
207277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
208277b19d0SLisandro Dalcin }
209277b19d0SLisandro Dalcin 
210*9371c9d4SSatish Balay PetscErrorCode TSDestroy_Sundials(TS ts) {
2117cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
2127cb94ee6SHong Zhang 
2137cb94ee6SHong Zhang   PetscFunctionBegin;
2149566063dSJacob Faibussowitsch   PetscCall(TSReset_Sundials(ts));
2159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free(&(cvode->comm_sundials)));
2169566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
2179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", NULL));
2189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", NULL));
2199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", NULL));
2219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", NULL));
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", NULL));
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", NULL));
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", NULL));
2269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", NULL));
2272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", NULL));
2287cb94ee6SHong Zhang   PetscFunctionReturn(0);
2297cb94ee6SHong Zhang }
2307cb94ee6SHong Zhang 
231*9371c9d4SSatish Balay PetscErrorCode TSSetUp_Sundials(TS ts) {
2327cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
23316016685SHong Zhang   PetscInt     glosize, locsize, i, flag;
2347cb94ee6SHong Zhang   PetscScalar *y_data, *parray;
23516016685SHong Zhang   void        *mem;
236dcbc6d53SSean Farley   PC           pc;
23719fd82e9SBarry Smith   PCType       pctype;
238ace3abfcSBarry Smith   PetscBool    pcnone;
2397cb94ee6SHong Zhang 
2407cb94ee6SHong Zhang   PetscFunctionBegin;
24108401ef6SPierre Jolivet   PetscCheck(ts->exact_final_time != TS_EXACTFINALTIME_MATCHSTEP, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for exact final time option 'MATCHSTEP' when using Sundials");
242236c45dcSShri Abhyankar 
2437cb94ee6SHong Zhang   /* get the vector size */
2449566063dSJacob Faibussowitsch   PetscCall(VecGetSize(ts->vec_sol, &glosize));
2459566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(ts->vec_sol, &locsize));
2467cb94ee6SHong Zhang 
2477cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
248918687b8SPatrick Sanan   if (cvode->use_dense) {
249918687b8SPatrick Sanan     PetscMPIInt size;
250918687b8SPatrick Sanan 
2519566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2523c633725SBarry Smith     PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case");
253918687b8SPatrick Sanan     cvode->y = N_VNew_Serial(locsize);
254918687b8SPatrick Sanan   } else {
255a07356b8SSatish Balay     cvode->y = N_VNew_Parallel(cvode->comm_sundials, locsize, glosize);
256918687b8SPatrick Sanan   }
2573c633725SBarry Smith   PetscCheck(cvode->y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cvode->y is not allocated");
2587cb94ee6SHong Zhang 
25928adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
2609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ts->vec_sol, &parray));
2617cb94ee6SHong Zhang   y_data = (PetscScalar *)N_VGetArrayPointer(cvode->y);
2627cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
2639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ts->vec_sol, NULL));
26442528757SHong Zhang 
2659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &cvode->update));
2669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &cvode->ydot));
2679566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)ts, (PetscObject)cvode->update));
2689566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)ts, (PetscObject)cvode->ydot));
2697cb94ee6SHong Zhang 
2707cb94ee6SHong Zhang   /*
2717cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
2727cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
2737cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
2747cb94ee6SHong Zhang   */
2759566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w1));
2769566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w2));
2779566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)ts, (PetscObject)cvode->w1));
2789566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)ts, (PetscObject)cvode->w2));
27916016685SHong Zhang 
28016016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
28116016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
2823c633725SBarry Smith   PetscCheck(mem, PETSC_COMM_SELF, PETSC_ERR_MEM, "CVodeCreate() fails");
28316016685SHong Zhang   cvode->mem = mem;
28416016685SHong Zhang 
28516016685SHong Zhang   /* Set the pointer to user-defined data */
28616016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
28728b400f6SJacob Faibussowitsch   PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetUserData() fails");
28816016685SHong Zhang 
289fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
290cdbf8f93SLisandro Dalcin   flag = CVodeSetInitStep(mem, (realtype)ts->time_step);
29128b400f6SJacob Faibussowitsch   PetscCheck(!flag, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetInitStep() failed");
292f1cd61daSJed Brown   if (cvode->mindt > 0) {
293f1cd61daSJed Brown     flag = CVodeSetMinStep(mem, (realtype)cvode->mindt);
2949f94935aSHong Zhang     if (flag) {
29508401ef6SPierre Jolivet       PetscCheck(flag != CV_MEM_NULL, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, cvode_mem pointer is NULL");
296f7d195e4SLawrence Mitchell       PetscCheck(flag != CV_ILL_INPUT, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
297f7d195e4SLawrence Mitchell       SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed");
2989f94935aSHong Zhang     }
299f1cd61daSJed Brown   }
300f1cd61daSJed Brown   if (cvode->maxdt > 0) {
301f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem, (realtype)cvode->maxdt);
30228b400f6SJacob Faibussowitsch     PetscCheck(!flag, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMaxStep() failed");
303f1cd61daSJed Brown   }
304f1cd61daSJed Brown 
30516016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
306a5b23f4aSJose E. Roman    * user's right hand side function in u'=f(t,u), the initial time T0, and
30716016685SHong Zhang    * the initial dependent variable vector cvode->y */
30816016685SHong Zhang   flag = CVodeInit(mem, TSFunction_Sundials, ts->ptime, cvode->y);
30928b400f6SJacob Faibussowitsch   PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeInit() fails, flag %d", flag);
31016016685SHong Zhang 
3119f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
31216016685SHong Zhang   flag = CVodeSStolerances(mem, cvode->reltol, cvode->abstol);
31328b400f6SJacob Faibussowitsch   PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSStolerances() fails, flag %d", flag);
31416016685SHong Zhang 
315c4e80e11SFlorian   /* Specify max order of BDF / ADAMS method */
316c4e80e11SFlorian   if (cvode->maxord != PETSC_DEFAULT) {
317c4e80e11SFlorian     flag = CVodeSetMaxOrd(mem, cvode->maxord);
31828b400f6SJacob Faibussowitsch     PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetMaxOrd() fails, flag %d", flag);
319c4e80e11SFlorian   }
320c4e80e11SFlorian 
3219f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3229f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem, ts->max_steps);
32328b400f6SJacob Faibussowitsch   PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVodeSetMaxNumSteps() fails, flag %d", flag);
32416016685SHong Zhang 
325918687b8SPatrick Sanan   if (cvode->use_dense) {
326918687b8SPatrick Sanan     /* call CVDense to use a dense linear solver. */
327918687b8SPatrick Sanan     PetscMPIInt size;
328918687b8SPatrick Sanan 
3299566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3303c633725SBarry Smith     PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case");
331918687b8SPatrick Sanan     flag = CVDense(mem, locsize);
33228b400f6SJacob Faibussowitsch     PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVDense() fails, flag %d", flag);
333918687b8SPatrick Sanan   } else {
33416016685SHong Zhang     /* call CVSpgmr to use GMRES as the linear solver.        */
33516016685SHong Zhang     /* setup the ode integrator with the given preconditioner */
3369566063dSJacob Faibussowitsch     PetscCall(TSSundialsGetPC(ts, &pc));
3379566063dSJacob Faibussowitsch     PetscCall(PCGetType(pc, &pctype));
3389566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCNONE, &pcnone));
33916016685SHong Zhang     if (pcnone) {
34016016685SHong Zhang       flag = CVSpgmr(mem, PREC_NONE, 0);
34128b400f6SJacob Faibussowitsch       PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpgmr() fails, flag %d", flag);
34216016685SHong Zhang     } else {
343f61b2b6aSHong Zhang       flag = CVSpgmr(mem, PREC_LEFT, cvode->maxl);
34428b400f6SJacob Faibussowitsch       PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpgmr() fails, flag %d", flag);
34516016685SHong Zhang 
34616016685SHong Zhang       /* Set preconditioner and solve routines Precond and PSolve,
34716016685SHong Zhang          and the pointer to the user-defined block data */
34816016685SHong Zhang       flag = CVSpilsSetPreconditioner(mem, TSPrecond_Sundials, TSPSolve_Sundials);
34928b400f6SJacob Faibussowitsch       PetscCheck(!flag, PETSC_COMM_SELF, PETSC_ERR_LIB, "CVSpilsSetPreconditioner() fails, flag %d", flag);
35016016685SHong Zhang     }
351918687b8SPatrick Sanan   }
3527cb94ee6SHong Zhang   PetscFunctionReturn(0);
3537cb94ee6SHong Zhang }
3547cb94ee6SHong Zhang 
3556fadb2cdSHong Zhang /* type of CVODE linear multistep method */
356c793f718SLisandro Dalcin const char *const TSSundialsLmmTypes[]         = {"", "ADAMS", "BDF", "TSSundialsLmmType", "SUNDIALS_", NULL};
3576fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
358c793f718SLisandro Dalcin const char *const TSSundialsGramSchmidtTypes[] = {"", "MODIFIED", "CLASSICAL", "TSSundialsGramSchmidtType", "SUNDIALS_", NULL};
359a04cf4d8SBarry Smith 
360*9371c9d4SSatish Balay PetscErrorCode TSSetFromOptions_Sundials(TS ts, PetscOptionItems *PetscOptionsObject) {
3617cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
3627cb94ee6SHong Zhang   int          indx;
363ace3abfcSBarry Smith   PetscBool    flag;
3647bda3a07SJed Brown   PC           pc;
3657cb94ee6SHong Zhang 
3667cb94ee6SHong Zhang   PetscFunctionBegin;
367d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SUNDIALS ODE solver options");
3689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-ts_sundials_type", "Scheme", "TSSundialsSetType", TSSundialsLmmTypes, 3, TSSundialsLmmTypes[cvode->cvode_type], &indx, &flag));
3691baa6e33SBarry Smith   if (flag) PetscCall(TSSundialsSetType(ts, (TSSundialsLmmType)indx));
3709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-ts_sundials_gramschmidt_type", "Type of orthogonalization", "TSSundialsSetGramSchmidtType", TSSundialsGramSchmidtTypes, 3, TSSundialsGramSchmidtTypes[cvode->gtype], &indx, &flag));
3711baa6e33SBarry Smith   if (flag) PetscCall(TSSundialsSetGramSchmidtType(ts, (TSSundialsGramSchmidtType)indx));
3729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_sundials_atol", "Absolute tolerance for convergence", "TSSundialsSetTolerance", cvode->abstol, &cvode->abstol, NULL));
3739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_sundials_rtol", "Relative tolerance for convergence", "TSSundialsSetTolerance", cvode->reltol, &cvode->reltol, NULL));
3749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_sundials_mindt", "Minimum step size", "TSSundialsSetMinTimeStep", cvode->mindt, &cvode->mindt, NULL));
3759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_sundials_maxdt", "Maximum step size", "TSSundialsSetMaxTimeStep", cvode->maxdt, &cvode->maxdt, NULL));
3769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_sundials_linear_tolerance", "Convergence tolerance for linear solve", "TSSundialsSetLinearTolerance", cvode->linear_tol, &cvode->linear_tol, NULL));
3779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ts_sundials_maxord", "Max Order for BDF/Adams method", "TSSundialsSetMaxOrd", cvode->maxord, &cvode->maxord, NULL));
3789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ts_sundials_maxl", "Max dimension of the Krylov subspace", "TSSundialsSetMaxl", cvode->maxl, &cvode->maxl, NULL));
3799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_sundials_monitor_steps", "Monitor SUNDIALS internal steps", "TSSundialsMonitorInternalSteps", cvode->monitorstep, &cvode->monitorstep, NULL));
3809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_sundials_use_dense", "Use dense internal solver in SUNDIALS (serial only)", "TSSundialsSetUseDense", cvode->use_dense, &cvode->use_dense, NULL));
381d0609cedSBarry Smith   PetscOptionsHeadEnd();
3829566063dSJacob Faibussowitsch   PetscCall(TSSundialsGetPC(ts, &pc));
3839566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
3847cb94ee6SHong Zhang   PetscFunctionReturn(0);
3857cb94ee6SHong Zhang }
3867cb94ee6SHong Zhang 
387*9371c9d4SSatish Balay PetscErrorCode TSView_Sundials(TS ts, PetscViewer viewer) {
3887cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
3897cb94ee6SHong Zhang   char        *type;
3907cb94ee6SHong Zhang   char         atype[] = "Adams";
3917cb94ee6SHong Zhang   char         btype[] = "BDF: backward differentiation formula";
392ace3abfcSBarry Smith   PetscBool    iascii, isstring;
3932c823083SHong Zhang   long int     nsteps, its, nfevals, nlinsetups, nfails, itmp;
3942c823083SHong Zhang   PetscInt     qlast, qcur;
3955d47aee6SHong Zhang   PetscReal    hinused, hlast, hcur, tcur, tolsfac;
39642528757SHong Zhang   PC           pc;
3977cb94ee6SHong Zhang 
3987cb94ee6SHong Zhang   PetscFunctionBegin;
399bbd56ea5SKarl Rupp   if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
400bbd56ea5SKarl Rupp   else type = btype;
4017cb94ee6SHong Zhang 
4029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
4039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
4047cb94ee6SHong Zhang   if (iascii) {
4059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator does not use SNES!\n"));
4069566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator type %s\n", type));
40763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials maxord %" PetscInt_FMT "\n", cvode->maxord));
4089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials abs tol %g rel tol %g\n", (double)cvode->abstol, (double)cvode->reltol));
409918687b8SPatrick Sanan     if (cvode->use_dense) {
4109566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials integrator using a dense linear solve\n"));
411918687b8SPatrick Sanan     } else {
4129566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials linear solver tolerance factor %g\n", (double)cvode->linear_tol));
41363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials max dimension of Krylov subspace %" PetscInt_FMT "\n", cvode->maxl));
4147cb94ee6SHong Zhang       if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4159566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n"));
4167cb94ee6SHong Zhang       } else {
4179566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n"));
4187cb94ee6SHong Zhang       }
419918687b8SPatrick Sanan     }
4209566063dSJacob Faibussowitsch     if (cvode->mindt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials minimum time step %g\n", (double)cvode->mindt));
4219566063dSJacob Faibussowitsch     if (cvode->maxdt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials maximum time step %g\n", (double)cvode->maxdt));
4222c823083SHong Zhang 
4235d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4249566063dSJacob Faibussowitsch     PetscCall(CVodeGetTolScaleFactor(cvode->mem, &tolsfac));
4259566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials suggested factor for tolerance scaling %g\n", tolsfac));
4269566063dSJacob Faibussowitsch     PetscCall(CVodeGetIntegratorStats(cvode->mem, &nsteps, &nfevals, &nlinsetups, &nfails, &qlast, &qcur, &hinused, &hlast, &hcur, &tcur));
42763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials cumulative number of internal steps %ld\n", nsteps));
42863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of calls to rhs function %ld\n", nfevals));
42963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of calls to linear solver setup function %ld\n", nlinsetups));
43063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of error test failures %ld\n", nfails));
4312c823083SHong Zhang 
4329566063dSJacob Faibussowitsch     PetscCall(CVodeGetNonlinSolvStats(cvode->mem, &its, &nfails));
43363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of nonlinear solver iterations %ld\n", its));
43463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of nonlinear convergence failure %ld\n", nfails));
435918687b8SPatrick Sanan     if (!cvode->use_dense) {
4369566063dSJacob Faibussowitsch       PetscCall(CVSpilsGetNumLinIters(cvode->mem, &its)); /* its = no. of calls to TSPrecond_Sundials() */
43763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of linear iterations %ld\n", its));
4389566063dSJacob Faibussowitsch       PetscCall(CVSpilsGetNumConvFails(cvode->mem, &itmp));
43963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of linear convergence failures %ld\n", itmp));
44042528757SHong Zhang 
4419566063dSJacob Faibussowitsch       PetscCall(TSSundialsGetPC(ts, &pc));
4429566063dSJacob Faibussowitsch       PetscCall(PCView(pc, viewer));
4439566063dSJacob Faibussowitsch       PetscCall(CVSpilsGetNumPrecEvals(cvode->mem, &itmp));
44463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of preconditioner evaluations %ld\n", itmp));
4459566063dSJacob Faibussowitsch       PetscCall(CVSpilsGetNumPrecSolves(cvode->mem, &itmp));
44663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of preconditioner solves %ld\n", itmp));
447918687b8SPatrick Sanan     }
4489566063dSJacob Faibussowitsch     PetscCall(CVSpilsGetNumJtimesEvals(cvode->mem, &itmp));
44963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of Jacobian-vector product evaluations %ld\n", itmp));
4509566063dSJacob Faibussowitsch     PetscCall(CVSpilsGetNumRhsEvals(cvode->mem, &itmp));
45163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Sundials no. of rhs calls for finite diff. Jacobian-vector evals %ld\n", itmp));
4527cb94ee6SHong Zhang   } else if (isstring) {
4539566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, "Sundials type %s", type));
4547cb94ee6SHong Zhang   }
4557cb94ee6SHong Zhang   PetscFunctionReturn(0);
4567cb94ee6SHong Zhang }
4577cb94ee6SHong Zhang 
4587cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
459*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetType_Sundials(TS ts, TSSundialsLmmType type) {
4607cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
4617cb94ee6SHong Zhang 
4627cb94ee6SHong Zhang   PetscFunctionBegin;
4637cb94ee6SHong Zhang   cvode->cvode_type = type;
4647cb94ee6SHong Zhang   PetscFunctionReturn(0);
4657cb94ee6SHong Zhang }
4667cb94ee6SHong Zhang 
467*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts, PetscInt maxl) {
4687cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
4697cb94ee6SHong Zhang 
4707cb94ee6SHong Zhang   PetscFunctionBegin;
471f61b2b6aSHong Zhang   cvode->maxl = maxl;
4727cb94ee6SHong Zhang   PetscFunctionReturn(0);
4737cb94ee6SHong Zhang }
4747cb94ee6SHong Zhang 
475*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts, PetscReal tol) {
4767cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
4777cb94ee6SHong Zhang 
4787cb94ee6SHong Zhang   PetscFunctionBegin;
4797cb94ee6SHong Zhang   cvode->linear_tol = tol;
4807cb94ee6SHong Zhang   PetscFunctionReturn(0);
4817cb94ee6SHong Zhang }
4827cb94ee6SHong Zhang 
483*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts, TSSundialsGramSchmidtType type) {
4847cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
4857cb94ee6SHong Zhang 
4867cb94ee6SHong Zhang   PetscFunctionBegin;
4877cb94ee6SHong Zhang   cvode->gtype = type;
4887cb94ee6SHong Zhang   PetscFunctionReturn(0);
4897cb94ee6SHong Zhang }
4907cb94ee6SHong Zhang 
491*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts, PetscReal aabs, PetscReal rel) {
4927cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
4937cb94ee6SHong Zhang 
4947cb94ee6SHong Zhang   PetscFunctionBegin;
4957cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
4967cb94ee6SHong Zhang   if (rel != PETSC_DECIDE) cvode->reltol = rel;
4977cb94ee6SHong Zhang   PetscFunctionReturn(0);
4987cb94ee6SHong Zhang }
4997cb94ee6SHong Zhang 
500*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts, PetscReal mindt) {
501f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials *)ts->data;
502f1cd61daSJed Brown 
503f1cd61daSJed Brown   PetscFunctionBegin;
504f1cd61daSJed Brown   cvode->mindt = mindt;
505f1cd61daSJed Brown   PetscFunctionReturn(0);
506f1cd61daSJed Brown }
507f1cd61daSJed Brown 
508*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts, PetscReal maxdt) {
509f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials *)ts->data;
510f1cd61daSJed Brown 
511f1cd61daSJed Brown   PetscFunctionBegin;
512f1cd61daSJed Brown   cvode->maxdt = maxdt;
513f1cd61daSJed Brown   PetscFunctionReturn(0);
514f1cd61daSJed Brown }
515918687b8SPatrick Sanan 
516*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts, PetscBool use_dense) {
517918687b8SPatrick Sanan   TS_Sundials *cvode = (TS_Sundials *)ts->data;
518918687b8SPatrick Sanan 
519918687b8SPatrick Sanan   PetscFunctionBegin;
520918687b8SPatrick Sanan   cvode->use_dense = use_dense;
521918687b8SPatrick Sanan   PetscFunctionReturn(0);
522918687b8SPatrick Sanan }
523918687b8SPatrick Sanan 
524*9371c9d4SSatish Balay PetscErrorCode TSSundialsGetPC_Sundials(TS ts, PC *pc) {
525dcbc6d53SSean Farley   SNES snes;
526dcbc6d53SSean Farley   KSP  ksp;
5277cb94ee6SHong Zhang 
5287cb94ee6SHong Zhang   PetscFunctionBegin;
5299566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
5309566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
5319566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, pc));
5327cb94ee6SHong Zhang   PetscFunctionReturn(0);
5337cb94ee6SHong Zhang }
5347cb94ee6SHong Zhang 
535*9371c9d4SSatish Balay PetscErrorCode TSSundialsGetIterations_Sundials(TS ts, int *nonlin, int *lin) {
5367cb94ee6SHong Zhang   PetscFunctionBegin;
5375ef26d82SJed Brown   if (nonlin) *nonlin = ts->snes_its;
5385ef26d82SJed Brown   if (lin) *lin = ts->ksp_its;
5397cb94ee6SHong Zhang   PetscFunctionReturn(0);
5407cb94ee6SHong Zhang }
5417cb94ee6SHong Zhang 
542*9371c9d4SSatish Balay PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts, PetscBool s) {
5432bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials *)ts->data;
5442bfc04deSHong Zhang 
5452bfc04deSHong Zhang   PetscFunctionBegin;
5462bfc04deSHong Zhang   cvode->monitorstep = s;
5472bfc04deSHong Zhang   PetscFunctionReturn(0);
5482bfc04deSHong Zhang }
5497cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
5507cb94ee6SHong Zhang 
5517cb94ee6SHong Zhang /*@C
5527cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
5537cb94ee6SHong Zhang 
5547cb94ee6SHong Zhang    Not Collective
5557cb94ee6SHong Zhang 
55697bb3fdcSJose E. Roman    Input Parameter:
5577cb94ee6SHong Zhang .    ts     - the time-step context
5587cb94ee6SHong Zhang 
5597cb94ee6SHong Zhang    Output Parameters:
5607cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
5617cb94ee6SHong Zhang -   lin    - number of linear iterations
5627cb94ee6SHong Zhang 
5637cb94ee6SHong Zhang    Level: advanced
5647cb94ee6SHong Zhang 
5657cb94ee6SHong Zhang    Notes:
5667cb94ee6SHong Zhang     These return the number since the creation of the TS object
5677cb94ee6SHong Zhang 
568db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
569db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
570db781477SPatrick Sanan           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
571db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()`
5727cb94ee6SHong Zhang 
5737cb94ee6SHong Zhang @*/
574*9371c9d4SSatish Balay PetscErrorCode TSSundialsGetIterations(TS ts, int *nonlin, int *lin) {
5757cb94ee6SHong Zhang   PetscFunctionBegin;
576cac4c232SBarry Smith   PetscUseMethod(ts, "TSSundialsGetIterations_C", (TS, int *, int *), (ts, nonlin, lin));
5777cb94ee6SHong Zhang   PetscFunctionReturn(0);
5787cb94ee6SHong Zhang }
5797cb94ee6SHong Zhang 
5807cb94ee6SHong Zhang /*@
5817cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
5827cb94ee6SHong Zhang 
5833f9fe445SBarry Smith    Logically Collective on TS
5847cb94ee6SHong Zhang 
5858f7d6fe5SPatrick Sanan    Input Parameters:
5867cb94ee6SHong Zhang +    ts     - the time-step context
5877cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
5887cb94ee6SHong Zhang 
5897cb94ee6SHong Zhang    Level: intermediate
5907cb94ee6SHong Zhang 
591db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`,
592db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
593db781477SPatrick Sanan           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
594db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
595db781477SPatrick Sanan           `TSSetExactFinalTime()`
5967cb94ee6SHong Zhang @*/
597*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetType(TS ts, TSSundialsLmmType type) {
5987cb94ee6SHong Zhang   PetscFunctionBegin;
599cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetType_C", (TS, TSSundialsLmmType), (ts, type));
6007cb94ee6SHong Zhang   PetscFunctionReturn(0);
6017cb94ee6SHong Zhang }
6027cb94ee6SHong Zhang 
6037cb94ee6SHong Zhang /*@
604c4e80e11SFlorian    TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS.
605c4e80e11SFlorian 
606c4e80e11SFlorian    Logically Collective on TS
607c4e80e11SFlorian 
6088f7d6fe5SPatrick Sanan    Input Parameters:
609c4e80e11SFlorian +    ts      - the time-step context
610c4e80e11SFlorian -    maxord  - maximum order of BDF / Adams method
611c4e80e11SFlorian 
612c4e80e11SFlorian    Level: advanced
613c4e80e11SFlorian 
614db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`,
615db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
616db781477SPatrick Sanan           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
617db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
618db781477SPatrick Sanan           `TSSetExactFinalTime()`
619c4e80e11SFlorian 
620c4e80e11SFlorian @*/
621*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxord(TS ts, PetscInt maxord) {
622c4e80e11SFlorian   PetscFunctionBegin;
623c4e80e11SFlorian   PetscValidLogicalCollectiveInt(ts, maxord, 2);
624cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetMaxOrd_C", (TS, PetscInt), (ts, maxord));
625c4e80e11SFlorian   PetscFunctionReturn(0);
626c4e80e11SFlorian }
627c4e80e11SFlorian 
628c4e80e11SFlorian /*@
629f61b2b6aSHong Zhang    TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
6307cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
631f61b2b6aSHong Zhang        this is the maximum number of GMRES steps that will be used.
6327cb94ee6SHong Zhang 
6333f9fe445SBarry Smith    Logically Collective on TS
6347cb94ee6SHong Zhang 
6358f7d6fe5SPatrick Sanan    Input Parameters:
6367cb94ee6SHong Zhang +    ts      - the time-step context
637f61b2b6aSHong Zhang -    maxl - number of direction vectors (the dimension of Krylov subspace).
6387cb94ee6SHong Zhang 
6397cb94ee6SHong Zhang    Level: advanced
6407cb94ee6SHong Zhang 
641db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`,
642db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
643db781477SPatrick Sanan           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
644db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
645db781477SPatrick Sanan           `TSSetExactFinalTime()`
6467cb94ee6SHong Zhang 
6477cb94ee6SHong Zhang @*/
648*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxl(TS ts, PetscInt maxl) {
6497cb94ee6SHong Zhang   PetscFunctionBegin;
650f61b2b6aSHong Zhang   PetscValidLogicalCollectiveInt(ts, maxl, 2);
651cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetMaxl_C", (TS, PetscInt), (ts, maxl));
6527cb94ee6SHong Zhang   PetscFunctionReturn(0);
6537cb94ee6SHong Zhang }
6547cb94ee6SHong Zhang 
6557cb94ee6SHong Zhang /*@
6567cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
6577cb94ee6SHong Zhang        system by SUNDIALS.
6587cb94ee6SHong Zhang 
6593f9fe445SBarry Smith    Logically Collective on TS
6607cb94ee6SHong Zhang 
6618f7d6fe5SPatrick Sanan    Input Parameters:
6627cb94ee6SHong Zhang +    ts     - the time-step context
6637cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
6647cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
6657cb94ee6SHong Zhang 
6667cb94ee6SHong Zhang    Level: advanced
6677cb94ee6SHong Zhang 
668db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
669db781477SPatrick Sanan           `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
670db781477SPatrick Sanan           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
671db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
672db781477SPatrick Sanan           `TSSetExactFinalTime()`
6737cb94ee6SHong Zhang 
6747cb94ee6SHong Zhang @*/
675*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetLinearTolerance(TS ts, PetscReal tol) {
6767cb94ee6SHong Zhang   PetscFunctionBegin;
677c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts, tol, 2);
678cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetLinearTolerance_C", (TS, PetscReal), (ts, tol));
6797cb94ee6SHong Zhang   PetscFunctionReturn(0);
6807cb94ee6SHong Zhang }
6817cb94ee6SHong Zhang 
6827cb94ee6SHong Zhang /*@
6837cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
6847cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
6857cb94ee6SHong Zhang 
6863f9fe445SBarry Smith    Logically Collective on TS
6877cb94ee6SHong Zhang 
6888f7d6fe5SPatrick Sanan    Input Parameters:
6897cb94ee6SHong Zhang +    ts  - the time-step context
6907cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
6917cb94ee6SHong Zhang 
6927cb94ee6SHong Zhang    Level: advanced
6937cb94ee6SHong Zhang 
694db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
695db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`,
696db781477SPatrick Sanan           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
697db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
698db781477SPatrick Sanan           `TSSetExactFinalTime()`
6997cb94ee6SHong Zhang 
7007cb94ee6SHong Zhang @*/
701*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetGramSchmidtType(TS ts, TSSundialsGramSchmidtType type) {
7027cb94ee6SHong Zhang   PetscFunctionBegin;
703cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetGramSchmidtType_C", (TS, TSSundialsGramSchmidtType), (ts, type));
7047cb94ee6SHong Zhang   PetscFunctionReturn(0);
7057cb94ee6SHong Zhang }
7067cb94ee6SHong Zhang 
7077cb94ee6SHong Zhang /*@
7087cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
7097cb94ee6SHong Zhang                          Sundials for error control.
7107cb94ee6SHong Zhang 
7113f9fe445SBarry Smith    Logically Collective on TS
7127cb94ee6SHong Zhang 
7138f7d6fe5SPatrick Sanan    Input Parameters:
7147cb94ee6SHong Zhang +    ts  - the time-step context
7157cb94ee6SHong Zhang .    aabs - the absolute tolerance
7167cb94ee6SHong Zhang -    rel - the relative tolerance
7177cb94ee6SHong Zhang 
7187cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
7197cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
7207cb94ee6SHong Zhang 
7217cb94ee6SHong Zhang    Level: intermediate
7227cb94ee6SHong Zhang 
723db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`,
724db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`,
725db781477SPatrick Sanan           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
726db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
727db781477SPatrick Sanan           `TSSetExactFinalTime()`
7287cb94ee6SHong Zhang 
7297cb94ee6SHong Zhang @*/
730*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetTolerance(TS ts, PetscReal aabs, PetscReal rel) {
7317cb94ee6SHong Zhang   PetscFunctionBegin;
732cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetTolerance_C", (TS, PetscReal, PetscReal), (ts, aabs, rel));
7337cb94ee6SHong Zhang   PetscFunctionReturn(0);
7347cb94ee6SHong Zhang }
7357cb94ee6SHong Zhang 
7367cb94ee6SHong Zhang /*@
7377cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
7387cb94ee6SHong Zhang 
7397cb94ee6SHong Zhang    Input Parameter:
7407cb94ee6SHong Zhang .    ts - the time-step context
7417cb94ee6SHong Zhang 
7427cb94ee6SHong Zhang    Output Parameter:
7437cb94ee6SHong Zhang .    pc - the preconditioner context
7447cb94ee6SHong Zhang 
7457cb94ee6SHong Zhang    Level: advanced
7467cb94ee6SHong Zhang 
747db781477SPatrick Sanan .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
748db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
749db781477SPatrick Sanan           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
750db781477SPatrick Sanan           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`
7517cb94ee6SHong Zhang @*/
752*9371c9d4SSatish Balay PetscErrorCode TSSundialsGetPC(TS ts, PC *pc) {
7537cb94ee6SHong Zhang   PetscFunctionBegin;
754cac4c232SBarry Smith   PetscUseMethod(ts, "TSSundialsGetPC_C", (TS, PC *), (ts, pc));
7557cb94ee6SHong Zhang   PetscFunctionReturn(0);
7567cb94ee6SHong Zhang }
7577cb94ee6SHong Zhang 
758f1cd61daSJed Brown /*@
759f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
760f1cd61daSJed Brown 
761d8d19677SJose E. Roman    Input Parameters:
762f1cd61daSJed Brown +   ts - the time-step context
763f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
764f1cd61daSJed Brown 
765fc6b9e64SJed Brown    Note:
766fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
767fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
768fc6b9e64SJed Brown 
769f1cd61daSJed Brown    Level: beginner
770f1cd61daSJed Brown 
771db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`,
772f1cd61daSJed Brown @*/
773*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetMinTimeStep(TS ts, PetscReal mindt) {
774f1cd61daSJed Brown   PetscFunctionBegin;
775cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetMinTimeStep_C", (TS, PetscReal), (ts, mindt));
776f1cd61daSJed Brown   PetscFunctionReturn(0);
777f1cd61daSJed Brown }
778f1cd61daSJed Brown 
779f1cd61daSJed Brown /*@
780f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
781f1cd61daSJed Brown 
782d8d19677SJose E. Roman    Input Parameters:
783f1cd61daSJed Brown +   ts - the time-step context
784f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
785f1cd61daSJed Brown 
786f1cd61daSJed Brown    Level: beginner
787f1cd61daSJed Brown 
788db781477SPatrick Sanan .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`,
789f1cd61daSJed Brown @*/
790*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetMaxTimeStep(TS ts, PetscReal maxdt) {
791f1cd61daSJed Brown   PetscFunctionBegin;
792cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt));
793f1cd61daSJed Brown   PetscFunctionReturn(0);
794f1cd61daSJed Brown }
795f1cd61daSJed Brown 
7962bfc04deSHong Zhang /*@
7972bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
7982bfc04deSHong Zhang 
799d8d19677SJose E. Roman    Input Parameters:
8002bfc04deSHong Zhang +   ts - the time-step context
8012bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
8022bfc04deSHong Zhang 
8032bfc04deSHong Zhang    Level: beginner
8042bfc04deSHong Zhang 
805f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
8062bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
807f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
8082bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
8092bfc04deSHong Zhang @*/
810*9371c9d4SSatish Balay PetscErrorCode TSSundialsMonitorInternalSteps(TS ts, PetscBool ft) {
8112bfc04deSHong Zhang   PetscFunctionBegin;
812cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsMonitorInternalSteps_C", (TS, PetscBool), (ts, ft));
8132bfc04deSHong Zhang   PetscFunctionReturn(0);
8142bfc04deSHong Zhang }
815918687b8SPatrick Sanan 
816918687b8SPatrick Sanan /*@
817918687b8SPatrick Sanan    TSSundialsSetUseDense - Set a flag to use a dense linear solver in SUNDIALS (serial only)
818918687b8SPatrick Sanan 
819918687b8SPatrick Sanan    Logically Collective
820918687b8SPatrick Sanan 
821918687b8SPatrick Sanan    Input Parameters:
822918687b8SPatrick Sanan +    ts         - the time-step context
823918687b8SPatrick Sanan -    use_dense  - PETSC_TRUE to use the dense solver
824918687b8SPatrick Sanan 
825918687b8SPatrick Sanan    Level: advanced
826918687b8SPatrick Sanan 
827db781477SPatrick Sanan .seealso: `TSSUNDIALS`
828918687b8SPatrick Sanan 
829918687b8SPatrick Sanan @*/
830*9371c9d4SSatish Balay PetscErrorCode TSSundialsSetUseDense(TS ts, PetscBool use_dense) {
831918687b8SPatrick Sanan   PetscFunctionBegin;
832918687b8SPatrick Sanan   PetscValidLogicalCollectiveInt(ts, use_dense, 2);
833cac4c232SBarry Smith   PetscTryMethod(ts, "TSSundialsSetUseDense_C", (TS, PetscBool), (ts, use_dense));
834918687b8SPatrick Sanan   PetscFunctionReturn(0);
835918687b8SPatrick Sanan }
836918687b8SPatrick Sanan 
8377cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
8387cb94ee6SHong Zhang /*MC
83996f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
8407cb94ee6SHong Zhang 
8417cb94ee6SHong Zhang    Options Database:
842897c9f78SHong Zhang +    -ts_sundials_type <bdf,adams> -
8437cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
8447cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
8457cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
846897c9f78SHong Zhang .    -ts_sundials_linear_tolerance <tol> -
847f61b2b6aSHong Zhang .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
848918687b8SPatrick Sanan .    -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps
849918687b8SPatrick Sanan -    -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only)
85016016685SHong Zhang 
85195452b02SPatrick Sanan     Notes:
85295452b02SPatrick Sanan     This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply,
853897c9f78SHong Zhang            only PETSc PC options.
8547cb94ee6SHong Zhang 
8557cb94ee6SHong Zhang     Level: beginner
8567cb94ee6SHong Zhang 
857db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`,
858db781477SPatrick Sanan           `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()`
8597cb94ee6SHong Zhang 
8607cb94ee6SHong Zhang M*/
861*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts) {
8627cb94ee6SHong Zhang   TS_Sundials *cvode;
86342528757SHong Zhang   PC           pc;
8647cb94ee6SHong Zhang 
8657cb94ee6SHong Zhang   PetscFunctionBegin;
866277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
86728adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
86828adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
869214bc6a2SJed Brown   ts->ops->setup          = TSSetUp_Sundials;
870b4eba00bSSean Farley   ts->ops->step           = TSStep_Sundials;
871b4eba00bSSean Farley   ts->ops->interpolate    = TSInterpolate_Sundials;
872214bc6a2SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
8732ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
8747cb94ee6SHong Zhang 
8759566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts, &cvode));
876bbd56ea5SKarl Rupp 
877efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
878efd4aadfSBarry Smith 
8797cb94ee6SHong Zhang   ts->data           = (void *)cvode;
8806fadb2cdSHong Zhang   cvode->cvode_type  = SUNDIALS_BDF;
8816fadb2cdSHong Zhang   cvode->gtype       = SUNDIALS_CLASSICAL_GS;
882f61b2b6aSHong Zhang   cvode->maxl        = 5;
883c4e80e11SFlorian   cvode->maxord      = PETSC_DEFAULT;
8847cb94ee6SHong Zhang   cvode->linear_tol  = .05;
885b4eba00bSSean Farley   cvode->monitorstep = PETSC_TRUE;
886918687b8SPatrick Sanan   cvode->use_dense   = PETSC_FALSE;
8877cb94ee6SHong Zhang 
8889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts), &(cvode->comm_sundials)));
889f1cd61daSJed Brown 
890f1cd61daSJed Brown   cvode->mindt = -1.;
891f1cd61daSJed Brown   cvode->maxdt = -1.;
892f1cd61daSJed Brown 
8937cb94ee6SHong Zhang   /* set tolerance for Sundials */
8947cb94ee6SHong Zhang   cvode->reltol = 1e-6;
8952c823083SHong Zhang   cvode->abstol = 1e-6;
8967cb94ee6SHong Zhang 
89742528757SHong Zhang   /* set PCNONE as default pctype */
8989566063dSJacob Faibussowitsch   PetscCall(TSSundialsGetPC_Sundials(ts, &pc));
8999566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
90042528757SHong Zhang 
9019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", TSSundialsSetType_Sundials));
9029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", TSSundialsSetMaxl_Sundials));
9039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", TSSundialsSetLinearTolerance_Sundials));
9049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", TSSundialsSetGramSchmidtType_Sundials));
9059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", TSSundialsSetTolerance_Sundials));
9069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", TSSundialsSetMinTimeStep_Sundials));
9079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", TSSundialsSetMaxTimeStep_Sundials));
9089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", TSSundialsGetPC_Sundials));
9099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", TSSundialsGetIterations_Sundials));
9109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", TSSundialsMonitorInternalSteps_Sundials));
9119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", TSSundialsSetUseDense_Sundials));
9127cb94ee6SHong Zhang   PetscFunctionReturn(0);
9137cb94ee6SHong Zhang }
914