xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 691b26d3a7ce3263bd9be9c446af0af2a46feecf)
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 */
14bbd56ea5SKarl Rupp PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,booleantype jok,booleantype *jcurPtr,
15bbd56ea5SKarl Rupp                                   realtype _gamma,void *P_data,N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
167cb94ee6SHong Zhang {
177cb94ee6SHong Zhang   TS             ts     = (TS) P_data;
187cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
19dcbc6d53SSean Farley   PC             pc;
207cb94ee6SHong Zhang   PetscErrorCode ierr;
210679a0aeSJed Brown   Mat            J,P;
220679a0aeSJed Brown   Vec            yy  = cvode->w1,yydot = cvode->ydot;
230679a0aeSJed Brown   PetscReal      gm  = (PetscReal)_gamma;
247cb94ee6SHong Zhang   PetscScalar    *y_data;
257cb94ee6SHong Zhang 
267cb94ee6SHong Zhang   PetscFunctionBegin;
270298fd71SBarry Smith   ierr   = TSGetIJacobian(ts,&J,&P,NULL,NULL);CHKERRQ(ierr);
287cb94ee6SHong Zhang   y_data = (PetscScalar*) N_VGetArrayPointer(y);
297cb94ee6SHong Zhang   ierr   = VecPlaceArray(yy,y_data);CHKERRQ(ierr);
300679a0aeSJed Brown   ierr   = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
310679a0aeSJed Brown   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
3209e82e2bSBarry Smith   ierr     = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,J,P,PETSC_FALSE);CHKERRQ(ierr);
337cb94ee6SHong Zhang   ierr     = VecResetArray(yy);CHKERRQ(ierr);
340679a0aeSJed Brown   ierr     = MatScale(P,gm);CHKERRQ(ierr); /* turn into I-gm*Jrest, J is not used by Sundials  */
357cb94ee6SHong Zhang   *jcurPtr = TRUE;
36dcbc6d53SSean Farley   ierr     = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
3709e82e2bSBarry Smith   ierr     = PCSetOperators(pc,J,P);CHKERRQ(ierr);
387cb94ee6SHong Zhang   PetscFunctionReturn(0);
397cb94ee6SHong Zhang }
407cb94ee6SHong Zhang 
417cb94ee6SHong Zhang /*
427cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
437cb94ee6SHong Zhang */
444ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
454ac7836bSHong Zhang                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
467cb94ee6SHong Zhang {
477cb94ee6SHong Zhang   TS             ts     = (TS) P_data;
487cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
49dcbc6d53SSean Farley   PC             pc;
507cb94ee6SHong Zhang   Vec            rr = cvode->w1,zz = cvode->w2;
517cb94ee6SHong Zhang   PetscErrorCode ierr;
527cb94ee6SHong Zhang   PetscScalar    *r_data,*z_data;
537cb94ee6SHong Zhang 
547cb94ee6SHong Zhang   PetscFunctionBegin;
557cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
567cb94ee6SHong Zhang   r_data = (PetscScalar*) N_VGetArrayPointer(r);
577cb94ee6SHong Zhang   z_data = (PetscScalar*) N_VGetArrayPointer(z);
587cb94ee6SHong Zhang   ierr   = VecPlaceArray(rr,r_data);CHKERRQ(ierr);
597cb94ee6SHong Zhang   ierr   = VecPlaceArray(zz,z_data);CHKERRQ(ierr);
604ac7836bSHong Zhang 
617cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
62dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
637cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz);CHKERRQ(ierr);
647cb94ee6SHong Zhang   ierr = VecResetArray(rr);CHKERRQ(ierr);
657cb94ee6SHong Zhang   ierr = VecResetArray(zz);CHKERRQ(ierr);
667cb94ee6SHong Zhang   PetscFunctionReturn(0);
677cb94ee6SHong Zhang }
687cb94ee6SHong Zhang 
697cb94ee6SHong Zhang /*
707cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
717cb94ee6SHong Zhang */
724ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
737cb94ee6SHong Zhang {
747cb94ee6SHong Zhang   TS             ts = (TS) ctx;
755fcef5e4SJed Brown   DM             dm;
76942e3340SBarry Smith   DMTS           tsdm;
775fcef5e4SJed Brown   TSIFunction    ifunction;
78ce94432eSBarry Smith   MPI_Comm       comm;
797cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
800679a0aeSJed Brown   Vec            yy     = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
817cb94ee6SHong Zhang   PetscScalar    *y_data,*ydot_data;
827cb94ee6SHong Zhang   PetscErrorCode ierr;
837cb94ee6SHong Zhang 
847cb94ee6SHong Zhang   PetscFunctionBegin;
85ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
865ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
877cb94ee6SHong Zhang   y_data    = (PetscScalar*) N_VGetArrayPointer(y);
887cb94ee6SHong Zhang   ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot);
894ac7836bSHong Zhang   ierr      = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
904ac7836bSHong Zhang   ierr      = VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr);
914ac7836bSHong Zhang 
925fcef5e4SJed Brown   /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
935fcef5e4SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
94942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
950298fd71SBarry Smith   ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr);
965fcef5e4SJed Brown   if (!ifunction) {
97bc0cc02bSJed Brown     ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr);
98bc0cc02bSJed Brown   } else {                      /* If rhsfunction is also set, this computes both parts and shifts them to the right */
99bc0cc02bSJed Brown     ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
1000679a0aeSJed Brown     ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr);
1010679a0aeSJed Brown     ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
102bc0cc02bSJed Brown   }
1037cb94ee6SHong Zhang   ierr = VecResetArray(yy);CHKERRABORT(comm,ierr);
1047cb94ee6SHong Zhang   ierr = VecResetArray(yyd);CHKERRABORT(comm,ierr);
1054ac7836bSHong Zhang   PetscFunctionReturn(0);
1067cb94ee6SHong Zhang }
1077cb94ee6SHong Zhang 
1087cb94ee6SHong Zhang /*
109b4eba00bSSean Farley        TSStep_Sundials - Calls Sundials to integrate the ODE.
1107cb94ee6SHong Zhang */
111b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts)
1127cb94ee6SHong Zhang {
1137cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1147cb94ee6SHong Zhang   PetscErrorCode ierr;
115b4eba00bSSean Farley   PetscInt       flag;
116be5899b3SLisandro Dalcin   long int       nits,lits,nsteps;
1177cb94ee6SHong Zhang   realtype       t,tout;
1187cb94ee6SHong Zhang   PetscScalar    *y_data;
1197cb94ee6SHong Zhang   void           *mem;
1207cb94ee6SHong Zhang 
1217cb94ee6SHong Zhang   PetscFunctionBegin;
12216016685SHong Zhang   mem  = cvode->mem;
1237cb94ee6SHong Zhang   tout = ts->max_time;
1247cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1257cb94ee6SHong Zhang   N_VSetArrayPointer((realtype*)y_data,cvode->y);
1260298fd71SBarry Smith   ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr);
127186e87acSLisandro Dalcin 
1289687d888SLisandro Dalcin   /* We would like to TSPreStage() and TSPostStage()
1299687d888SLisandro Dalcin    * before each stage solve but CVode does not appear to support this. */
130be5899b3SLisandro Dalcin   if (cvode->monitorstep)
131be5899b3SLisandro Dalcin     flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
132be5899b3SLisandro Dalcin   else
133be5899b3SLisandro Dalcin     flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
1349f94935aSHong Zhang 
1359f94935aSHong Zhang   if (flag) { /* display error message */
1369f94935aSHong Zhang     switch (flag) {
1379f94935aSHong Zhang       case CV_ILL_INPUT:
1389f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
1399f94935aSHong Zhang         break;
1409f94935aSHong Zhang       case CV_TOO_CLOSE:
1419f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
1429f94935aSHong Zhang         break;
1433c7fefeeSJed Brown       case CV_TOO_MUCH_WORK: {
1449f94935aSHong Zhang         PetscReal tcur;
1459f94935aSHong Zhang         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
1469f94935aSHong Zhang         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
14719eac22cSLisandro Dalcin         SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%g, nsteps %D exceeds maxstep %D. Increase '-ts_max_steps <>' or modify TSSetMaxSteps()",(double)tcur,nsteps,ts->max_steps);
1483c7fefeeSJed Brown       } break;
1499f94935aSHong Zhang       case CV_TOO_MUCH_ACC:
1509f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
1519f94935aSHong Zhang         break;
1529f94935aSHong Zhang       case CV_ERR_FAILURE:
1539f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
1549f94935aSHong Zhang         break;
1559f94935aSHong Zhang       case CV_CONV_FAILURE:
1569f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
1579f94935aSHong Zhang         break;
1589f94935aSHong Zhang       case CV_LINIT_FAIL:
1599f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
1609f94935aSHong Zhang         break;
1619f94935aSHong Zhang       case CV_LSETUP_FAIL:
1629f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
1639f94935aSHong Zhang         break;
1649f94935aSHong Zhang       case CV_LSOLVE_FAIL:
1659f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
1669f94935aSHong Zhang         break;
1679f94935aSHong Zhang       case CV_RHSFUNC_FAIL:
1689f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
1699f94935aSHong Zhang         break;
1709f94935aSHong Zhang       case CV_FIRST_RHSFUNC_ERR:
1719f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
1729f94935aSHong Zhang         break;
1739f94935aSHong Zhang       case CV_REPTD_RHSFUNC_ERR:
1749f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
1759f94935aSHong Zhang         break;
1769f94935aSHong Zhang       case CV_UNREC_RHSFUNC_ERR:
1779f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
1789f94935aSHong Zhang         break;
1799f94935aSHong Zhang       case CV_RTFUNC_FAIL:
1809f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
1819f94935aSHong Zhang         break;
1829f94935aSHong Zhang       default:
1839f94935aSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
1849f94935aSHong Zhang     }
1859f94935aSHong Zhang   }
1869f94935aSHong Zhang 
187be5899b3SLisandro Dalcin   /* log inner nonlinear and linear iterations */
188be5899b3SLisandro Dalcin   ierr = CVodeGetNumNonlinSolvIters(mem,&nits);CHKERRQ(ierr);
189be5899b3SLisandro Dalcin   ierr = CVSpilsGetNumLinIters(mem,&lits);CHKERRQ(ierr);
190be5899b3SLisandro Dalcin   ts->snes_its += nits; ts->ksp_its = lits;
191be5899b3SLisandro Dalcin 
1927cb94ee6SHong Zhang   /* copy the solution from cvode->y to cvode->update and sol */
1937cb94ee6SHong Zhang   ierr = VecPlaceArray(cvode->w1,y_data);CHKERRQ(ierr);
1947cb94ee6SHong Zhang   ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
1957cb94ee6SHong Zhang   ierr = VecResetArray(cvode->w1);CHKERRQ(ierr);
196bc0cc02bSJed Brown   ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr);
197186e87acSLisandro Dalcin 
198186e87acSLisandro Dalcin   ts->time_step = t - ts->ptime;
199186e87acSLisandro Dalcin   ts->ptime = t;
200186e87acSLisandro Dalcin 
2019f94935aSHong Zhang   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
2022808aa04SLisandro Dalcin   if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */
203b4eba00bSSean Farley   PetscFunctionReturn(0);
204b4eba00bSSean Farley }
205b4eba00bSSean Farley 
206b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
207b4eba00bSSean Farley {
208b4eba00bSSean Farley   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
209b4eba00bSSean Farley   N_Vector       y;
210b4eba00bSSean Farley   PetscErrorCode ierr;
211b4eba00bSSean Farley   PetscScalar    *x_data;
212b4eba00bSSean Farley   PetscInt       glosize,locsize;
213b4eba00bSSean Farley 
214b4eba00bSSean Farley   PetscFunctionBegin;
215b4eba00bSSean Farley   /* get the vector size */
216b4eba00bSSean Farley   ierr = VecGetSize(X,&glosize);CHKERRQ(ierr);
217b4eba00bSSean Farley   ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr);
218b4eba00bSSean Farley 
219b4eba00bSSean Farley   /* allocate the memory for N_Vec y */
220b4eba00bSSean Farley   y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
221*691b26d3SBarry Smith   if (!y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Interpolated y is not allocated");
222b4eba00bSSean Farley 
223b4eba00bSSean Farley   ierr = VecGetArray(X,&x_data);CHKERRQ(ierr);
224b4eba00bSSean Farley   N_VSetArrayPointer((realtype*)x_data,y);
225b4eba00bSSean Farley   ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr);
226b4eba00bSSean Farley   ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr);
2277cb94ee6SHong Zhang   PetscFunctionReturn(0);
2287cb94ee6SHong Zhang }
2297cb94ee6SHong Zhang 
230277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts)
231277b19d0SLisandro Dalcin {
232277b19d0SLisandro Dalcin   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
233277b19d0SLisandro Dalcin   PetscErrorCode ierr;
234277b19d0SLisandro Dalcin 
235277b19d0SLisandro Dalcin   PetscFunctionBegin;
2365c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
2370679a0aeSJed Brown   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
2385c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
2395c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
240bbd56ea5SKarl Rupp   if (cvode->mem) CVodeFree(&cvode->mem);
241277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
242277b19d0SLisandro Dalcin }
243277b19d0SLisandro Dalcin 
2447cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
2457cb94ee6SHong Zhang {
2467cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2477cb94ee6SHong Zhang   PetscErrorCode ierr;
2487cb94ee6SHong Zhang 
2497cb94ee6SHong Zhang   PetscFunctionBegin;
250277b19d0SLisandro Dalcin   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
251a07356b8SSatish Balay   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
252277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
253bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL);CHKERRQ(ierr);
254bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL);CHKERRQ(ierr);
255bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL);CHKERRQ(ierr);
256bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL);CHKERRQ(ierr);
257bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL);CHKERRQ(ierr);
258bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL);CHKERRQ(ierr);
259bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL);CHKERRQ(ierr);
260bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL);CHKERRQ(ierr);
261bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL);CHKERRQ(ierr);
262bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL);CHKERRQ(ierr);
2637cb94ee6SHong Zhang   PetscFunctionReturn(0);
2647cb94ee6SHong Zhang }
2657cb94ee6SHong Zhang 
266214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts)
2677cb94ee6SHong Zhang {
2687cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2697cb94ee6SHong Zhang   PetscErrorCode ierr;
27016016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2717cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
27216016685SHong Zhang   void           *mem;
273dcbc6d53SSean Farley   PC             pc;
27419fd82e9SBarry Smith   PCType         pctype;
275ace3abfcSBarry Smith   PetscBool      pcnone;
2767cb94ee6SHong Zhang 
2777cb94ee6SHong Zhang   PetscFunctionBegin;
278236c45dcSShri Abhyankar   if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for exact final time option 'MATCHSTEP' when using Sundials");
279236c45dcSShri Abhyankar 
2807cb94ee6SHong Zhang   /* get the vector size */
2817cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2827cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
2837cb94ee6SHong Zhang 
2847cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
285a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
286*691b26d3SBarry Smith   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cvode->y is not allocated");
2877cb94ee6SHong Zhang 
28828adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
2897cb94ee6SHong Zhang   ierr   = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
2907cb94ee6SHong Zhang   y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y);
2917cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
2920298fd71SBarry Smith   ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr);
29342528757SHong Zhang 
2947cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
2950679a0aeSJed Brown   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
2963bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update);CHKERRQ(ierr);
2973bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot);CHKERRQ(ierr);
2987cb94ee6SHong Zhang 
2997cb94ee6SHong Zhang   /*
3007cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
3017cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
3027cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
3037cb94ee6SHong Zhang   */
304ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
305ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
3063bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1);CHKERRQ(ierr);
3073bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2);CHKERRQ(ierr);
30816016685SHong Zhang 
30916016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
31016016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
311e32f2f54SBarry Smith   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
31216016685SHong Zhang   cvode->mem = mem;
31316016685SHong Zhang 
31416016685SHong Zhang   /* Set the pointer to user-defined data */
31516016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
316e32f2f54SBarry Smith   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
31716016685SHong Zhang 
318fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
319cdbf8f93SLisandro Dalcin   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
320ce94432eSBarry Smith   if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed");
321f1cd61daSJed Brown   if (cvode->mindt > 0) {
322f1cd61daSJed Brown     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
3239f94935aSHong Zhang     if (flag) {
324ce94432eSBarry Smith       if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
325ce94432eSBarry Smith       else if (flag == CV_ILL_INPUT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
326ce94432eSBarry Smith       else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed");
3279f94935aSHong Zhang     }
328f1cd61daSJed Brown   }
329f1cd61daSJed Brown   if (cvode->maxdt > 0) {
330f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
331ce94432eSBarry Smith     if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
332f1cd61daSJed Brown   }
333f1cd61daSJed Brown 
33416016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
33516016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
33616016685SHong Zhang    * the initial dependent variable vector cvode->y */
33716016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
338bbd56ea5SKarl Rupp   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
33916016685SHong Zhang 
3409f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
34116016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
342bbd56ea5SKarl Rupp   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
34316016685SHong Zhang 
3449f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3459f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
346be5899b3SLisandro Dalcin   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag);
34716016685SHong Zhang 
34816016685SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
34916016685SHong Zhang   /* setup the ode integrator with the given preconditioner */
350dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
351dcbc6d53SSean Farley   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
352251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
35316016685SHong Zhang   if (pcnone) {
35416016685SHong Zhang     flag = CVSpgmr(mem,PREC_NONE,0);
355e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
35616016685SHong Zhang   } else {
357f61b2b6aSHong Zhang     flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
358e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
35916016685SHong Zhang 
36016016685SHong Zhang     /* Set preconditioner and solve routines Precond and PSolve,
36116016685SHong Zhang      and the pointer to the user-defined block data */
36216016685SHong Zhang     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
363e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
36416016685SHong Zhang   }
36516016685SHong Zhang 
36616016685SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
367e32f2f54SBarry Smith   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
3687cb94ee6SHong Zhang   PetscFunctionReturn(0);
3697cb94ee6SHong Zhang }
3707cb94ee6SHong Zhang 
3716fadb2cdSHong Zhang /* type of CVODE linear multistep method */
3726a6fc655SJed Brown const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0};
3736fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
3746a6fc655SJed Brown const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0};
375a04cf4d8SBarry Smith 
3764416b707SBarry Smith PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts)
3777cb94ee6SHong Zhang {
3787cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3797cb94ee6SHong Zhang   PetscErrorCode ierr;
3807cb94ee6SHong Zhang   int            indx;
381ace3abfcSBarry Smith   PetscBool      flag;
3827bda3a07SJed Brown   PC             pc;
3837cb94ee6SHong Zhang 
3847cb94ee6SHong Zhang   PetscFunctionBegin;
385e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SUNDIALS ODE solver options");CHKERRQ(ierr);
3866fadb2cdSHong Zhang   ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
3877cb94ee6SHong Zhang   if (flag) {
3886fadb2cdSHong Zhang     ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
3897cb94ee6SHong Zhang   }
3906fadb2cdSHong Zhang   ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
3917cb94ee6SHong Zhang   if (flag) {
3927cb94ee6SHong Zhang     ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
3937cb94ee6SHong Zhang   }
3940298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);CHKERRQ(ierr);
3950298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);CHKERRQ(ierr);
3960298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);CHKERRQ(ierr);
3970298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);CHKERRQ(ierr);
39894ae4db5SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL);CHKERRQ(ierr);
39994ae4db5SBarry Smith   ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL);CHKERRQ(ierr);
400be5899b3SLisandro Dalcin   ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);CHKERRQ(ierr);
4017cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4027bda3a07SJed Brown   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
4037bda3a07SJed Brown   ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
4047cb94ee6SHong Zhang   PetscFunctionReturn(0);
4057cb94ee6SHong Zhang }
4067cb94ee6SHong Zhang 
4077cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
4087cb94ee6SHong Zhang {
4097cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4107cb94ee6SHong Zhang   PetscErrorCode ierr;
4117cb94ee6SHong Zhang   char           *type;
4127cb94ee6SHong Zhang   char           atype[] = "Adams";
4137cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
414ace3abfcSBarry Smith   PetscBool      iascii,isstring;
4152c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
4162c823083SHong Zhang   PetscInt       qlast,qcur;
4175d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
41842528757SHong Zhang   PC             pc;
4197cb94ee6SHong Zhang 
4207cb94ee6SHong Zhang   PetscFunctionBegin;
421bbd56ea5SKarl Rupp   if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
422bbd56ea5SKarl Rupp   else                                     type = btype;
4237cb94ee6SHong Zhang 
424251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
425251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
4267cb94ee6SHong Zhang   if (iascii) {
4277cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
4287cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
4297cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
4307cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
431f61b2b6aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr);
4327cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4337cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4347cb94ee6SHong Zhang     } else {
4357cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4367cb94ee6SHong Zhang     }
437f1cd61daSJed Brown     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
438f1cd61daSJed Brown     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
4392c823083SHong Zhang 
4405d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4415d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
4425d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
4432c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
4442c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
4452c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
4462c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
4472c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
4482c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
4492c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
4502c823083SHong Zhang 
4512c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
4522c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
4532c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
4542c823083SHong Zhang 
4552c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
4562c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
4572c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
4582c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
45942528757SHong Zhang 
46042528757SHong Zhang     ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
46142528757SHong Zhang     ierr = PCView(pc,viewer);CHKERRQ(ierr);
4622c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4632c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
4642c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
4652c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
46642528757SHong Zhang 
4672c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4682c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
4692c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4702c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
4717cb94ee6SHong Zhang   } else if (isstring) {
4727cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
4737cb94ee6SHong Zhang   }
4747cb94ee6SHong Zhang   PetscFunctionReturn(0);
4757cb94ee6SHong Zhang }
4767cb94ee6SHong Zhang 
4777cb94ee6SHong Zhang 
4787cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
4797087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
4807cb94ee6SHong Zhang {
4817cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4827cb94ee6SHong Zhang 
4837cb94ee6SHong Zhang   PetscFunctionBegin;
4847cb94ee6SHong Zhang   cvode->cvode_type = type;
4857cb94ee6SHong Zhang   PetscFunctionReturn(0);
4867cb94ee6SHong Zhang }
4877cb94ee6SHong Zhang 
488f61b2b6aSHong Zhang PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
4897cb94ee6SHong Zhang {
4907cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4917cb94ee6SHong Zhang 
4927cb94ee6SHong Zhang   PetscFunctionBegin;
493f61b2b6aSHong Zhang   cvode->maxl = maxl;
4947cb94ee6SHong Zhang   PetscFunctionReturn(0);
4957cb94ee6SHong Zhang }
4967cb94ee6SHong Zhang 
4977087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
4987cb94ee6SHong Zhang {
4997cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5007cb94ee6SHong Zhang 
5017cb94ee6SHong Zhang   PetscFunctionBegin;
5027cb94ee6SHong Zhang   cvode->linear_tol = tol;
5037cb94ee6SHong Zhang   PetscFunctionReturn(0);
5047cb94ee6SHong Zhang }
5057cb94ee6SHong Zhang 
5067087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
5077cb94ee6SHong Zhang {
5087cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5097cb94ee6SHong Zhang 
5107cb94ee6SHong Zhang   PetscFunctionBegin;
5117cb94ee6SHong Zhang   cvode->gtype = type;
5127cb94ee6SHong Zhang   PetscFunctionReturn(0);
5137cb94ee6SHong Zhang }
5147cb94ee6SHong Zhang 
5157087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
5167cb94ee6SHong Zhang {
5177cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5187cb94ee6SHong Zhang 
5197cb94ee6SHong Zhang   PetscFunctionBegin;
5207cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5217cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
5227cb94ee6SHong Zhang   PetscFunctionReturn(0);
5237cb94ee6SHong Zhang }
5247cb94ee6SHong Zhang 
5257087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
526f1cd61daSJed Brown {
527f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
528f1cd61daSJed Brown 
529f1cd61daSJed Brown   PetscFunctionBegin;
530f1cd61daSJed Brown   cvode->mindt = mindt;
531f1cd61daSJed Brown   PetscFunctionReturn(0);
532f1cd61daSJed Brown }
533f1cd61daSJed Brown 
5347087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
535f1cd61daSJed Brown {
536f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
537f1cd61daSJed Brown 
538f1cd61daSJed Brown   PetscFunctionBegin;
539f1cd61daSJed Brown   cvode->maxdt = maxdt;
540f1cd61daSJed Brown   PetscFunctionReturn(0);
541f1cd61daSJed Brown }
5427087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
5437cb94ee6SHong Zhang {
544dcbc6d53SSean Farley   SNES           snes;
545dcbc6d53SSean Farley   KSP            ksp;
546dcbc6d53SSean Farley   PetscErrorCode ierr;
5477cb94ee6SHong Zhang 
5487cb94ee6SHong Zhang   PetscFunctionBegin;
549dcbc6d53SSean Farley   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
550dcbc6d53SSean Farley   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
551dcbc6d53SSean Farley   ierr = KSPGetPC(ksp,pc);CHKERRQ(ierr);
5527cb94ee6SHong Zhang   PetscFunctionReturn(0);
5537cb94ee6SHong Zhang }
5547cb94ee6SHong Zhang 
5557087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
5567cb94ee6SHong Zhang {
5577cb94ee6SHong Zhang   PetscFunctionBegin;
5585ef26d82SJed Brown   if (nonlin) *nonlin = ts->snes_its;
5595ef26d82SJed Brown   if (lin)    *lin    = ts->ksp_its;
5607cb94ee6SHong Zhang   PetscFunctionReturn(0);
5617cb94ee6SHong Zhang }
5627cb94ee6SHong Zhang 
5637087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)
5642bfc04deSHong Zhang {
5652bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5662bfc04deSHong Zhang 
5672bfc04deSHong Zhang   PetscFunctionBegin;
5682bfc04deSHong Zhang   cvode->monitorstep = s;
5692bfc04deSHong Zhang   PetscFunctionReturn(0);
5702bfc04deSHong Zhang }
5717cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
5727cb94ee6SHong Zhang 
5737cb94ee6SHong Zhang /*@C
5747cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
5757cb94ee6SHong Zhang 
5767cb94ee6SHong Zhang    Not Collective
5777cb94ee6SHong Zhang 
5787cb94ee6SHong Zhang    Input parameters:
5797cb94ee6SHong Zhang .    ts     - the time-step context
5807cb94ee6SHong Zhang 
5817cb94ee6SHong Zhang    Output Parameters:
5827cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
5837cb94ee6SHong Zhang -   lin    - number of linear iterations
5847cb94ee6SHong Zhang 
5857cb94ee6SHong Zhang    Level: advanced
5867cb94ee6SHong Zhang 
5877cb94ee6SHong Zhang    Notes:
5887cb94ee6SHong Zhang     These return the number since the creation of the TS object
5897cb94ee6SHong Zhang 
590f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
5917cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
592f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
593a43b19c4SJed Brown           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
5947cb94ee6SHong Zhang 
5957cb94ee6SHong Zhang @*/
5967087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
5977cb94ee6SHong Zhang {
5984ac538c5SBarry Smith   PetscErrorCode ierr;
5997cb94ee6SHong Zhang 
6007cb94ee6SHong Zhang   PetscFunctionBegin;
6014ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
6027cb94ee6SHong Zhang   PetscFunctionReturn(0);
6037cb94ee6SHong Zhang }
6047cb94ee6SHong Zhang 
6057cb94ee6SHong Zhang /*@
6067cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
6077cb94ee6SHong Zhang 
6083f9fe445SBarry Smith    Logically Collective on TS
6097cb94ee6SHong Zhang 
6107cb94ee6SHong Zhang    Input parameters:
6117cb94ee6SHong Zhang +    ts     - the time-step context
6127cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
6137cb94ee6SHong Zhang 
6147cb94ee6SHong Zhang    Level: intermediate
6157cb94ee6SHong Zhang 
616f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetMaxl(),
6177cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
618f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
6197cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
620a43b19c4SJed Brown           TSSetExactFinalTime()
6217cb94ee6SHong Zhang @*/
6227087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
6237cb94ee6SHong Zhang {
6244ac538c5SBarry Smith   PetscErrorCode ierr;
6257cb94ee6SHong Zhang 
6267cb94ee6SHong Zhang   PetscFunctionBegin;
6274ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
6287cb94ee6SHong Zhang   PetscFunctionReturn(0);
6297cb94ee6SHong Zhang }
6307cb94ee6SHong Zhang 
6317cb94ee6SHong Zhang /*@
632f61b2b6aSHong Zhang    TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
6337cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
634f61b2b6aSHong Zhang        this is the maximum number of GMRES steps that will be used.
6357cb94ee6SHong Zhang 
6363f9fe445SBarry Smith    Logically Collective on TS
6377cb94ee6SHong Zhang 
6387cb94ee6SHong Zhang    Input parameters:
6397cb94ee6SHong Zhang +    ts      - the time-step context
640f61b2b6aSHong Zhang -    maxl - number of direction vectors (the dimension of Krylov subspace).
6417cb94ee6SHong Zhang 
6427cb94ee6SHong Zhang    Level: advanced
6437cb94ee6SHong Zhang 
6447cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
6457cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
646f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
6477cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
648a43b19c4SJed Brown           TSSetExactFinalTime()
6497cb94ee6SHong Zhang 
6507cb94ee6SHong Zhang @*/
651f61b2b6aSHong Zhang PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
6527cb94ee6SHong Zhang {
6534ac538c5SBarry Smith   PetscErrorCode ierr;
6547cb94ee6SHong Zhang 
6557cb94ee6SHong Zhang   PetscFunctionBegin;
656f61b2b6aSHong Zhang   PetscValidLogicalCollectiveInt(ts,maxl,2);
657f61b2b6aSHong Zhang   ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr);
6587cb94ee6SHong Zhang   PetscFunctionReturn(0);
6597cb94ee6SHong Zhang }
6607cb94ee6SHong Zhang 
6617cb94ee6SHong Zhang /*@
6627cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
6637cb94ee6SHong Zhang        system by SUNDIALS.
6647cb94ee6SHong Zhang 
6653f9fe445SBarry Smith    Logically Collective on TS
6667cb94ee6SHong Zhang 
6677cb94ee6SHong Zhang    Input parameters:
6687cb94ee6SHong Zhang +    ts     - the time-step context
6697cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
6707cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
6717cb94ee6SHong Zhang 
6727cb94ee6SHong Zhang    Level: advanced
6737cb94ee6SHong Zhang 
674f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
6757cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
676f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
6777cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
678a43b19c4SJed Brown           TSSetExactFinalTime()
6797cb94ee6SHong Zhang 
6807cb94ee6SHong Zhang @*/
6817087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
6827cb94ee6SHong Zhang {
6834ac538c5SBarry Smith   PetscErrorCode ierr;
6847cb94ee6SHong Zhang 
6857cb94ee6SHong Zhang   PetscFunctionBegin;
686c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,tol,2);
6874ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
6887cb94ee6SHong Zhang   PetscFunctionReturn(0);
6897cb94ee6SHong Zhang }
6907cb94ee6SHong Zhang 
6917cb94ee6SHong Zhang /*@
6927cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
6937cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
6947cb94ee6SHong Zhang 
6953f9fe445SBarry Smith    Logically Collective on TS
6967cb94ee6SHong Zhang 
6977cb94ee6SHong Zhang    Input parameters:
6987cb94ee6SHong Zhang +    ts  - the time-step context
6997cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
7007cb94ee6SHong Zhang 
7017cb94ee6SHong Zhang    Level: advanced
7027cb94ee6SHong Zhang 
703f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
7047cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
705f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7067cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
707a43b19c4SJed Brown           TSSetExactFinalTime()
7087cb94ee6SHong Zhang 
7097cb94ee6SHong Zhang @*/
7107087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
7117cb94ee6SHong Zhang {
7124ac538c5SBarry Smith   PetscErrorCode ierr;
7137cb94ee6SHong Zhang 
7147cb94ee6SHong Zhang   PetscFunctionBegin;
7154ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
7167cb94ee6SHong Zhang   PetscFunctionReturn(0);
7177cb94ee6SHong Zhang }
7187cb94ee6SHong Zhang 
7197cb94ee6SHong Zhang /*@
7207cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
7217cb94ee6SHong Zhang                          Sundials for error control.
7227cb94ee6SHong Zhang 
7233f9fe445SBarry Smith    Logically Collective on TS
7247cb94ee6SHong Zhang 
7257cb94ee6SHong Zhang    Input parameters:
7267cb94ee6SHong Zhang +    ts  - the time-step context
7277cb94ee6SHong Zhang .    aabs - the absolute tolerance
7287cb94ee6SHong Zhang -    rel - the relative tolerance
7297cb94ee6SHong Zhang 
7307cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
7317cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
7327cb94ee6SHong Zhang 
7337cb94ee6SHong Zhang    Level: intermediate
7347cb94ee6SHong Zhang 
735f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
7367cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
737f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7387cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
739a43b19c4SJed Brown           TSSetExactFinalTime()
7407cb94ee6SHong Zhang 
7417cb94ee6SHong Zhang @*/
7427087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
7437cb94ee6SHong Zhang {
7444ac538c5SBarry Smith   PetscErrorCode ierr;
7457cb94ee6SHong Zhang 
7467cb94ee6SHong Zhang   PetscFunctionBegin;
7474ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
7487cb94ee6SHong Zhang   PetscFunctionReturn(0);
7497cb94ee6SHong Zhang }
7507cb94ee6SHong Zhang 
7517cb94ee6SHong Zhang /*@
7527cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
7537cb94ee6SHong Zhang 
7547cb94ee6SHong Zhang    Input Parameter:
7557cb94ee6SHong Zhang .    ts - the time-step context
7567cb94ee6SHong Zhang 
7577cb94ee6SHong Zhang    Output Parameter:
7587cb94ee6SHong Zhang .    pc - the preconditioner context
7597cb94ee6SHong Zhang 
7607cb94ee6SHong Zhang    Level: advanced
7617cb94ee6SHong Zhang 
762f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
7637cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
764f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7657cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
7667cb94ee6SHong Zhang @*/
7677087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
7687cb94ee6SHong Zhang {
7694ac538c5SBarry Smith   PetscErrorCode ierr;
7707cb94ee6SHong Zhang 
7717cb94ee6SHong Zhang   PetscFunctionBegin;
7724ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));CHKERRQ(ierr);
7737cb94ee6SHong Zhang   PetscFunctionReturn(0);
7747cb94ee6SHong Zhang }
7757cb94ee6SHong Zhang 
776f1cd61daSJed Brown /*@
777f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
778f1cd61daSJed Brown 
779f1cd61daSJed Brown    Input Parameter:
780f1cd61daSJed Brown +   ts - the time-step context
781f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
782f1cd61daSJed Brown 
783fc6b9e64SJed Brown    Note:
784fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
785fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
786fc6b9e64SJed Brown 
787f1cd61daSJed Brown    Level: beginner
788f1cd61daSJed Brown 
789f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
790f1cd61daSJed Brown @*/
7917087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
792f1cd61daSJed Brown {
7934ac538c5SBarry Smith   PetscErrorCode ierr;
794f1cd61daSJed Brown 
795f1cd61daSJed Brown   PetscFunctionBegin;
7964ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
797f1cd61daSJed Brown   PetscFunctionReturn(0);
798f1cd61daSJed Brown }
799f1cd61daSJed Brown 
800f1cd61daSJed Brown /*@
801f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
802f1cd61daSJed Brown 
803f1cd61daSJed Brown    Input Parameter:
804f1cd61daSJed Brown +   ts - the time-step context
805f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
806f1cd61daSJed Brown 
807f1cd61daSJed Brown    Level: beginner
808f1cd61daSJed Brown 
809f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
810f1cd61daSJed Brown @*/
8117087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
812f1cd61daSJed Brown {
8134ac538c5SBarry Smith   PetscErrorCode ierr;
814f1cd61daSJed Brown 
815f1cd61daSJed Brown   PetscFunctionBegin;
8164ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
817f1cd61daSJed Brown   PetscFunctionReturn(0);
818f1cd61daSJed Brown }
819f1cd61daSJed Brown 
8202bfc04deSHong Zhang /*@
8212bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
8222bfc04deSHong Zhang 
8232bfc04deSHong Zhang    Input Parameter:
8242bfc04deSHong Zhang +   ts - the time-step context
8252bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
8262bfc04deSHong Zhang 
8272bfc04deSHong Zhang    Level: beginner
8282bfc04deSHong Zhang 
829f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
8302bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
831f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
8322bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
8332bfc04deSHong Zhang @*/
8347087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool ft)
8352bfc04deSHong Zhang {
8364ac538c5SBarry Smith   PetscErrorCode ierr;
8372bfc04deSHong Zhang 
8382bfc04deSHong Zhang   PetscFunctionBegin;
8394ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
8402bfc04deSHong Zhang   PetscFunctionReturn(0);
8412bfc04deSHong Zhang }
8427cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
8437cb94ee6SHong Zhang /*MC
84496f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
8457cb94ee6SHong Zhang 
8467cb94ee6SHong Zhang    Options Database:
847897c9f78SHong Zhang +    -ts_sundials_type <bdf,adams> -
8487cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
8497cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
8507cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
851897c9f78SHong Zhang .    -ts_sundials_linear_tolerance <tol> -
852f61b2b6aSHong Zhang .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
853897c9f78SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps
85416016685SHong Zhang 
8557cb94ee6SHong Zhang 
85695452b02SPatrick Sanan     Notes:
85795452b02SPatrick Sanan     This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply,
858897c9f78SHong Zhang            only PETSc PC options.
8597cb94ee6SHong Zhang 
8607cb94ee6SHong Zhang     Level: beginner
8617cb94ee6SHong Zhang 
862f61b2b6aSHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
863a43b19c4SJed Brown            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
8647cb94ee6SHong Zhang 
8657cb94ee6SHong Zhang M*/
8668cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts)
8677cb94ee6SHong Zhang {
8687cb94ee6SHong Zhang   TS_Sundials    *cvode;
8697cb94ee6SHong Zhang   PetscErrorCode ierr;
87042528757SHong Zhang   PC             pc;
8717cb94ee6SHong Zhang 
8727cb94ee6SHong Zhang   PetscFunctionBegin;
873277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
87428adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
87528adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
876214bc6a2SJed Brown   ts->ops->setup          = TSSetUp_Sundials;
877b4eba00bSSean Farley   ts->ops->step           = TSStep_Sundials;
878b4eba00bSSean Farley   ts->ops->interpolate    = TSInterpolate_Sundials;
879214bc6a2SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
8802ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
8817cb94ee6SHong Zhang 
882b00a9115SJed Brown   ierr = PetscNewLog(ts,&cvode);CHKERRQ(ierr);
883bbd56ea5SKarl Rupp 
884efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
885efd4aadfSBarry Smith 
8867cb94ee6SHong Zhang   ts->data           = (void*)cvode;
8876fadb2cdSHong Zhang   cvode->cvode_type  = SUNDIALS_BDF;
8886fadb2cdSHong Zhang   cvode->gtype       = SUNDIALS_CLASSICAL_GS;
889f61b2b6aSHong Zhang   cvode->maxl        = 5;
8907cb94ee6SHong Zhang   cvode->linear_tol  = .05;
891b4eba00bSSean Farley   cvode->monitorstep = PETSC_TRUE;
8927cb94ee6SHong Zhang 
893ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));CHKERRQ(ierr);
894f1cd61daSJed Brown 
895f1cd61daSJed Brown   cvode->mindt = -1.;
896f1cd61daSJed Brown   cvode->maxdt = -1.;
897f1cd61daSJed Brown 
8987cb94ee6SHong Zhang   /* set tolerance for Sundials */
8997cb94ee6SHong Zhang   cvode->reltol = 1e-6;
9002c823083SHong Zhang   cvode->abstol = 1e-6;
9017cb94ee6SHong Zhang 
90242528757SHong Zhang   /* set PCNONE as default pctype */
90342528757SHong Zhang   ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr);
90442528757SHong Zhang   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
90542528757SHong Zhang 
906bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);CHKERRQ(ierr);
907bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);CHKERRQ(ierr);
908bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
909bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
910bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
911bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
912bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
913bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);CHKERRQ(ierr);
914bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
915bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
9167cb94ee6SHong Zhang   PetscFunctionReturn(0);
9177cb94ee6SHong Zhang }
918