17cb94ee6SHong Zhang /*
2bcf0153eSBarry Smith Provides a PETSc interface to version 2.5 of SUNDIALS/CVODE solver (a very old version)
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 */
TSPrecond_Sundials_Petsc(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)143ba16761SJacob Faibussowitsch static PetscErrorCode TSPrecond_Sundials_Petsc(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)
15d71ae5a4SJacob Faibussowitsch {
167cb94ee6SHong Zhang TS ts = (TS)P_data;
177cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
18dcbc6d53SSean Farley PC pc;
190679a0aeSJed Brown Mat J, P;
200679a0aeSJed Brown Vec yy = cvode->w1, yydot = cvode->ydot;
210679a0aeSJed Brown PetscReal gm = (PetscReal)_gamma;
227cb94ee6SHong Zhang PetscScalar *y_data;
237cb94ee6SHong Zhang
247cb94ee6SHong Zhang PetscFunctionBegin;
259566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &P, NULL, NULL));
267cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(y);
279566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(yy, y_data));
289566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot)); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
290679a0aeSJed Brown /* compute the shifted Jacobian (1/gm)*I + Jrest */
309566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ts->ptime, yy, yydot, 1 / gm, J, P, PETSC_FALSE));
319566063dSJacob Faibussowitsch PetscCall(VecResetArray(yy));
32bcf0153eSBarry Smith PetscCall(MatScale(P, gm)); /* turn into I-gm*Jrest, J is not used by SUNDIALS */
337cb94ee6SHong Zhang *jcurPtr = TRUE;
349566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc));
359566063dSJacob Faibussowitsch PetscCall(PCSetOperators(pc, J, P));
363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
373ba16761SJacob Faibussowitsch }
383ba16761SJacob Faibussowitsch
393ba16761SJacob Faibussowitsch /* Sundial expects an int (*)(args...) but PetscErrorCode is an enum. Instead of switching out
403ba16761SJacob Faibussowitsch all the PetscCalls in TSPrecond_Sundials_Petsc we just wrap it */
TSPrecond_Sundials_Private(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)413ba16761SJacob Faibussowitsch static int TSPrecond_Sundials_Private(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)
423ba16761SJacob Faibussowitsch {
433ba16761SJacob Faibussowitsch return (int)TSPrecond_Sundials_Petsc(tn, y, fy, jok, jcurPtr, _gamma, P_data, vtemp1, vtemp2, vtemp3);
447cb94ee6SHong Zhang }
457cb94ee6SHong Zhang
467cb94ee6SHong Zhang /*
47bcf0153eSBarry Smith TSPSolve_Sundials - routine that we provide to SUNDIALS that applies the preconditioner.
487cb94ee6SHong Zhang */
TSPSolve_Sundials_Petsc(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)493ba16761SJacob Faibussowitsch static PetscErrorCode TSPSolve_Sundials_Petsc(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)
50d71ae5a4SJacob Faibussowitsch {
517cb94ee6SHong Zhang TS ts = (TS)P_data;
527cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
53dcbc6d53SSean Farley PC pc;
547cb94ee6SHong Zhang Vec rr = cvode->w1, zz = cvode->w2;
557cb94ee6SHong Zhang PetscScalar *r_data, *z_data;
567cb94ee6SHong Zhang
577cb94ee6SHong Zhang PetscFunctionBegin;
587cb94ee6SHong Zhang /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
597cb94ee6SHong Zhang r_data = (PetscScalar *)N_VGetArrayPointer(r);
607cb94ee6SHong Zhang z_data = (PetscScalar *)N_VGetArrayPointer(z);
619566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(rr, r_data));
629566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(zz, z_data));
634ac7836bSHong Zhang
647cb94ee6SHong Zhang /* Solve the Px=r and put the result in zz */
659566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc));
669566063dSJacob Faibussowitsch PetscCall(PCApply(pc, rr, zz));
679566063dSJacob Faibussowitsch PetscCall(VecResetArray(rr));
689566063dSJacob Faibussowitsch PetscCall(VecResetArray(zz));
693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
703ba16761SJacob Faibussowitsch }
713ba16761SJacob Faibussowitsch
723ba16761SJacob Faibussowitsch /* See TSPrecond_Sundials_Private() */
TSPSolve_Sundials_Private(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)733ba16761SJacob Faibussowitsch static int TSPSolve_Sundials_Private(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)
743ba16761SJacob Faibussowitsch {
753ba16761SJacob Faibussowitsch return (int)TSPSolve_Sundials_Petsc(tn, y, fy, r, z, _gamma, delta, lr, P_data, vtemp);
767cb94ee6SHong Zhang }
777cb94ee6SHong Zhang
787cb94ee6SHong Zhang /*
79dd8e379bSPierre Jolivet TSFunction_Sundials - routine that we provide to SUNDIALS that applies the right-hand side.
807cb94ee6SHong Zhang */
TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,PetscCtx ctx)81*2a8381b2SBarry Smith static int TSFunction_Sundials(realtype t, N_Vector y, N_Vector ydot, PetscCtx ctx)
82d71ae5a4SJacob Faibussowitsch {
837cb94ee6SHong Zhang TS ts = (TS)ctx;
845fcef5e4SJed Brown DM dm;
85942e3340SBarry Smith DMTS tsdm;
868434afd1SBarry Smith TSIFunctionFn *ifunction;
87ce94432eSBarry Smith MPI_Comm comm;
887cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
890679a0aeSJed Brown Vec yy = cvode->w1, yyd = cvode->w2, yydot = cvode->ydot;
907cb94ee6SHong Zhang PetscScalar *y_data, *ydot_data;
917cb94ee6SHong Zhang
927cb94ee6SHong Zhang PetscFunctionBegin;
939566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
945ff03cf7SDinesh Kaushik /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
957cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(y);
967cb94ee6SHong Zhang ydot_data = (PetscScalar *)N_VGetArrayPointer(ydot);
979566063dSJacob Faibussowitsch PetscCallAbort(comm, VecPlaceArray(yy, y_data));
989566063dSJacob Faibussowitsch PetscCallAbort(comm, VecPlaceArray(yyd, ydot_data));
994ac7836bSHong Zhang
100dd8e379bSPierre Jolivet /* Now compute the right-hand side function, via IFunction unless only the more efficient RHSFunction is set */
1019566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
1029566063dSJacob Faibussowitsch PetscCall(DMGetDMTS(dm, &tsdm));
1039566063dSJacob Faibussowitsch PetscCall(DMTSGetIFunction(dm, &ifunction, NULL));
1045fcef5e4SJed Brown if (!ifunction) {
1059566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t, yy, yyd));
106bc0cc02bSJed Brown } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */
1079566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(yydot));
1089566063dSJacob Faibussowitsch PetscCallAbort(comm, TSComputeIFunction(ts, t, yy, yydot, yyd, PETSC_FALSE));
1099566063dSJacob Faibussowitsch PetscCall(VecScale(yyd, -1.));
110bc0cc02bSJed Brown }
1119566063dSJacob Faibussowitsch PetscCallAbort(comm, VecResetArray(yy));
1129566063dSJacob Faibussowitsch PetscCallAbort(comm, VecResetArray(yyd));
1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1147cb94ee6SHong Zhang }
1157cb94ee6SHong Zhang
1167cb94ee6SHong Zhang /*
117bcf0153eSBarry Smith TSStep_Sundials - Calls SUNDIALS to integrate the ODE.
1187cb94ee6SHong Zhang */
TSStep_Sundials(TS ts)11966976f2fSJacob Faibussowitsch static PetscErrorCode TSStep_Sundials(TS ts)
120d71ae5a4SJacob Faibussowitsch {
1217cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
122b4eba00bSSean Farley PetscInt flag;
123be5899b3SLisandro Dalcin long int nits, lits, nsteps;
1247cb94ee6SHong Zhang realtype t, tout;
1257cb94ee6SHong Zhang PetscScalar *y_data;
1267cb94ee6SHong Zhang void *mem;
1277cb94ee6SHong Zhang
1287cb94ee6SHong Zhang PetscFunctionBegin;
12916016685SHong Zhang mem = cvode->mem;
1307cb94ee6SHong Zhang tout = ts->max_time;
1319566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol, &y_data));
1327cb94ee6SHong Zhang N_VSetArrayPointer((realtype *)y_data, cvode->y);
1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol, NULL));
134186e87acSLisandro Dalcin
1359687d888SLisandro Dalcin /* We would like to TSPreStage() and TSPostStage()
1369687d888SLisandro Dalcin * before each stage solve but CVode does not appear to support this. */
1379371c9d4SSatish Balay if (cvode->monitorstep) flag = CVode(mem, tout, cvode->y, &t, CV_ONE_STEP);
1389371c9d4SSatish Balay else flag = CVode(mem, tout, cvode->y, &t, CV_NORMAL);
1399f94935aSHong Zhang
1409f94935aSHong Zhang if (flag) { /* display error message */
1419f94935aSHong Zhang switch (flag) {
142d71ae5a4SJacob Faibussowitsch case CV_ILL_INPUT:
143d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ILL_INPUT");
144d71ae5a4SJacob Faibussowitsch break;
145d71ae5a4SJacob Faibussowitsch case CV_TOO_CLOSE:
146d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_CLOSE");
147d71ae5a4SJacob Faibussowitsch break;
1483c7fefeeSJed Brown case CV_TOO_MUCH_WORK: {
1499f94935aSHong Zhang PetscReal tcur;
1503ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetNumSteps, mem, &nsteps);
1513ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetCurrentTime, mem, &tcur);
15263a3b9bcSJacob 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);
1533c7fefeeSJed Brown } break;
154d71ae5a4SJacob Faibussowitsch case CV_TOO_MUCH_ACC:
155d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_TOO_MUCH_ACC");
156d71ae5a4SJacob Faibussowitsch break;
157d71ae5a4SJacob Faibussowitsch case CV_ERR_FAILURE:
158d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_ERR_FAILURE");
159d71ae5a4SJacob Faibussowitsch break;
160d71ae5a4SJacob Faibussowitsch case CV_CONV_FAILURE:
161d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_CONV_FAILURE");
162d71ae5a4SJacob Faibussowitsch break;
163d71ae5a4SJacob Faibussowitsch case CV_LINIT_FAIL:
164d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LINIT_FAIL");
165d71ae5a4SJacob Faibussowitsch break;
166d71ae5a4SJacob Faibussowitsch case CV_LSETUP_FAIL:
167d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSETUP_FAIL");
168d71ae5a4SJacob Faibussowitsch break;
169d71ae5a4SJacob Faibussowitsch case CV_LSOLVE_FAIL:
170d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_LSOLVE_FAIL");
171d71ae5a4SJacob Faibussowitsch break;
172d71ae5a4SJacob Faibussowitsch case CV_RHSFUNC_FAIL:
173d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RHSFUNC_FAIL");
174d71ae5a4SJacob Faibussowitsch break;
175d71ae5a4SJacob Faibussowitsch case CV_FIRST_RHSFUNC_ERR:
176d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_FIRST_RHSFUNC_ERR");
177d71ae5a4SJacob Faibussowitsch break;
178d71ae5a4SJacob Faibussowitsch case CV_REPTD_RHSFUNC_ERR:
179d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_REPTD_RHSFUNC_ERR");
180d71ae5a4SJacob Faibussowitsch break;
181d71ae5a4SJacob Faibussowitsch case CV_UNREC_RHSFUNC_ERR:
182d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_UNREC_RHSFUNC_ERR");
183d71ae5a4SJacob Faibussowitsch break;
184d71ae5a4SJacob Faibussowitsch case CV_RTFUNC_FAIL:
185d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, CV_RTFUNC_FAIL");
186d71ae5a4SJacob Faibussowitsch break;
187d71ae5a4SJacob Faibussowitsch default:
188d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "CVode() fails, flag %d", flag);
1899f94935aSHong Zhang }
1909f94935aSHong Zhang }
1919f94935aSHong Zhang
192be5899b3SLisandro Dalcin /* log inner nonlinear and linear iterations */
1933ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetNumNonlinSolvIters, mem, &nits);
1943ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumLinIters, mem, &lits);
1959371c9d4SSatish Balay ts->snes_its += nits;
1969371c9d4SSatish Balay ts->ksp_its = lits;
197be5899b3SLisandro Dalcin
1987cb94ee6SHong Zhang /* copy the solution from cvode->y to cvode->update and sol */
1999566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(cvode->w1, y_data));
2009566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->w1, cvode->update));
2019566063dSJacob Faibussowitsch PetscCall(VecResetArray(cvode->w1));
2029566063dSJacob Faibussowitsch PetscCall(VecCopy(cvode->update, ts->vec_sol));
203186e87acSLisandro Dalcin
204186e87acSLisandro Dalcin ts->time_step = t - ts->ptime;
205186e87acSLisandro Dalcin ts->ptime = t;
206186e87acSLisandro Dalcin
2073ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetNumSteps, mem, &nsteps);
2082808aa04SLisandro Dalcin if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */
2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
210b4eba00bSSean Farley }
211b4eba00bSSean Farley
TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)212d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Sundials(TS ts, PetscReal t, Vec X)
213d71ae5a4SJacob Faibussowitsch {
214b4eba00bSSean Farley TS_Sundials *cvode = (TS_Sundials *)ts->data;
215b4eba00bSSean Farley N_Vector y;
216b4eba00bSSean Farley PetscScalar *x_data;
217b4eba00bSSean Farley PetscInt glosize, locsize;
218b4eba00bSSean Farley
219b4eba00bSSean Farley PetscFunctionBegin;
220b4eba00bSSean Farley /* get the vector size */
2219566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &glosize));
2229566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &locsize));
2239566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x_data));
224b4eba00bSSean Farley
2254d78c9e0SClaas Abert /* Initialize N_Vec y with x_data */
226918687b8SPatrick Sanan if (cvode->use_dense) {
227918687b8SPatrick Sanan PetscMPIInt size;
228918687b8SPatrick Sanan
2299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2303c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case");
231918687b8SPatrick Sanan y = N_VMake_Serial(locsize, (realtype *)x_data);
232918687b8SPatrick Sanan } else {
2334d78c9e0SClaas Abert y = N_VMake_Parallel(cvode->comm_sundials, locsize, glosize, (realtype *)x_data);
234918687b8SPatrick Sanan }
235918687b8SPatrick Sanan
2363c633725SBarry Smith PetscCheck(y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Interpolated y is not allocated");
237b4eba00bSSean Farley
2383ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetDky, cvode->mem, t, 0, y);
2399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x_data));
2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2417cb94ee6SHong Zhang }
2427cb94ee6SHong Zhang
TSReset_Sundials(TS ts)24366976f2fSJacob Faibussowitsch static PetscErrorCode TSReset_Sundials(TS ts)
244d71ae5a4SJacob Faibussowitsch {
245277b19d0SLisandro Dalcin TS_Sundials *cvode = (TS_Sundials *)ts->data;
246277b19d0SLisandro Dalcin
247277b19d0SLisandro Dalcin PetscFunctionBegin;
2489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->update));
2499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->ydot));
2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w1));
2519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cvode->w2));
252bbd56ea5SKarl Rupp if (cvode->mem) CVodeFree(&cvode->mem);
2533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
254277b19d0SLisandro Dalcin }
255277b19d0SLisandro Dalcin
TSDestroy_Sundials(TS ts)25666976f2fSJacob Faibussowitsch static PetscErrorCode TSDestroy_Sundials(TS ts)
257d71ae5a4SJacob Faibussowitsch {
2587cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
2597cb94ee6SHong Zhang
2607cb94ee6SHong Zhang PetscFunctionBegin;
2619566063dSJacob Faibussowitsch PetscCall(TSReset_Sundials(ts));
262f4f49eeaSPierre Jolivet PetscCallMPI(MPI_Comm_free(&cvode->comm_sundials));
2639566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data));
2649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", NULL));
2659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", NULL));
2669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", NULL));
2679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", NULL));
2689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", NULL));
2699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", NULL));
2709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", NULL));
2719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", NULL));
2729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", NULL));
2739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", NULL));
2742e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", NULL));
2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2767cb94ee6SHong Zhang }
2777cb94ee6SHong Zhang
TSSetUp_Sundials(TS ts)27866976f2fSJacob Faibussowitsch static PetscErrorCode TSSetUp_Sundials(TS ts)
279d71ae5a4SJacob Faibussowitsch {
2807cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
2813ba16761SJacob Faibussowitsch PetscInt glosize, locsize, i;
2827cb94ee6SHong Zhang PetscScalar *y_data, *parray;
283dcbc6d53SSean Farley PC pc;
28419fd82e9SBarry Smith PCType pctype;
285ace3abfcSBarry Smith PetscBool pcnone;
2867cb94ee6SHong Zhang
2877cb94ee6SHong Zhang PetscFunctionBegin;
288bcf0153eSBarry Smith 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");
289236c45dcSShri Abhyankar
2907cb94ee6SHong Zhang /* get the vector size */
2919566063dSJacob Faibussowitsch PetscCall(VecGetSize(ts->vec_sol, &glosize));
2929566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol, &locsize));
2937cb94ee6SHong Zhang
2947cb94ee6SHong Zhang /* allocate the memory for N_Vec y */
295918687b8SPatrick Sanan if (cvode->use_dense) {
296918687b8SPatrick Sanan PetscMPIInt size;
297918687b8SPatrick Sanan
2989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2993c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "TSSUNDIALS only supports a dense solve in the serial case");
300918687b8SPatrick Sanan cvode->y = N_VNew_Serial(locsize);
301918687b8SPatrick Sanan } else {
302a07356b8SSatish Balay cvode->y = N_VNew_Parallel(cvode->comm_sundials, locsize, glosize);
303918687b8SPatrick Sanan }
3043c633725SBarry Smith PetscCheck(cvode->y, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cvode->y is not allocated");
3057cb94ee6SHong Zhang
30628adb3f7SHong Zhang /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
3079566063dSJacob Faibussowitsch PetscCall(VecGetArray(ts->vec_sol, &parray));
3087cb94ee6SHong Zhang y_data = (PetscScalar *)N_VGetArrayPointer(cvode->y);
3097cb94ee6SHong Zhang for (i = 0; i < locsize; i++) y_data[i] = parray[i];
3109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ts->vec_sol, NULL));
31142528757SHong Zhang
3129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &cvode->update));
3139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &cvode->ydot));
3147cb94ee6SHong Zhang
3157cb94ee6SHong Zhang /*
3167cb94ee6SHong Zhang Create work vectors for the TSPSolve_Sundials() routine. Note these are
3177cb94ee6SHong Zhang allocated with zero space arrays because the actual array space is provided
318bcf0153eSBarry Smith by SUNDIALS and set using VecPlaceArray().
3197cb94ee6SHong Zhang */
3209566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w1));
3219566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts), 1, locsize, PETSC_DECIDE, NULL, &cvode->w2));
32216016685SHong Zhang
32316016685SHong Zhang /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
3243ba16761SJacob Faibussowitsch cvode->mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
3253ba16761SJacob Faibussowitsch PetscCheck(cvode->mem, PETSC_COMM_SELF, PETSC_ERR_MEM, "CVodeCreate() fails");
32616016685SHong Zhang
32716016685SHong Zhang /* Set the pointer to user-defined data */
3283ba16761SJacob Faibussowitsch PetscCallExternal(CVodeSetUserData, cvode->mem, ts);
32916016685SHong Zhang
330bcf0153eSBarry Smith /* SUNDIALS may choose to use a smaller initial step, but will never use a larger step. */
3313ba16761SJacob Faibussowitsch PetscCallExternal(CVodeSetInitStep, cvode->mem, (realtype)ts->time_step);
332f1cd61daSJed Brown if (cvode->mindt > 0) {
3333ba16761SJacob Faibussowitsch int flag = CVodeSetMinStep(cvode->mem, (realtype)cvode->mindt);
3349f94935aSHong Zhang if (flag) {
33508401ef6SPierre Jolivet PetscCheck(flag != CV_MEM_NULL, PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed, cvode_mem pointer is NULL");
336f7d195e4SLawrence 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");
337f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_LIB, "CVodeSetMinStep() failed");
3389f94935aSHong Zhang }
339f1cd61daSJed Brown }
3403ba16761SJacob Faibussowitsch if (cvode->maxdt > 0) PetscCallExternal(CVodeSetMaxStep, cvode->mem, (realtype)cvode->maxdt);
341f1cd61daSJed Brown
34216016685SHong Zhang /* Call CVodeInit to initialize the integrator memory and specify the
343dd8e379bSPierre Jolivet * user's right-hand side function in u'=f(t,u), the initial time T0, and
34416016685SHong Zhang * the initial dependent variable vector cvode->y */
3453ba16761SJacob Faibussowitsch PetscCallExternal(CVodeInit, cvode->mem, TSFunction_Sundials, ts->ptime, cvode->y);
34616016685SHong Zhang
3479f94935aSHong Zhang /* specifies scalar relative and absolute tolerances */
3483ba16761SJacob Faibussowitsch PetscCallExternal(CVodeSStolerances, cvode->mem, cvode->reltol, cvode->abstol);
34916016685SHong Zhang
350c4e80e11SFlorian /* Specify max order of BDF / ADAMS method */
3513ba16761SJacob Faibussowitsch if (cvode->maxord != PETSC_DEFAULT) PetscCallExternal(CVodeSetMaxOrd, cvode->mem, cvode->maxord);
352c4e80e11SFlorian
3539f94935aSHong Zhang /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3543ba16761SJacob Faibussowitsch PetscCallExternal(CVodeSetMaxNumSteps, cvode->mem, ts->max_steps);
35516016685SHong Zhang
356918687b8SPatrick Sanan if (cvode->use_dense) {
3573ba16761SJacob Faibussowitsch PetscCallExternal(CVDense, cvode->mem, locsize);
358918687b8SPatrick Sanan } else {
35916016685SHong Zhang /* call CVSpgmr to use GMRES as the linear solver. */
36016016685SHong Zhang /* setup the ode integrator with the given preconditioner */
3619566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc));
3629566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &pctype));
3639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCNONE, &pcnone));
36416016685SHong Zhang if (pcnone) {
3653ba16761SJacob Faibussowitsch PetscCallExternal(CVSpgmr, cvode->mem, PREC_NONE, 0);
36616016685SHong Zhang } else {
3673ba16761SJacob Faibussowitsch PetscCallExternal(CVSpgmr, cvode->mem, PREC_LEFT, cvode->maxl);
36816016685SHong Zhang
36916016685SHong Zhang /* Set preconditioner and solve routines Precond and PSolve,
37016016685SHong Zhang and the pointer to the user-defined block data */
3713ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsSetPreconditioner, cvode->mem, TSPrecond_Sundials_Private, TSPSolve_Sundials_Private);
37216016685SHong Zhang }
373918687b8SPatrick Sanan }
3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3757cb94ee6SHong Zhang }
3767cb94ee6SHong Zhang
3776fadb2cdSHong Zhang /* type of CVODE linear multistep method */
378c793f718SLisandro Dalcin const char *const TSSundialsLmmTypes[] = {"", "ADAMS", "BDF", "TSSundialsLmmType", "SUNDIALS_", NULL};
3796fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
380c793f718SLisandro Dalcin const char *const TSSundialsGramSchmidtTypes[] = {"", "MODIFIED", "CLASSICAL", "TSSundialsGramSchmidtType", "SUNDIALS_", NULL};
381a04cf4d8SBarry Smith
TSSetFromOptions_Sundials(TS ts,PetscOptionItems PetscOptionsObject)382ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Sundials(TS ts, PetscOptionItems PetscOptionsObject)
383d71ae5a4SJacob Faibussowitsch {
3847cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
3857cb94ee6SHong Zhang int indx;
386ace3abfcSBarry Smith PetscBool flag;
3877bda3a07SJed Brown PC pc;
3887cb94ee6SHong Zhang
3897cb94ee6SHong Zhang PetscFunctionBegin;
390d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SUNDIALS ODE solver options");
3919566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_type", "Scheme", "TSSundialsSetType", TSSundialsLmmTypes, 3, TSSundialsLmmTypes[cvode->cvode_type], &indx, &flag));
3921baa6e33SBarry Smith if (flag) PetscCall(TSSundialsSetType(ts, (TSSundialsLmmType)indx));
3939566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_sundials_gramschmidt_type", "Type of orthogonalization", "TSSundialsSetGramSchmidtType", TSSundialsGramSchmidtTypes, 3, TSSundialsGramSchmidtTypes[cvode->gtype], &indx, &flag));
3941baa6e33SBarry Smith if (flag) PetscCall(TSSundialsSetGramSchmidtType(ts, (TSSundialsGramSchmidtType)indx));
3959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_atol", "Absolute tolerance for convergence", "TSSundialsSetTolerance", cvode->abstol, &cvode->abstol, NULL));
3969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_rtol", "Relative tolerance for convergence", "TSSundialsSetTolerance", cvode->reltol, &cvode->reltol, NULL));
3979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_mindt", "Minimum step size", "TSSundialsSetMinTimeStep", cvode->mindt, &cvode->mindt, NULL));
3989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_maxdt", "Maximum step size", "TSSundialsSetMaxTimeStep", cvode->maxdt, &cvode->maxdt, NULL));
3999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_sundials_linear_tolerance", "Convergence tolerance for linear solve", "TSSundialsSetLinearTolerance", cvode->linear_tol, &cvode->linear_tol, NULL));
4009566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxord", "Max Order for BDF/Adams method", "TSSundialsSetMaxOrd", cvode->maxord, &cvode->maxord, NULL));
4019566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_sundials_maxl", "Max dimension of the Krylov subspace", "TSSundialsSetMaxl", cvode->maxl, &cvode->maxl, NULL));
4029566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_monitor_steps", "Monitor SUNDIALS internal steps", "TSSundialsMonitorInternalSteps", cvode->monitorstep, &cvode->monitorstep, NULL));
4039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_sundials_use_dense", "Use dense internal solver in SUNDIALS (serial only)", "TSSundialsSetUseDense", cvode->use_dense, &cvode->use_dense, NULL));
404d0609cedSBarry Smith PetscOptionsHeadEnd();
4059566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc));
4069566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(pc));
4073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4087cb94ee6SHong Zhang }
4097cb94ee6SHong Zhang
TSView_Sundials(TS ts,PetscViewer viewer)41066976f2fSJacob Faibussowitsch static PetscErrorCode TSView_Sundials(TS ts, PetscViewer viewer)
411d71ae5a4SJacob Faibussowitsch {
4127cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
4137cb94ee6SHong Zhang char *type;
4147cb94ee6SHong Zhang char atype[] = "Adams";
4157cb94ee6SHong Zhang char btype[] = "BDF: backward differentiation formula";
4169f196a02SMartin Diehl PetscBool isascii, isstring;
4172c823083SHong Zhang long int nsteps, its, nfevals, nlinsetups, nfails, itmp;
4182c823083SHong Zhang PetscInt qlast, qcur;
4195d47aee6SHong Zhang PetscReal hinused, hlast, hcur, tcur, tolsfac;
42042528757SHong Zhang PC pc;
4217cb94ee6SHong Zhang
4227cb94ee6SHong Zhang PetscFunctionBegin;
423bbd56ea5SKarl Rupp if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
424bbd56ea5SKarl Rupp else type = btype;
4257cb94ee6SHong Zhang
4269f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
4279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
4289f196a02SMartin Diehl if (isascii) {
429bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS integrator does not use SNES!\n"));
430bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS integrator type %s\n", type));
431bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS maxord %" PetscInt_FMT "\n", cvode->maxord));
432bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS abs tol %g rel tol %g\n", (double)cvode->abstol, (double)cvode->reltol));
433918687b8SPatrick Sanan if (cvode->use_dense) {
434bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS integrator using a dense linear solve\n"));
435918687b8SPatrick Sanan } else {
436bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS linear solver tolerance factor %g\n", (double)cvode->linear_tol));
437bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS max dimension of Krylov subspace %" PetscInt_FMT "\n", cvode->maxl));
4387cb94ee6SHong Zhang if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
439bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS using modified Gram-Schmidt for orthogonalization in GMRES\n"));
4407cb94ee6SHong Zhang } else {
441bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n"));
4427cb94ee6SHong Zhang }
443918687b8SPatrick Sanan }
444bcf0153eSBarry Smith if (cvode->mindt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS minimum time step %g\n", (double)cvode->mindt));
445bcf0153eSBarry Smith if (cvode->maxdt > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS maximum time step %g\n", (double)cvode->maxdt));
4462c823083SHong Zhang
4475d47aee6SHong Zhang /* Outputs from CVODE, CVSPILS */
4483ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetTolScaleFactor, cvode->mem, &tolsfac);
449bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS suggested factor for tolerance scaling %g\n", tolsfac));
4503ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetIntegratorStats, cvode->mem, &nsteps, &nfevals, &nlinsetups, &nfails, &qlast, &qcur, &hinused, &hlast, &hcur, &tcur);
451bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS cumulative number of internal steps %ld\n", nsteps));
452bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of calls to rhs function %ld\n", nfevals));
453bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of calls to linear solver setup function %ld\n", nlinsetups));
454bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of error test failures %ld\n", nfails));
4552c823083SHong Zhang
4563ba16761SJacob Faibussowitsch PetscCallExternal(CVodeGetNonlinSolvStats, cvode->mem, &its, &nfails);
457bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of nonlinear solver iterations %ld\n", its));
458bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of nonlinear convergence failure %ld\n", nfails));
459918687b8SPatrick Sanan if (!cvode->use_dense) {
4603ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumLinIters, cvode->mem, &its); /* its = no. of calls to TSPrecond_Sundials() */
461bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of linear iterations %ld\n", its));
4623ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumConvFails, cvode->mem, &itmp);
463bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of linear convergence failures %ld\n", itmp));
46442528757SHong Zhang
4659566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC(ts, &pc));
4669566063dSJacob Faibussowitsch PetscCall(PCView(pc, viewer));
4673ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumPrecEvals, cvode->mem, &itmp);
468bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of preconditioner evaluations %ld\n", itmp));
4693ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumPrecSolves, cvode->mem, &itmp);
470bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of preconditioner solves %ld\n", itmp));
471918687b8SPatrick Sanan }
4723ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumJtimesEvals, cvode->mem, &itmp);
473bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of Jacobian-vector product evaluations %ld\n", itmp));
4743ba16761SJacob Faibussowitsch PetscCallExternal(CVSpilsGetNumRhsEvals, cvode->mem, &itmp);
475bcf0153eSBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "SUNDIALS no. of rhs calls for finite diff. Jacobian-vector evals %ld\n", itmp));
4767cb94ee6SHong Zhang } else if (isstring) {
477bcf0153eSBarry Smith PetscCall(PetscViewerStringSPrintf(viewer, "SUNDIALS type %s", type));
4787cb94ee6SHong Zhang }
4793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4807cb94ee6SHong Zhang }
4817cb94ee6SHong Zhang
TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)48266976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetType_Sundials(TS ts, TSSundialsLmmType type)
483d71ae5a4SJacob Faibussowitsch {
4847cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
4857cb94ee6SHong Zhang
4867cb94ee6SHong Zhang PetscFunctionBegin;
4877cb94ee6SHong Zhang cvode->cvode_type = type;
4883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4897cb94ee6SHong Zhang }
4907cb94ee6SHong Zhang
TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)49166976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts, PetscInt maxl)
492d71ae5a4SJacob Faibussowitsch {
4937cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
4947cb94ee6SHong Zhang
4957cb94ee6SHong Zhang PetscFunctionBegin;
496f61b2b6aSHong Zhang cvode->maxl = maxl;
4973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4987cb94ee6SHong Zhang }
4997cb94ee6SHong Zhang
TSSundialsSetLinearTolerance_Sundials(TS ts,PetscReal tol)50066976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts, PetscReal tol)
501d71ae5a4SJacob Faibussowitsch {
5027cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
5037cb94ee6SHong Zhang
5047cb94ee6SHong Zhang PetscFunctionBegin;
5057cb94ee6SHong Zhang cvode->linear_tol = tol;
5063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5077cb94ee6SHong Zhang }
5087cb94ee6SHong Zhang
TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)50966976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts, TSSundialsGramSchmidtType type)
510d71ae5a4SJacob Faibussowitsch {
5117cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
5127cb94ee6SHong Zhang
5137cb94ee6SHong Zhang PetscFunctionBegin;
5147cb94ee6SHong Zhang cvode->gtype = type;
5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5167cb94ee6SHong Zhang }
5177cb94ee6SHong Zhang
TSSundialsSetTolerance_Sundials(TS ts,PetscReal aabs,PetscReal rel)51866976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts, PetscReal aabs, PetscReal rel)
519d71ae5a4SJacob Faibussowitsch {
5207cb94ee6SHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
5217cb94ee6SHong Zhang
5227cb94ee6SHong Zhang PetscFunctionBegin;
5237cb94ee6SHong Zhang if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5247cb94ee6SHong Zhang if (rel != PETSC_DECIDE) cvode->reltol = rel;
5253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5267cb94ee6SHong Zhang }
5277cb94ee6SHong Zhang
TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)52866976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts, PetscReal mindt)
529d71ae5a4SJacob Faibussowitsch {
530f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials *)ts->data;
531f1cd61daSJed Brown
532f1cd61daSJed Brown PetscFunctionBegin;
533f1cd61daSJed Brown cvode->mindt = mindt;
5343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
535f1cd61daSJed Brown }
536f1cd61daSJed Brown
TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)53766976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts, PetscReal maxdt)
538d71ae5a4SJacob Faibussowitsch {
539f1cd61daSJed Brown TS_Sundials *cvode = (TS_Sundials *)ts->data;
540f1cd61daSJed Brown
541f1cd61daSJed Brown PetscFunctionBegin;
542f1cd61daSJed Brown cvode->maxdt = maxdt;
5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
544f1cd61daSJed Brown }
545918687b8SPatrick Sanan
TSSundialsSetUseDense_Sundials(TS ts,PetscBool use_dense)54666976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsSetUseDense_Sundials(TS ts, PetscBool use_dense)
547d71ae5a4SJacob Faibussowitsch {
548918687b8SPatrick Sanan TS_Sundials *cvode = (TS_Sundials *)ts->data;
549918687b8SPatrick Sanan
550918687b8SPatrick Sanan PetscFunctionBegin;
551918687b8SPatrick Sanan cvode->use_dense = use_dense;
5523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
553918687b8SPatrick Sanan }
554918687b8SPatrick Sanan
TSSundialsGetPC_Sundials(TS ts,PC * pc)55566976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsGetPC_Sundials(TS ts, PC *pc)
556d71ae5a4SJacob Faibussowitsch {
557dcbc6d53SSean Farley SNES snes;
558dcbc6d53SSean Farley KSP ksp;
5597cb94ee6SHong Zhang
5607cb94ee6SHong Zhang PetscFunctionBegin;
5619566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
5629566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp));
5639566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, pc));
5643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5657cb94ee6SHong Zhang }
5667cb94ee6SHong Zhang
TSSundialsGetIterations_Sundials(TS ts,int * nonlin,int * lin)56766976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsGetIterations_Sundials(TS ts, int *nonlin, int *lin)
568d71ae5a4SJacob Faibussowitsch {
5697cb94ee6SHong Zhang PetscFunctionBegin;
5705ef26d82SJed Brown if (nonlin) *nonlin = ts->snes_its;
5715ef26d82SJed Brown if (lin) *lin = ts->ksp_its;
5723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5737cb94ee6SHong Zhang }
5747cb94ee6SHong Zhang
TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)57566976f2fSJacob Faibussowitsch static PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts, PetscBool s)
576d71ae5a4SJacob Faibussowitsch {
5772bfc04deSHong Zhang TS_Sundials *cvode = (TS_Sundials *)ts->data;
5782bfc04deSHong Zhang
5792bfc04deSHong Zhang PetscFunctionBegin;
5802bfc04deSHong Zhang cvode->monitorstep = s;
5813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5822bfc04deSHong Zhang }
5837cb94ee6SHong Zhang
584cc4c1da9SBarry Smith /*@
585bcf0153eSBarry Smith TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by `TSSUNDIALS`.
5867cb94ee6SHong Zhang
5877cb94ee6SHong Zhang Not Collective
5887cb94ee6SHong Zhang
58997bb3fdcSJose E. Roman Input Parameter:
5907cb94ee6SHong Zhang . ts - the time-step context
5917cb94ee6SHong Zhang
5927cb94ee6SHong Zhang Output Parameters:
5937cb94ee6SHong Zhang + nonlin - number of nonlinear iterations
5947cb94ee6SHong Zhang - lin - number of linear iterations
5957cb94ee6SHong Zhang
5967cb94ee6SHong Zhang Level: advanced
5977cb94ee6SHong Zhang
598bcf0153eSBarry Smith Note:
599bcf0153eSBarry Smith These return the number since the creation of the `TS` object
6007cb94ee6SHong Zhang
6011cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
602db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
60342747ad1SJacob Faibussowitsch `TSSundialsGetPC()`, `TSSetExactFinalTime()`
6047cb94ee6SHong Zhang @*/
TSSundialsGetIterations(TS ts,int * nonlin,int * lin)605d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsGetIterations(TS ts, int *nonlin, int *lin)
606d71ae5a4SJacob Faibussowitsch {
6077cb94ee6SHong Zhang PetscFunctionBegin;
608cac4c232SBarry Smith PetscUseMethod(ts, "TSSundialsGetIterations_C", (TS, int *, int *), (ts, nonlin, lin));
6093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6107cb94ee6SHong Zhang }
6117cb94ee6SHong Zhang
6127cb94ee6SHong Zhang /*@
613bcf0153eSBarry Smith TSSundialsSetType - Sets the method that `TSSUNDIALS` will use for integration.
6147cb94ee6SHong Zhang
615c3339decSBarry Smith Logically Collective
6167cb94ee6SHong Zhang
6178f7d6fe5SPatrick Sanan Input Parameters:
6187cb94ee6SHong Zhang + ts - the time-step context
619bcf0153eSBarry Smith - type - one of `SUNDIALS_ADAMS` or `SUNDIALS_BDF`
6207cb94ee6SHong Zhang
6217cb94ee6SHong Zhang Level: intermediate
6227cb94ee6SHong Zhang
6231cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`,
62442747ad1SJacob Faibussowitsch `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`,
62542747ad1SJacob Faibussowitsch `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()`
6267cb94ee6SHong Zhang @*/
TSSundialsSetType(TS ts,TSSundialsLmmType type)627d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetType(TS ts, TSSundialsLmmType type)
628d71ae5a4SJacob Faibussowitsch {
6297cb94ee6SHong Zhang PetscFunctionBegin;
630cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetType_C", (TS, TSSundialsLmmType), (ts, type));
6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6327cb94ee6SHong Zhang }
6337cb94ee6SHong Zhang
6347cb94ee6SHong Zhang /*@
635bcf0153eSBarry Smith TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by `TSSUNDIALS`.
636c4e80e11SFlorian
637c3339decSBarry Smith Logically Collective
638c4e80e11SFlorian
6398f7d6fe5SPatrick Sanan Input Parameters:
640c4e80e11SFlorian + ts - the time-step context
641c4e80e11SFlorian - maxord - maximum order of BDF / Adams method
642c4e80e11SFlorian
643c4e80e11SFlorian Level: advanced
644c4e80e11SFlorian
6451cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`,
64642747ad1SJacob Faibussowitsch `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`,
64742747ad1SJacob Faibussowitsch `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()`
648c4e80e11SFlorian @*/
TSSundialsSetMaxord(TS ts,PetscInt maxord)649d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxord(TS ts, PetscInt maxord)
650d71ae5a4SJacob Faibussowitsch {
651c4e80e11SFlorian PetscFunctionBegin;
652c4e80e11SFlorian PetscValidLogicalCollectiveInt(ts, maxord, 2);
653cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxOrd_C", (TS, PetscInt), (ts, maxord));
6543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
655c4e80e11SFlorian }
656c4e80e11SFlorian
657c4e80e11SFlorian /*@
658f61b2b6aSHong Zhang TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
659bcf0153eSBarry Smith GMRES in the linear solver in `TSSUNDIALS`. `TSSUNDIALS` DOES NOT use restarted GMRES so
660f61b2b6aSHong Zhang this is the maximum number of GMRES steps that will be used.
6617cb94ee6SHong Zhang
662c3339decSBarry Smith Logically Collective
6637cb94ee6SHong Zhang
6648f7d6fe5SPatrick Sanan Input Parameters:
6657cb94ee6SHong Zhang + ts - the time-step context
666f61b2b6aSHong Zhang - maxl - number of direction vectors (the dimension of Krylov subspace).
6677cb94ee6SHong Zhang
6687cb94ee6SHong Zhang Level: advanced
6697cb94ee6SHong Zhang
6701cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`,
671db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
672b43aa488SJacob Faibussowitsch `TSSundialsGetPC()`, `TSSetExactFinalTime()`
6737cb94ee6SHong Zhang @*/
TSSundialsSetMaxl(TS ts,PetscInt maxl)674d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxl(TS ts, PetscInt maxl)
675d71ae5a4SJacob Faibussowitsch {
6767cb94ee6SHong Zhang PetscFunctionBegin;
677f61b2b6aSHong Zhang PetscValidLogicalCollectiveInt(ts, maxl, 2);
678cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxl_C", (TS, PetscInt), (ts, maxl));
6793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6807cb94ee6SHong Zhang }
6817cb94ee6SHong Zhang
6827cb94ee6SHong Zhang /*@
6837cb94ee6SHong Zhang TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
684bcf0153eSBarry Smith system by `TSSUNDIALS`.
6857cb94ee6SHong Zhang
686c3339decSBarry Smith Logically Collective
6877cb94ee6SHong Zhang
6888f7d6fe5SPatrick Sanan Input Parameters:
6897cb94ee6SHong Zhang + ts - the time-step context
6907cb94ee6SHong Zhang - tol - the factor by which the tolerance on the nonlinear solver is
6917cb94ee6SHong Zhang multiplied to get the tolerance on the linear solver, .05 by default.
6927cb94ee6SHong Zhang
6937cb94ee6SHong Zhang Level: advanced
6947cb94ee6SHong Zhang
6951cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
696db781477SPatrick Sanan `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
69742747ad1SJacob Faibussowitsch `TSSundialsGetPC()`,
698db781477SPatrick Sanan `TSSetExactFinalTime()`
6997cb94ee6SHong Zhang @*/
TSSundialsSetLinearTolerance(TS ts,PetscReal tol)700d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetLinearTolerance(TS ts, PetscReal tol)
701d71ae5a4SJacob Faibussowitsch {
7027cb94ee6SHong Zhang PetscFunctionBegin;
703c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts, tol, 2);
704cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetLinearTolerance_C", (TS, PetscReal), (ts, tol));
7053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7067cb94ee6SHong Zhang }
7077cb94ee6SHong Zhang
7087cb94ee6SHong Zhang /*@
7097cb94ee6SHong Zhang TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
710bcf0153eSBarry Smith in GMRES method by `TSSUNDIALS` linear solver.
7117cb94ee6SHong Zhang
712c3339decSBarry Smith Logically Collective
7137cb94ee6SHong Zhang
7148f7d6fe5SPatrick Sanan Input Parameters:
7157cb94ee6SHong Zhang + ts - the time-step context
716bcf0153eSBarry Smith - type - either `SUNDIALS_MODIFIED_GS` or `SUNDIALS_CLASSICAL_GS`
7177cb94ee6SHong Zhang
7187cb94ee6SHong Zhang Level: advanced
7197cb94ee6SHong Zhang
7201cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
721db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`,
722b43aa488SJacob Faibussowitsch `TSSundialsGetPC()`,
723db781477SPatrick Sanan `TSSetExactFinalTime()`
7247cb94ee6SHong Zhang @*/
TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)725d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetGramSchmidtType(TS ts, TSSundialsGramSchmidtType type)
726d71ae5a4SJacob Faibussowitsch {
7277cb94ee6SHong Zhang PetscFunctionBegin;
728cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetGramSchmidtType_C", (TS, TSSundialsGramSchmidtType), (ts, type));
7293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7307cb94ee6SHong Zhang }
7317cb94ee6SHong Zhang
7327cb94ee6SHong Zhang /*@
7337cb94ee6SHong Zhang TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
734bcf0153eSBarry Smith `TSSUNDIALS` for error control.
7357cb94ee6SHong Zhang
736c3339decSBarry Smith Logically Collective
7377cb94ee6SHong Zhang
7388f7d6fe5SPatrick Sanan Input Parameters:
7397cb94ee6SHong Zhang + ts - the time-step context
7407cb94ee6SHong Zhang . aabs - the absolute tolerance
7417cb94ee6SHong Zhang - rel - the relative tolerance
7427cb94ee6SHong Zhang
743bcf0153eSBarry Smith See the CVODE/SUNDIALS users manual for exact details on these parameters. Essentially
7447cb94ee6SHong Zhang these regulate the size of the error for a SINGLE timestep.
7457cb94ee6SHong Zhang
7467cb94ee6SHong Zhang Level: intermediate
7477cb94ee6SHong Zhang
7481cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`,
749db781477SPatrick Sanan `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`,
75042747ad1SJacob Faibussowitsch `TSSundialsGetPC()`,
751db781477SPatrick Sanan `TSSetExactFinalTime()`
7527cb94ee6SHong Zhang @*/
TSSundialsSetTolerance(TS ts,PetscReal aabs,PetscReal rel)753d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetTolerance(TS ts, PetscReal aabs, PetscReal rel)
754d71ae5a4SJacob Faibussowitsch {
7557cb94ee6SHong Zhang PetscFunctionBegin;
756cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetTolerance_C", (TS, PetscReal, PetscReal), (ts, aabs, rel));
7573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7587cb94ee6SHong Zhang }
7597cb94ee6SHong Zhang
7607cb94ee6SHong Zhang /*@
761bcf0153eSBarry Smith TSSundialsGetPC - Extract the PC context from a time-step context for `TSSUNDIALS`.
7627cb94ee6SHong Zhang
7637cb94ee6SHong Zhang Input Parameter:
7647cb94ee6SHong Zhang . ts - the time-step context
7657cb94ee6SHong Zhang
7667cb94ee6SHong Zhang Output Parameter:
7677cb94ee6SHong Zhang . pc - the preconditioner context
7687cb94ee6SHong Zhang
7697cb94ee6SHong Zhang Level: advanced
7707cb94ee6SHong Zhang
7711cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
772b43aa488SJacob Faibussowitsch `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`
7737cb94ee6SHong Zhang @*/
TSSundialsGetPC(TS ts,PC * pc)774d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsGetPC(TS ts, PC *pc)
775d71ae5a4SJacob Faibussowitsch {
7767cb94ee6SHong Zhang PetscFunctionBegin;
777cac4c232SBarry Smith PetscUseMethod(ts, "TSSundialsGetPC_C", (TS, PC *), (ts, pc));
7783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7797cb94ee6SHong Zhang }
7807cb94ee6SHong Zhang
781f1cd61daSJed Brown /*@
782f1cd61daSJed Brown TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
783f1cd61daSJed Brown
784d8d19677SJose E. Roman Input Parameters:
785f1cd61daSJed Brown + ts - the time-step context
786f1cd61daSJed Brown - mindt - lowest time step if positive, negative to deactivate
787f1cd61daSJed Brown
788fc6b9e64SJed Brown Note:
789bcf0153eSBarry Smith `TSSUNDIALS` will error if it is not possible to keep the estimated truncation error below
790bcf0153eSBarry Smith the tolerance set with `TSSundialsSetTolerance()` without going below this step size.
791fc6b9e64SJed Brown
792f1cd61daSJed Brown Level: beginner
793f1cd61daSJed Brown
7941cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsSetType()`, `TSSundialsSetTolerance()`,
795f1cd61daSJed Brown @*/
TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)796d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMinTimeStep(TS ts, PetscReal mindt)
797d71ae5a4SJacob Faibussowitsch {
798f1cd61daSJed Brown PetscFunctionBegin;
799cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMinTimeStep_C", (TS, PetscReal), (ts, mindt));
8003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
801f1cd61daSJed Brown }
802f1cd61daSJed Brown
803f1cd61daSJed Brown /*@
804f1cd61daSJed Brown TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
805f1cd61daSJed Brown
806d8d19677SJose E. Roman Input Parameters:
807f1cd61daSJed Brown + ts - the time-step context
808f1cd61daSJed Brown - maxdt - lowest time step if positive, negative to deactivate
809f1cd61daSJed Brown
810f1cd61daSJed Brown Level: beginner
811f1cd61daSJed Brown
8121cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsSetType()`, `TSSundialsSetTolerance()`,
813f1cd61daSJed Brown @*/
TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)814d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetMaxTimeStep(TS ts, PetscReal maxdt)
815d71ae5a4SJacob Faibussowitsch {
816f1cd61daSJed Brown PetscFunctionBegin;
817cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetMaxTimeStep_C", (TS, PetscReal), (ts, maxdt));
8183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
819f1cd61daSJed Brown }
820f1cd61daSJed Brown
8212bfc04deSHong Zhang /*@
822bcf0153eSBarry Smith TSSundialsMonitorInternalSteps - Monitor `TSSUNDIALS` internal steps (Defaults to false).
8232bfc04deSHong Zhang
824d8d19677SJose E. Roman Input Parameters:
8252bfc04deSHong Zhang + ts - the time-step context
826bcf0153eSBarry Smith - ft - `PETSC_TRUE` if monitor, else `PETSC_FALSE`
8272bfc04deSHong Zhang
8282bfc04deSHong Zhang Level: beginner
8292bfc04deSHong Zhang
8301cc06b55SBarry Smith .seealso: [](ch_ts), `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
831bcf0153eSBarry Smith `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
832b43aa488SJacob Faibussowitsch `TSSundialsGetPC()`
8332bfc04deSHong Zhang @*/
TSSundialsMonitorInternalSteps(TS ts,PetscBool ft)834d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsMonitorInternalSteps(TS ts, PetscBool ft)
835d71ae5a4SJacob Faibussowitsch {
8362bfc04deSHong Zhang PetscFunctionBegin;
837cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsMonitorInternalSteps_C", (TS, PetscBool), (ts, ft));
8383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8392bfc04deSHong Zhang }
840918687b8SPatrick Sanan
841918687b8SPatrick Sanan /*@
842bcf0153eSBarry Smith TSSundialsSetUseDense - Set a flag to use a dense linear solver in `TSSUNDIALS` (serial only)
843918687b8SPatrick Sanan
844918687b8SPatrick Sanan Logically Collective
845918687b8SPatrick Sanan
846918687b8SPatrick Sanan Input Parameters:
847918687b8SPatrick Sanan + ts - the time-step context
848bcf0153eSBarry Smith - use_dense - `PETSC_TRUE` to use the dense solver
849918687b8SPatrick Sanan
850918687b8SPatrick Sanan Level: advanced
851918687b8SPatrick Sanan
8521cc06b55SBarry Smith .seealso: [](ch_ts), `TSSUNDIALS`
853918687b8SPatrick Sanan @*/
TSSundialsSetUseDense(TS ts,PetscBool use_dense)854d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSundialsSetUseDense(TS ts, PetscBool use_dense)
855d71ae5a4SJacob Faibussowitsch {
856918687b8SPatrick Sanan PetscFunctionBegin;
857918687b8SPatrick Sanan PetscValidLogicalCollectiveInt(ts, use_dense, 2);
858cac4c232SBarry Smith PetscTryMethod(ts, "TSSundialsSetUseDense_C", (TS, PetscBool), (ts, use_dense));
8593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
860918687b8SPatrick Sanan }
861918687b8SPatrick Sanan
8627cb94ee6SHong Zhang /*MC
863bcf0153eSBarry Smith TSSUNDIALS - ODE solver using a very old version of the LLNL CVODE/SUNDIALS package, version 2.5 (now called SUNDIALS). Requires ./configure --download-sundials
8647cb94ee6SHong Zhang
865bcf0153eSBarry Smith Options Database Keys:
866897c9f78SHong Zhang + -ts_sundials_type <bdf,adams> -
8677cb94ee6SHong Zhang . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
8687cb94ee6SHong Zhang . -ts_sundials_atol <tol> - Absolute tolerance for convergence
8697cb94ee6SHong Zhang . -ts_sundials_rtol <tol> - Relative tolerance for convergence
870897c9f78SHong Zhang . -ts_sundials_linear_tolerance <tol> -
871f61b2b6aSHong Zhang . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
872918687b8SPatrick Sanan . -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps
873918687b8SPatrick Sanan - -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only)
87416016685SHong Zhang
8757cb94ee6SHong Zhang Level: beginner
8767cb94ee6SHong Zhang
877bcf0153eSBarry Smith Note:
878bcf0153eSBarry Smith This uses its own nonlinear solver and Krylov method so PETSc `SNES` and `KSP` options do not apply,
879bcf0153eSBarry Smith only PETSc `PC` options.
8807cb94ee6SHong Zhang
8811cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`, `TSType`,
882bcf0153eSBarry Smith `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()`
8837cb94ee6SHong Zhang M*/
TSCreate_Sundials(TS ts)884d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts)
885d71ae5a4SJacob Faibussowitsch {
8867cb94ee6SHong Zhang TS_Sundials *cvode;
88742528757SHong Zhang PC pc;
8887cb94ee6SHong Zhang
8897cb94ee6SHong Zhang PetscFunctionBegin;
890277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Sundials;
89128adb3f7SHong Zhang ts->ops->destroy = TSDestroy_Sundials;
89228adb3f7SHong Zhang ts->ops->view = TSView_Sundials;
893214bc6a2SJed Brown ts->ops->setup = TSSetUp_Sundials;
894b4eba00bSSean Farley ts->ops->step = TSStep_Sundials;
895b4eba00bSSean Farley ts->ops->interpolate = TSInterpolate_Sundials;
896214bc6a2SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Sundials;
8972ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE;
8987cb94ee6SHong Zhang
8994dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cvode));
900bbd56ea5SKarl Rupp
901efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE;
902efd4aadfSBarry Smith
9037cb94ee6SHong Zhang ts->data = (void *)cvode;
9046fadb2cdSHong Zhang cvode->cvode_type = SUNDIALS_BDF;
9056fadb2cdSHong Zhang cvode->gtype = SUNDIALS_CLASSICAL_GS;
906f61b2b6aSHong Zhang cvode->maxl = 5;
907c4e80e11SFlorian cvode->maxord = PETSC_DEFAULT;
9087cb94ee6SHong Zhang cvode->linear_tol = .05;
909b4eba00bSSean Farley cvode->monitorstep = PETSC_TRUE;
910918687b8SPatrick Sanan cvode->use_dense = PETSC_FALSE;
9117cb94ee6SHong Zhang
912f4f49eeaSPierre Jolivet PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts), &cvode->comm_sundials));
913f1cd61daSJed Brown
914f1cd61daSJed Brown cvode->mindt = -1.;
915f1cd61daSJed Brown cvode->maxdt = -1.;
916f1cd61daSJed Brown
917bcf0153eSBarry Smith /* set tolerance for SUNDIALS */
9187cb94ee6SHong Zhang cvode->reltol = 1e-6;
9192c823083SHong Zhang cvode->abstol = 1e-6;
9207cb94ee6SHong Zhang
92142528757SHong Zhang /* set PCNONE as default pctype */
9229566063dSJacob Faibussowitsch PetscCall(TSSundialsGetPC_Sundials(ts, &pc));
9239566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE));
92442528757SHong Zhang
9259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetType_C", TSSundialsSetType_Sundials));
9269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxl_C", TSSundialsSetMaxl_Sundials));
9279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetLinearTolerance_C", TSSundialsSetLinearTolerance_Sundials));
9289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetGramSchmidtType_C", TSSundialsSetGramSchmidtType_Sundials));
9299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetTolerance_C", TSSundialsSetTolerance_Sundials));
9309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMinTimeStep_C", TSSundialsSetMinTimeStep_Sundials));
9319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetMaxTimeStep_C", TSSundialsSetMaxTimeStep_Sundials));
9329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetPC_C", TSSundialsGetPC_Sundials));
9339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsGetIterations_C", TSSundialsGetIterations_Sundials));
9349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsMonitorInternalSteps_C", TSSundialsMonitorInternalSteps_Sundials));
9359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSundialsSetUseDense_C", TSSundialsSetUseDense_Sundials));
9363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9377cb94ee6SHong Zhang }
938