xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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 */
856a740aaSHong Zhang #include "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 */
147cb94ee6SHong Zhang #undef __FUNCT__
157cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials"
16bbd56ea5SKarl Rupp PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,booleantype jok,booleantype *jcurPtr,
17bbd56ea5SKarl Rupp                                   realtype _gamma,void *P_data,N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
187cb94ee6SHong Zhang {
197cb94ee6SHong Zhang   TS             ts     = (TS) P_data;
207cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
21dcbc6d53SSean Farley   PC             pc;
227cb94ee6SHong Zhang   PetscErrorCode ierr;
230679a0aeSJed Brown   Mat            J,P;
240679a0aeSJed Brown   Vec            yy  = cvode->w1,yydot = cvode->ydot;
250679a0aeSJed Brown   PetscReal      gm  = (PetscReal)_gamma;
267cb94ee6SHong Zhang   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
277cb94ee6SHong Zhang   PetscScalar    *y_data;
287cb94ee6SHong Zhang 
297cb94ee6SHong Zhang   PetscFunctionBegin;
300298fd71SBarry Smith   ierr   = TSGetIJacobian(ts,&J,&P,NULL,NULL);CHKERRQ(ierr);
317cb94ee6SHong Zhang   y_data = (PetscScalar*) N_VGetArrayPointer(y);
327cb94ee6SHong Zhang   ierr   = VecPlaceArray(yy,y_data);CHKERRQ(ierr);
330679a0aeSJed Brown   ierr   = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
340679a0aeSJed Brown   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
350679a0aeSJed Brown   ierr     = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr);
367cb94ee6SHong Zhang   ierr     = VecResetArray(yy);CHKERRQ(ierr);
370679a0aeSJed Brown   ierr     = MatScale(P,gm);CHKERRQ(ierr); /* turn into I-gm*Jrest, J is not used by Sundials  */
387cb94ee6SHong Zhang   *jcurPtr = TRUE;
39dcbc6d53SSean Farley   ierr     = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
400679a0aeSJed Brown   ierr     = PCSetOperators(pc,J,P,str);CHKERRQ(ierr);
417cb94ee6SHong Zhang   PetscFunctionReturn(0);
427cb94ee6SHong Zhang }
437cb94ee6SHong Zhang 
447cb94ee6SHong Zhang /*
457cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
467cb94ee6SHong Zhang */
477cb94ee6SHong Zhang #undef __FUNCT__
487cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials"
494ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
504ac7836bSHong Zhang                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
517cb94ee6SHong Zhang {
527cb94ee6SHong Zhang   TS             ts     = (TS) P_data;
537cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
54dcbc6d53SSean Farley   PC             pc;
557cb94ee6SHong Zhang   Vec            rr = cvode->w1,zz = cvode->w2;
567cb94ee6SHong Zhang   PetscErrorCode ierr;
577cb94ee6SHong Zhang   PetscScalar    *r_data,*z_data;
587cb94ee6SHong Zhang 
597cb94ee6SHong Zhang   PetscFunctionBegin;
607cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
617cb94ee6SHong Zhang   r_data = (PetscScalar*) N_VGetArrayPointer(r);
627cb94ee6SHong Zhang   z_data = (PetscScalar*) N_VGetArrayPointer(z);
637cb94ee6SHong Zhang   ierr   = VecPlaceArray(rr,r_data);CHKERRQ(ierr);
647cb94ee6SHong Zhang   ierr   = VecPlaceArray(zz,z_data);CHKERRQ(ierr);
654ac7836bSHong Zhang 
667cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
67dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
687cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz);CHKERRQ(ierr);
697cb94ee6SHong Zhang   ierr = VecResetArray(rr);CHKERRQ(ierr);
707cb94ee6SHong Zhang   ierr = VecResetArray(zz);CHKERRQ(ierr);
717cb94ee6SHong Zhang   PetscFunctionReturn(0);
727cb94ee6SHong Zhang }
737cb94ee6SHong Zhang 
747cb94ee6SHong Zhang /*
757cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
767cb94ee6SHong Zhang */
777cb94ee6SHong Zhang #undef __FUNCT__
787cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials"
794ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
807cb94ee6SHong Zhang {
817cb94ee6SHong Zhang   TS             ts = (TS) ctx;
825fcef5e4SJed Brown   DM             dm;
83942e3340SBarry Smith   DMTS           tsdm;
845fcef5e4SJed Brown   TSIFunction    ifunction;
85*ce94432eSBarry Smith   MPI_Comm       comm;
867cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
870679a0aeSJed Brown   Vec            yy     = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
887cb94ee6SHong Zhang   PetscScalar    *y_data,*ydot_data;
897cb94ee6SHong Zhang   PetscErrorCode ierr;
907cb94ee6SHong Zhang 
917cb94ee6SHong Zhang   PetscFunctionBegin;
92*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
935ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
947cb94ee6SHong Zhang   y_data    = (PetscScalar*) N_VGetArrayPointer(y);
957cb94ee6SHong Zhang   ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot);
964ac7836bSHong Zhang   ierr      = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
974ac7836bSHong Zhang   ierr      = VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr);
984ac7836bSHong Zhang 
995fcef5e4SJed Brown   /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
1005fcef5e4SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
101942e3340SBarry Smith   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
1020298fd71SBarry Smith   ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr);
1035fcef5e4SJed Brown   if (!ifunction) {
104bc0cc02bSJed Brown     ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr);
105bc0cc02bSJed Brown   } else {                      /* If rhsfunction is also set, this computes both parts and shifts them to the right */
106bc0cc02bSJed Brown     ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
1070679a0aeSJed Brown     ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr);
1080679a0aeSJed Brown     ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
109bc0cc02bSJed Brown   }
1107cb94ee6SHong Zhang   ierr = VecResetArray(yy);CHKERRABORT(comm,ierr);
1117cb94ee6SHong Zhang   ierr = VecResetArray(yyd);CHKERRABORT(comm,ierr);
1124ac7836bSHong Zhang   PetscFunctionReturn(0);
1137cb94ee6SHong Zhang }
1147cb94ee6SHong Zhang 
1157cb94ee6SHong Zhang /*
116b4eba00bSSean Farley        TSStep_Sundials - Calls Sundials to integrate the ODE.
1177cb94ee6SHong Zhang */
1187cb94ee6SHong Zhang #undef __FUNCT__
119b4eba00bSSean Farley #define __FUNCT__ "TSStep_Sundials"
120b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts)
1217cb94ee6SHong Zhang {
1227cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1237cb94ee6SHong Zhang   PetscErrorCode ierr;
124b4eba00bSSean Farley   PetscInt       flag;
1259f94935aSHong Zhang   long int       its,nsteps;
1267cb94ee6SHong Zhang   realtype       t,tout;
1277cb94ee6SHong Zhang   PetscScalar    *y_data;
1287cb94ee6SHong Zhang   void           *mem;
1297cb94ee6SHong Zhang 
1307cb94ee6SHong Zhang   PetscFunctionBegin;
13116016685SHong Zhang   mem  = cvode->mem;
1327cb94ee6SHong Zhang   tout = ts->max_time;
1337cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1347cb94ee6SHong Zhang   N_VSetArrayPointer((realtype*)y_data,cvode->y);
1350298fd71SBarry Smith   ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr);
136186e87acSLisandro Dalcin 
1373f2090d5SJed Brown   ierr = TSPreStep(ts);CHKERRQ(ierr);
138186e87acSLisandro Dalcin 
139b8123daeSJed Brown   /* We would like to call TSPreStep() when starting each step (including rejections) and TSPreStage() before each
140b8123daeSJed Brown    * stage solve, but CVode does not appear to support this. */
141bbd56ea5SKarl Rupp   if (cvode->monitorstep) flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
142bbd56ea5SKarl Rupp   else flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
1439f94935aSHong Zhang 
1449f94935aSHong Zhang   if (flag) { /* display error message */
1459f94935aSHong Zhang     switch (flag) {
1469f94935aSHong Zhang       case CV_ILL_INPUT:
1479f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
1489f94935aSHong Zhang         break;
1499f94935aSHong Zhang       case CV_TOO_CLOSE:
1509f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
1519f94935aSHong Zhang         break;
1523c7fefeeSJed Brown       case CV_TOO_MUCH_WORK: {
1539f94935aSHong Zhang         PetscReal      tcur;
1549f94935aSHong Zhang         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
1559f94935aSHong Zhang         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
1569f94935aSHong Zhang         SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%G, nsteps %D exceeds mxstep %D. Increase '-ts_max_steps <>' or modify TSSetDuration()",tcur,nsteps,ts->max_steps);
1573c7fefeeSJed Brown       } break;
1589f94935aSHong Zhang       case CV_TOO_MUCH_ACC:
1599f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
1609f94935aSHong Zhang         break;
1619f94935aSHong Zhang       case CV_ERR_FAILURE:
1629f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
1639f94935aSHong Zhang         break;
1649f94935aSHong Zhang       case CV_CONV_FAILURE:
1659f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
1669f94935aSHong Zhang         break;
1679f94935aSHong Zhang       case CV_LINIT_FAIL:
1689f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
1699f94935aSHong Zhang         break;
1709f94935aSHong Zhang       case CV_LSETUP_FAIL:
1719f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
1729f94935aSHong Zhang         break;
1739f94935aSHong Zhang       case CV_LSOLVE_FAIL:
1749f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
1759f94935aSHong Zhang         break;
1769f94935aSHong Zhang       case CV_RHSFUNC_FAIL:
1779f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
1789f94935aSHong Zhang         break;
1799f94935aSHong Zhang       case CV_FIRST_RHSFUNC_ERR:
1809f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
1819f94935aSHong Zhang         break;
1829f94935aSHong Zhang       case CV_REPTD_RHSFUNC_ERR:
1839f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
1849f94935aSHong Zhang         break;
1859f94935aSHong Zhang       case CV_UNREC_RHSFUNC_ERR:
1869f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
1879f94935aSHong Zhang         break;
1889f94935aSHong Zhang       case CV_RTFUNC_FAIL:
1899f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
1909f94935aSHong Zhang         break;
1919f94935aSHong Zhang       default:
1929f94935aSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
1939f94935aSHong Zhang     }
1949f94935aSHong Zhang   }
1959f94935aSHong Zhang 
1967cb94ee6SHong Zhang   /* copy the solution from cvode->y to cvode->update and sol */
1977cb94ee6SHong Zhang   ierr = VecPlaceArray(cvode->w1,y_data);CHKERRQ(ierr);
1987cb94ee6SHong Zhang   ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
1997cb94ee6SHong Zhang   ierr = VecResetArray(cvode->w1);CHKERRQ(ierr);
200bc0cc02bSJed Brown   ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr);
2017cb94ee6SHong Zhang   ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
2024ac7836bSHong Zhang   ierr = CVSpilsGetNumLinIters(mem, &its);
2035ef26d82SJed Brown   ts->snes_its = its; ts->ksp_its = its;
204186e87acSLisandro Dalcin 
205186e87acSLisandro Dalcin   ts->time_step = t - ts->ptime;
206186e87acSLisandro Dalcin   ts->ptime     = t;
2077cb94ee6SHong Zhang   ts->steps++;
208186e87acSLisandro Dalcin 
2099f94935aSHong Zhang   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
210b4eba00bSSean Farley   if (!cvode->monitorstep) ts->steps = nsteps;
211b4eba00bSSean Farley   PetscFunctionReturn(0);
212b4eba00bSSean Farley }
213b4eba00bSSean Farley 
214b4eba00bSSean Farley #undef __FUNCT__
215b4eba00bSSean Farley #define __FUNCT__ "TSInterpolate_Sundials"
216b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
217b4eba00bSSean Farley {
218b4eba00bSSean Farley   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
219b4eba00bSSean Farley   N_Vector       y;
220b4eba00bSSean Farley   PetscErrorCode ierr;
221b4eba00bSSean Farley   PetscScalar    *x_data;
222b4eba00bSSean Farley   PetscInt       glosize,locsize;
223b4eba00bSSean Farley 
224b4eba00bSSean Farley   PetscFunctionBegin;
225b4eba00bSSean Farley   /* get the vector size */
226b4eba00bSSean Farley   ierr = VecGetSize(X,&glosize);CHKERRQ(ierr);
227b4eba00bSSean Farley   ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr);
228b4eba00bSSean Farley 
229b4eba00bSSean Farley   /* allocate the memory for N_Vec y */
230b4eba00bSSean Farley   y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
231b4eba00bSSean Farley   if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");
232b4eba00bSSean Farley 
233b4eba00bSSean Farley   ierr = VecGetArray(X,&x_data);CHKERRQ(ierr);
234b4eba00bSSean Farley   N_VSetArrayPointer((realtype*)x_data,y);
235b4eba00bSSean Farley   ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr);
236b4eba00bSSean Farley   ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr);
2377cb94ee6SHong Zhang   PetscFunctionReturn(0);
2387cb94ee6SHong Zhang }
2397cb94ee6SHong Zhang 
2407cb94ee6SHong Zhang #undef __FUNCT__
241277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials"
242277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts)
243277b19d0SLisandro Dalcin {
244277b19d0SLisandro Dalcin   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
245277b19d0SLisandro Dalcin   PetscErrorCode ierr;
246277b19d0SLisandro Dalcin 
247277b19d0SLisandro Dalcin   PetscFunctionBegin;
2485c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
2490679a0aeSJed Brown   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
2505c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
2515c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
252bbd56ea5SKarl Rupp   if (cvode->mem) CVodeFree(&cvode->mem);
253277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
254277b19d0SLisandro Dalcin }
255277b19d0SLisandro Dalcin 
256277b19d0SLisandro Dalcin #undef __FUNCT__
2577cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials"
2587cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
2597cb94ee6SHong Zhang {
2607cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2617cb94ee6SHong Zhang   PetscErrorCode ierr;
2627cb94ee6SHong Zhang 
2637cb94ee6SHong Zhang   PetscFunctionBegin;
264277b19d0SLisandro Dalcin   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
265a07356b8SSatish Balay   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
266277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
2670298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",NULL);CHKERRQ(ierr);
2680298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C","",NULL);CHKERRQ(ierr);
2690298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",NULL);CHKERRQ(ierr);
2700298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",NULL);CHKERRQ(ierr);
2710298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",NULL);CHKERRQ(ierr);
2720298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",NULL);CHKERRQ(ierr);
2730298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",NULL);CHKERRQ(ierr);
2740298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",NULL);CHKERRQ(ierr);
2750298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",NULL);CHKERRQ(ierr);
2760298fd71SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",NULL);CHKERRQ(ierr);
2777cb94ee6SHong Zhang   PetscFunctionReturn(0);
2787cb94ee6SHong Zhang }
2797cb94ee6SHong Zhang 
2807cb94ee6SHong Zhang #undef __FUNCT__
281214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials"
282214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts)
2837cb94ee6SHong Zhang {
2847cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2857cb94ee6SHong Zhang   PetscErrorCode ierr;
28616016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2877cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
28816016685SHong Zhang   void           *mem;
289dcbc6d53SSean Farley   PC             pc;
29019fd82e9SBarry Smith   PCType         pctype;
291ace3abfcSBarry Smith   PetscBool      pcnone;
2927cb94ee6SHong Zhang 
2937cb94ee6SHong Zhang   PetscFunctionBegin;
2947cb94ee6SHong Zhang   /* get the vector size */
2957cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2967cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
2977cb94ee6SHong Zhang 
2987cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
299a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
300e32f2f54SBarry Smith   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
3017cb94ee6SHong Zhang 
30228adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
3037cb94ee6SHong Zhang   ierr   = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
3047cb94ee6SHong Zhang   y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y);
3057cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
3060298fd71SBarry Smith   ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr);
30742528757SHong Zhang 
3087cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
3090679a0aeSJed Brown   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
3107cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
3110679a0aeSJed Brown   ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr);
3127cb94ee6SHong Zhang 
3137cb94ee6SHong Zhang   /*
3147cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
3157cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
3167cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
3177cb94ee6SHong Zhang   */
318*ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
319*ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
3207cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
3217cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
32216016685SHong Zhang 
32316016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
32416016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
325e32f2f54SBarry Smith   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
32616016685SHong Zhang   cvode->mem = mem;
32716016685SHong Zhang 
32816016685SHong Zhang   /* Set the pointer to user-defined data */
32916016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
330e32f2f54SBarry Smith   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
33116016685SHong Zhang 
332fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
333cdbf8f93SLisandro Dalcin   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
334*ce94432eSBarry Smith   if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed");
335f1cd61daSJed Brown   if (cvode->mindt > 0) {
336f1cd61daSJed Brown     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
3379f94935aSHong Zhang     if (flag) {
338*ce94432eSBarry Smith       if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
339*ce94432eSBarry 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");
340*ce94432eSBarry Smith       else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed");
3419f94935aSHong Zhang     }
342f1cd61daSJed Brown   }
343f1cd61daSJed Brown   if (cvode->maxdt > 0) {
344f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
345*ce94432eSBarry Smith     if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
346f1cd61daSJed Brown   }
347f1cd61daSJed Brown 
34816016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
34916016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
35016016685SHong Zhang    * the initial dependent variable vector cvode->y */
35116016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
352bbd56ea5SKarl Rupp   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
35316016685SHong Zhang 
3549f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
35516016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
356bbd56ea5SKarl Rupp   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
35716016685SHong Zhang 
3589f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3599f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
36016016685SHong Zhang 
36116016685SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
36216016685SHong Zhang   /* setup the ode integrator with the given preconditioner */
363dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
364dcbc6d53SSean Farley   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
365251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
36616016685SHong Zhang   if (pcnone) {
36716016685SHong Zhang     flag = CVSpgmr(mem,PREC_NONE,0);
368e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
36916016685SHong Zhang   } else {
370f61b2b6aSHong Zhang     flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
371e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
37216016685SHong Zhang 
37316016685SHong Zhang     /* Set preconditioner and solve routines Precond and PSolve,
37416016685SHong Zhang      and the pointer to the user-defined block data */
37516016685SHong Zhang     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
376e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
37716016685SHong Zhang   }
37816016685SHong Zhang 
37916016685SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
380e32f2f54SBarry Smith   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
3817cb94ee6SHong Zhang   PetscFunctionReturn(0);
3827cb94ee6SHong Zhang }
3837cb94ee6SHong Zhang 
3846fadb2cdSHong Zhang /* type of CVODE linear multistep method */
3856a6fc655SJed Brown const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0};
3866fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
3876a6fc655SJed Brown const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0};
388a04cf4d8SBarry Smith 
3897cb94ee6SHong Zhang #undef __FUNCT__
390214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials"
391214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts)
3927cb94ee6SHong Zhang {
3937cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3947cb94ee6SHong Zhang   PetscErrorCode ierr;
3957cb94ee6SHong Zhang   int            indx;
396ace3abfcSBarry Smith   PetscBool      flag;
3977bda3a07SJed Brown   PC             pc;
3987cb94ee6SHong Zhang 
3997cb94ee6SHong Zhang   PetscFunctionBegin;
4007cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
4016fadb2cdSHong Zhang   ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
4027cb94ee6SHong Zhang   if (flag) {
4036fadb2cdSHong Zhang     ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
4047cb94ee6SHong Zhang   }
4056fadb2cdSHong Zhang   ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
4067cb94ee6SHong Zhang   if (flag) {
4077cb94ee6SHong Zhang     ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
4087cb94ee6SHong Zhang   }
4090298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);CHKERRQ(ierr);
4100298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);CHKERRQ(ierr);
4110298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);CHKERRQ(ierr);
4120298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);CHKERRQ(ierr);
4137cb94ee6SHong Zhang   ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
414f61b2b6aSHong Zhang   ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);CHKERRQ(ierr);
4150298fd71SBarry Smith   ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);CHKERRQ(ierr);
4167cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4177bda3a07SJed Brown   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
4187bda3a07SJed Brown   ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
4197cb94ee6SHong Zhang   PetscFunctionReturn(0);
4207cb94ee6SHong Zhang }
4217cb94ee6SHong Zhang 
4227cb94ee6SHong Zhang #undef __FUNCT__
4237cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
4247cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
4257cb94ee6SHong Zhang {
4267cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4277cb94ee6SHong Zhang   PetscErrorCode ierr;
4287cb94ee6SHong Zhang   char           *type;
4297cb94ee6SHong Zhang   char           atype[] = "Adams";
4307cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
431ace3abfcSBarry Smith   PetscBool      iascii,isstring;
4322c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
4332c823083SHong Zhang   PetscInt       qlast,qcur;
4345d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
43542528757SHong Zhang   PC             pc;
4367cb94ee6SHong Zhang 
4377cb94ee6SHong Zhang   PetscFunctionBegin;
438bbd56ea5SKarl Rupp   if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
439bbd56ea5SKarl Rupp   else                                     type = btype;
4407cb94ee6SHong Zhang 
441251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
442251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
4437cb94ee6SHong Zhang   if (iascii) {
4447cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
4457cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
4467cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
4477cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
448f61b2b6aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr);
4497cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4507cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4517cb94ee6SHong Zhang     } else {
4527cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4537cb94ee6SHong Zhang     }
454f1cd61daSJed Brown     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
455f1cd61daSJed Brown     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
4562c823083SHong Zhang 
4575d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4585d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
4595d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
4602c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
4612c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
4622c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
4632c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
4642c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
4652c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
4662c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
4672c823083SHong Zhang 
4682c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
4692c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
4702c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
4712c823083SHong Zhang 
4722c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
4732c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
4742c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
4752c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
47642528757SHong Zhang 
47742528757SHong Zhang     ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
47842528757SHong Zhang     ierr = PCView(pc,viewer);CHKERRQ(ierr);
4792c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4802c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
4812c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
4822c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
48342528757SHong Zhang 
4842c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4852c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
4862c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4872c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
4887cb94ee6SHong Zhang   } else if (isstring) {
4897cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
4907cb94ee6SHong Zhang   }
4917cb94ee6SHong Zhang   PetscFunctionReturn(0);
4927cb94ee6SHong Zhang }
4937cb94ee6SHong Zhang 
4947cb94ee6SHong Zhang 
4957cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
4967cb94ee6SHong Zhang EXTERN_C_BEGIN
4977cb94ee6SHong Zhang #undef __FUNCT__
4987cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
4997087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
5007cb94ee6SHong Zhang {
5017cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5027cb94ee6SHong Zhang 
5037cb94ee6SHong Zhang   PetscFunctionBegin;
5047cb94ee6SHong Zhang   cvode->cvode_type = type;
5057cb94ee6SHong Zhang   PetscFunctionReturn(0);
5067cb94ee6SHong Zhang }
5077cb94ee6SHong Zhang EXTERN_C_END
5087cb94ee6SHong Zhang 
5097cb94ee6SHong Zhang EXTERN_C_BEGIN
5107cb94ee6SHong Zhang #undef __FUNCT__
511f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl_Sundials"
512f61b2b6aSHong Zhang PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
5137cb94ee6SHong Zhang {
5147cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5157cb94ee6SHong Zhang 
5167cb94ee6SHong Zhang   PetscFunctionBegin;
517f61b2b6aSHong Zhang   cvode->maxl = maxl;
5187cb94ee6SHong Zhang   PetscFunctionReturn(0);
5197cb94ee6SHong Zhang }
5207cb94ee6SHong Zhang EXTERN_C_END
5217cb94ee6SHong Zhang 
5227cb94ee6SHong Zhang EXTERN_C_BEGIN
5237cb94ee6SHong Zhang #undef __FUNCT__
5247cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
5257087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
5267cb94ee6SHong Zhang {
5277cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5287cb94ee6SHong Zhang 
5297cb94ee6SHong Zhang   PetscFunctionBegin;
5307cb94ee6SHong Zhang   cvode->linear_tol = tol;
5317cb94ee6SHong Zhang   PetscFunctionReturn(0);
5327cb94ee6SHong Zhang }
5337cb94ee6SHong Zhang EXTERN_C_END
5347cb94ee6SHong Zhang 
5357cb94ee6SHong Zhang EXTERN_C_BEGIN
5367cb94ee6SHong Zhang #undef __FUNCT__
5377cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
5387087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
5397cb94ee6SHong Zhang {
5407cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5417cb94ee6SHong Zhang 
5427cb94ee6SHong Zhang   PetscFunctionBegin;
5437cb94ee6SHong Zhang   cvode->gtype = type;
5447cb94ee6SHong Zhang   PetscFunctionReturn(0);
5457cb94ee6SHong Zhang }
5467cb94ee6SHong Zhang EXTERN_C_END
5477cb94ee6SHong Zhang 
5487cb94ee6SHong Zhang EXTERN_C_BEGIN
5497cb94ee6SHong Zhang #undef __FUNCT__
5507cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
5517087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
5527cb94ee6SHong Zhang {
5537cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5547cb94ee6SHong Zhang 
5557cb94ee6SHong Zhang   PetscFunctionBegin;
5567cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5577cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
5587cb94ee6SHong Zhang   PetscFunctionReturn(0);
5597cb94ee6SHong Zhang }
5607cb94ee6SHong Zhang EXTERN_C_END
5617cb94ee6SHong Zhang 
5627cb94ee6SHong Zhang EXTERN_C_BEGIN
5637cb94ee6SHong Zhang #undef __FUNCT__
564f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
5657087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
566f1cd61daSJed Brown {
567f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
568f1cd61daSJed Brown 
569f1cd61daSJed Brown   PetscFunctionBegin;
570f1cd61daSJed Brown   cvode->mindt = mindt;
571f1cd61daSJed Brown   PetscFunctionReturn(0);
572f1cd61daSJed Brown }
573f1cd61daSJed Brown EXTERN_C_END
574f1cd61daSJed Brown 
575f1cd61daSJed Brown EXTERN_C_BEGIN
576f1cd61daSJed Brown #undef __FUNCT__
577f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
5787087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
579f1cd61daSJed Brown {
580f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
581f1cd61daSJed Brown 
582f1cd61daSJed Brown   PetscFunctionBegin;
583f1cd61daSJed Brown   cvode->maxdt = maxdt;
584f1cd61daSJed Brown   PetscFunctionReturn(0);
585f1cd61daSJed Brown }
586f1cd61daSJed Brown EXTERN_C_END
587f1cd61daSJed Brown 
588f1cd61daSJed Brown EXTERN_C_BEGIN
589f1cd61daSJed Brown #undef __FUNCT__
5907cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
5917087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
5927cb94ee6SHong Zhang {
593dcbc6d53SSean Farley   SNES           snes;
594dcbc6d53SSean Farley   KSP            ksp;
595dcbc6d53SSean Farley   PetscErrorCode ierr;
5967cb94ee6SHong Zhang 
5977cb94ee6SHong Zhang   PetscFunctionBegin;
598dcbc6d53SSean Farley   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
599dcbc6d53SSean Farley   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
600dcbc6d53SSean Farley   ierr = KSPGetPC(ksp,pc);CHKERRQ(ierr);
6017cb94ee6SHong Zhang   PetscFunctionReturn(0);
6027cb94ee6SHong Zhang }
6037cb94ee6SHong Zhang EXTERN_C_END
6047cb94ee6SHong Zhang 
6057cb94ee6SHong Zhang EXTERN_C_BEGIN
6067cb94ee6SHong Zhang #undef __FUNCT__
6077cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
6087087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
6097cb94ee6SHong Zhang {
6107cb94ee6SHong Zhang   PetscFunctionBegin;
6115ef26d82SJed Brown   if (nonlin) *nonlin = ts->snes_its;
6125ef26d82SJed Brown   if (lin)    *lin    = ts->ksp_its;
6137cb94ee6SHong Zhang   PetscFunctionReturn(0);
6147cb94ee6SHong Zhang }
6157cb94ee6SHong Zhang EXTERN_C_END
6167cb94ee6SHong Zhang 
6177cb94ee6SHong Zhang EXTERN_C_BEGIN
6187cb94ee6SHong Zhang #undef __FUNCT__
6192bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
6207087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)
6212bfc04deSHong Zhang {
6222bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6232bfc04deSHong Zhang 
6242bfc04deSHong Zhang   PetscFunctionBegin;
6252bfc04deSHong Zhang   cvode->monitorstep = s;
6262bfc04deSHong Zhang   PetscFunctionReturn(0);
6272bfc04deSHong Zhang }
6282bfc04deSHong Zhang EXTERN_C_END
6297cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
6307cb94ee6SHong Zhang 
6317cb94ee6SHong Zhang #undef __FUNCT__
6327cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
6337cb94ee6SHong Zhang /*@C
6347cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
6357cb94ee6SHong Zhang 
6367cb94ee6SHong Zhang    Not Collective
6377cb94ee6SHong Zhang 
6387cb94ee6SHong Zhang    Input parameters:
6397cb94ee6SHong Zhang .    ts     - the time-step context
6407cb94ee6SHong Zhang 
6417cb94ee6SHong Zhang    Output Parameters:
6427cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
6437cb94ee6SHong Zhang -   lin    - number of linear iterations
6447cb94ee6SHong Zhang 
6457cb94ee6SHong Zhang    Level: advanced
6467cb94ee6SHong Zhang 
6477cb94ee6SHong Zhang    Notes:
6487cb94ee6SHong Zhang     These return the number since the creation of the TS object
6497cb94ee6SHong Zhang 
6507cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
6517cb94ee6SHong Zhang 
652f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
6537cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
654f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
655a43b19c4SJed Brown           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
6567cb94ee6SHong Zhang 
6577cb94ee6SHong Zhang @*/
6587087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
6597cb94ee6SHong Zhang {
6604ac538c5SBarry Smith   PetscErrorCode ierr;
6617cb94ee6SHong Zhang 
6627cb94ee6SHong Zhang   PetscFunctionBegin;
6634ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
6647cb94ee6SHong Zhang   PetscFunctionReturn(0);
6657cb94ee6SHong Zhang }
6667cb94ee6SHong Zhang 
6677cb94ee6SHong Zhang #undef __FUNCT__
6687cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
6697cb94ee6SHong Zhang /*@
6707cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
6717cb94ee6SHong Zhang 
6723f9fe445SBarry Smith    Logically Collective on TS
6737cb94ee6SHong Zhang 
6747cb94ee6SHong Zhang    Input parameters:
6757cb94ee6SHong Zhang +    ts     - the time-step context
6767cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
6777cb94ee6SHong Zhang 
6787cb94ee6SHong Zhang    Level: intermediate
6797cb94ee6SHong Zhang 
6807cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
6817cb94ee6SHong Zhang 
682f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetMaxl(),
6837cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
684f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
6857cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
686a43b19c4SJed Brown           TSSetExactFinalTime()
6877cb94ee6SHong Zhang @*/
6887087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
6897cb94ee6SHong Zhang {
6904ac538c5SBarry Smith   PetscErrorCode ierr;
6917cb94ee6SHong Zhang 
6927cb94ee6SHong Zhang   PetscFunctionBegin;
6934ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
6947cb94ee6SHong Zhang   PetscFunctionReturn(0);
6957cb94ee6SHong Zhang }
6967cb94ee6SHong Zhang 
6977cb94ee6SHong Zhang #undef __FUNCT__
698f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl"
6997cb94ee6SHong Zhang /*@
700f61b2b6aSHong Zhang    TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
7017cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
702f61b2b6aSHong Zhang        this is the maximum number of GMRES steps that will be used.
7037cb94ee6SHong Zhang 
7043f9fe445SBarry Smith    Logically Collective on TS
7057cb94ee6SHong Zhang 
7067cb94ee6SHong Zhang    Input parameters:
7077cb94ee6SHong Zhang +    ts      - the time-step context
708f61b2b6aSHong Zhang -    maxl - number of direction vectors (the dimension of Krylov subspace).
7097cb94ee6SHong Zhang 
7107cb94ee6SHong Zhang    Level: advanced
7117cb94ee6SHong Zhang 
712f61b2b6aSHong Zhang .keywords: GMRES
7137cb94ee6SHong Zhang 
7147cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
7157cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
716f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7177cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
718a43b19c4SJed Brown           TSSetExactFinalTime()
7197cb94ee6SHong Zhang 
7207cb94ee6SHong Zhang @*/
721f61b2b6aSHong Zhang PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
7227cb94ee6SHong Zhang {
7234ac538c5SBarry Smith   PetscErrorCode ierr;
7247cb94ee6SHong Zhang 
7257cb94ee6SHong Zhang   PetscFunctionBegin;
726f61b2b6aSHong Zhang   PetscValidLogicalCollectiveInt(ts,maxl,2);
727f61b2b6aSHong Zhang   ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr);
7287cb94ee6SHong Zhang   PetscFunctionReturn(0);
7297cb94ee6SHong Zhang }
7307cb94ee6SHong Zhang 
7317cb94ee6SHong Zhang #undef __FUNCT__
7327cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
7337cb94ee6SHong Zhang /*@
7347cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
7357cb94ee6SHong Zhang        system by SUNDIALS.
7367cb94ee6SHong Zhang 
7373f9fe445SBarry Smith    Logically Collective on TS
7387cb94ee6SHong Zhang 
7397cb94ee6SHong Zhang    Input parameters:
7407cb94ee6SHong Zhang +    ts     - the time-step context
7417cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
7427cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
7437cb94ee6SHong Zhang 
7447cb94ee6SHong Zhang    Level: advanced
7457cb94ee6SHong Zhang 
7467cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
7477cb94ee6SHong Zhang 
748f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
7497cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
750f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7517cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
752a43b19c4SJed Brown           TSSetExactFinalTime()
7537cb94ee6SHong Zhang 
7547cb94ee6SHong Zhang @*/
7557087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
7567cb94ee6SHong Zhang {
7574ac538c5SBarry Smith   PetscErrorCode ierr;
7587cb94ee6SHong Zhang 
7597cb94ee6SHong Zhang   PetscFunctionBegin;
760c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,tol,2);
7614ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
7627cb94ee6SHong Zhang   PetscFunctionReturn(0);
7637cb94ee6SHong Zhang }
7647cb94ee6SHong Zhang 
7657cb94ee6SHong Zhang #undef __FUNCT__
7667cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
7677cb94ee6SHong Zhang /*@
7687cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
7697cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
7707cb94ee6SHong Zhang 
7713f9fe445SBarry Smith    Logically Collective on TS
7727cb94ee6SHong Zhang 
7737cb94ee6SHong Zhang    Input parameters:
7747cb94ee6SHong Zhang +    ts  - the time-step context
7757cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
7767cb94ee6SHong Zhang 
7777cb94ee6SHong Zhang    Level: advanced
7787cb94ee6SHong Zhang 
7797cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
7807cb94ee6SHong Zhang 
781f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
7827cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
783f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7847cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
785a43b19c4SJed Brown           TSSetExactFinalTime()
7867cb94ee6SHong Zhang 
7877cb94ee6SHong Zhang @*/
7887087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
7897cb94ee6SHong Zhang {
7904ac538c5SBarry Smith   PetscErrorCode ierr;
7917cb94ee6SHong Zhang 
7927cb94ee6SHong Zhang   PetscFunctionBegin;
7934ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
7947cb94ee6SHong Zhang   PetscFunctionReturn(0);
7957cb94ee6SHong Zhang }
7967cb94ee6SHong Zhang 
7977cb94ee6SHong Zhang #undef __FUNCT__
7987cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
7997cb94ee6SHong Zhang /*@
8007cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
8017cb94ee6SHong Zhang                          Sundials for error control.
8027cb94ee6SHong Zhang 
8033f9fe445SBarry Smith    Logically Collective on TS
8047cb94ee6SHong Zhang 
8057cb94ee6SHong Zhang    Input parameters:
8067cb94ee6SHong Zhang +    ts  - the time-step context
8077cb94ee6SHong Zhang .    aabs - the absolute tolerance
8087cb94ee6SHong Zhang -    rel - the relative tolerance
8097cb94ee6SHong Zhang 
8107cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
8117cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
8127cb94ee6SHong Zhang 
8137cb94ee6SHong Zhang    Level: intermediate
8147cb94ee6SHong Zhang 
8157cb94ee6SHong Zhang .keywords: Sundials, tolerance
8167cb94ee6SHong Zhang 
817f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
8187cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
819f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
8207cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
821a43b19c4SJed Brown           TSSetExactFinalTime()
8227cb94ee6SHong Zhang 
8237cb94ee6SHong Zhang @*/
8247087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
8257cb94ee6SHong Zhang {
8264ac538c5SBarry Smith   PetscErrorCode ierr;
8277cb94ee6SHong Zhang 
8287cb94ee6SHong Zhang   PetscFunctionBegin;
8294ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
8307cb94ee6SHong Zhang   PetscFunctionReturn(0);
8317cb94ee6SHong Zhang }
8327cb94ee6SHong Zhang 
8337cb94ee6SHong Zhang #undef __FUNCT__
8347cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
8357cb94ee6SHong Zhang /*@
8367cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
8377cb94ee6SHong Zhang 
8387cb94ee6SHong Zhang    Input Parameter:
8397cb94ee6SHong Zhang .    ts - the time-step context
8407cb94ee6SHong Zhang 
8417cb94ee6SHong Zhang    Output Parameter:
8427cb94ee6SHong Zhang .    pc - the preconditioner context
8437cb94ee6SHong Zhang 
8447cb94ee6SHong Zhang    Level: advanced
8457cb94ee6SHong Zhang 
846f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
8477cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
848f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
8497cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
8507cb94ee6SHong Zhang @*/
8517087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
8527cb94ee6SHong Zhang {
8534ac538c5SBarry Smith   PetscErrorCode ierr;
8547cb94ee6SHong Zhang 
8557cb94ee6SHong Zhang   PetscFunctionBegin;
8564ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));CHKERRQ(ierr);
8577cb94ee6SHong Zhang   PetscFunctionReturn(0);
8587cb94ee6SHong Zhang }
8597cb94ee6SHong Zhang 
8607cb94ee6SHong Zhang #undef __FUNCT__
861f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep"
862f1cd61daSJed Brown /*@
863f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
864f1cd61daSJed Brown 
865f1cd61daSJed Brown    Input Parameter:
866f1cd61daSJed Brown +   ts - the time-step context
867f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
868f1cd61daSJed Brown 
869fc6b9e64SJed Brown    Note:
870fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
871fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
872fc6b9e64SJed Brown 
873f1cd61daSJed Brown    Level: beginner
874f1cd61daSJed Brown 
875f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
876f1cd61daSJed Brown @*/
8777087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
878f1cd61daSJed Brown {
8794ac538c5SBarry Smith   PetscErrorCode ierr;
880f1cd61daSJed Brown 
881f1cd61daSJed Brown   PetscFunctionBegin;
8824ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
883f1cd61daSJed Brown   PetscFunctionReturn(0);
884f1cd61daSJed Brown }
885f1cd61daSJed Brown 
886f1cd61daSJed Brown #undef __FUNCT__
887f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep"
888f1cd61daSJed Brown /*@
889f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
890f1cd61daSJed Brown 
891f1cd61daSJed Brown    Input Parameter:
892f1cd61daSJed Brown +   ts - the time-step context
893f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
894f1cd61daSJed Brown 
895f1cd61daSJed Brown    Level: beginner
896f1cd61daSJed Brown 
897f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
898f1cd61daSJed Brown @*/
8997087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
900f1cd61daSJed Brown {
9014ac538c5SBarry Smith   PetscErrorCode ierr;
902f1cd61daSJed Brown 
903f1cd61daSJed Brown   PetscFunctionBegin;
9044ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
905f1cd61daSJed Brown   PetscFunctionReturn(0);
906f1cd61daSJed Brown }
907f1cd61daSJed Brown 
908f1cd61daSJed Brown #undef __FUNCT__
9092bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps"
9102bfc04deSHong Zhang /*@
9112bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
9122bfc04deSHong Zhang 
9132bfc04deSHong Zhang    Input Parameter:
9142bfc04deSHong Zhang +   ts - the time-step context
9152bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
9162bfc04deSHong Zhang 
9172bfc04deSHong Zhang    Level: beginner
9182bfc04deSHong Zhang 
919f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
9202bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
921f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
9222bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
9232bfc04deSHong Zhang @*/
9247087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool ft)
9252bfc04deSHong Zhang {
9264ac538c5SBarry Smith   PetscErrorCode ierr;
9272bfc04deSHong Zhang 
9282bfc04deSHong Zhang   PetscFunctionBegin;
9294ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
9302bfc04deSHong Zhang   PetscFunctionReturn(0);
9312bfc04deSHong Zhang }
9327cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
9337cb94ee6SHong Zhang /*MC
93496f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
9357cb94ee6SHong Zhang 
9367cb94ee6SHong Zhang    Options Database:
9377cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
9387cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
9397cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
9407cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
9417cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
942f61b2b6aSHong Zhang .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
94316016685SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
94416016685SHong Zhang 
9457cb94ee6SHong Zhang 
9467cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
9477cb94ee6SHong Zhang            only PETSc PC options
9487cb94ee6SHong Zhang 
9497cb94ee6SHong Zhang     Level: beginner
9507cb94ee6SHong Zhang 
951f61b2b6aSHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
952a43b19c4SJed Brown            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
9537cb94ee6SHong Zhang 
9547cb94ee6SHong Zhang M*/
9557cb94ee6SHong Zhang EXTERN_C_BEGIN
9567cb94ee6SHong Zhang #undef __FUNCT__
9577cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
9587087cfbeSBarry Smith PetscErrorCode  TSCreate_Sundials(TS ts)
9597cb94ee6SHong Zhang {
9607cb94ee6SHong Zhang   TS_Sundials    *cvode;
9617cb94ee6SHong Zhang   PetscErrorCode ierr;
96242528757SHong Zhang   PC             pc;
9637cb94ee6SHong Zhang 
9647cb94ee6SHong Zhang   PetscFunctionBegin;
965277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
96628adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
96728adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
968214bc6a2SJed Brown   ts->ops->setup          = TSSetUp_Sundials;
969b4eba00bSSean Farley   ts->ops->step           = TSStep_Sundials;
970b4eba00bSSean Farley   ts->ops->interpolate    = TSInterpolate_Sundials;
971214bc6a2SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
9727cb94ee6SHong Zhang 
97338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
974bbd56ea5SKarl Rupp 
9757cb94ee6SHong Zhang   ts->data           = (void*)cvode;
9766fadb2cdSHong Zhang   cvode->cvode_type  = SUNDIALS_BDF;
9776fadb2cdSHong Zhang   cvode->gtype       = SUNDIALS_CLASSICAL_GS;
978f61b2b6aSHong Zhang   cvode->maxl        = 5;
9797cb94ee6SHong Zhang   cvode->linear_tol  = .05;
980b4eba00bSSean Farley   cvode->monitorstep = PETSC_TRUE;
9817cb94ee6SHong Zhang 
982*ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));CHKERRQ(ierr);
983f1cd61daSJed Brown 
984f1cd61daSJed Brown   cvode->mindt = -1.;
985f1cd61daSJed Brown   cvode->maxdt = -1.;
986f1cd61daSJed Brown 
9877cb94ee6SHong Zhang   /* set tolerance for Sundials */
9887cb94ee6SHong Zhang   cvode->reltol = 1e-6;
9892c823083SHong Zhang   cvode->abstol = 1e-6;
9907cb94ee6SHong Zhang 
99142528757SHong Zhang   /* set PCNONE as default pctype */
99242528757SHong Zhang   ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr);
99342528757SHong Zhang   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
99442528757SHong Zhang 
99549354f04SShri Abhyankar   if (ts->exact_final_time != TS_EXACTFINALTIME_STEPOVER) ts->exact_final_time = TS_EXACTFINALTIME_STEPOVER;
996a43b19c4SJed Brown 
9977cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
9987cb94ee6SHong Zhang                                            TSSundialsSetType_Sundials);CHKERRQ(ierr);
999f61b2b6aSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C",
1000f61b2b6aSHong Zhang                                            "TSSundialsSetMaxl_Sundials",
1001f61b2b6aSHong Zhang                                            TSSundialsSetMaxl_Sundials);CHKERRQ(ierr);
10027cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
10037cb94ee6SHong Zhang                                            "TSSundialsSetLinearTolerance_Sundials",
10047cb94ee6SHong Zhang                                            TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
10057cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
10067cb94ee6SHong Zhang                                            "TSSundialsSetGramSchmidtType_Sundials",
10077cb94ee6SHong Zhang                                            TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
10087cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
10097cb94ee6SHong Zhang                                            "TSSundialsSetTolerance_Sundials",
10107cb94ee6SHong Zhang                                            TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1011f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1012f1cd61daSJed Brown                                            "TSSundialsSetMinTimeStep_Sundials",
1013f1cd61daSJed Brown                                            TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1014f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1015f1cd61daSJed Brown                                            "TSSundialsSetMaxTimeStep_Sundials",
1016f1cd61daSJed Brown                                            TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
10177cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
10187cb94ee6SHong Zhang                                            "TSSundialsGetPC_Sundials",
10197cb94ee6SHong Zhang                                            TSSundialsGetPC_Sundials);CHKERRQ(ierr);
10207cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
10217cb94ee6SHong Zhang                                            "TSSundialsGetIterations_Sundials",
10227cb94ee6SHong Zhang                                            TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
10232bfc04deSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
10242bfc04deSHong Zhang                                            "TSSundialsMonitorInternalSteps_Sundials",
10252bfc04deSHong Zhang                                            TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
10267cb94ee6SHong Zhang   PetscFunctionReturn(0);
10277cb94ee6SHong Zhang }
10287cb94ee6SHong Zhang EXTERN_C_END
1029