xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 5c6b4a3d39520986e9d877fb6aedc988c8b7badf)
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"
167cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
177cb94ee6SHong Zhang                     booleantype jok,booleantype *jcurPtr,
187cb94ee6SHong Zhang                     realtype _gamma,void *P_data,
197cb94ee6SHong Zhang                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
207cb94ee6SHong Zhang {
217cb94ee6SHong Zhang   TS             ts = (TS) P_data;
227cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
23dcbc6d53SSean Farley   PC             pc;
247cb94ee6SHong Zhang   PetscErrorCode ierr;
25089b2837SJed Brown   Mat            Jac;
267cb94ee6SHong Zhang   Vec            yy = cvode->w1;
277cb94ee6SHong Zhang   PetscScalar    one = 1.0,gm;
287cb94ee6SHong Zhang   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
297cb94ee6SHong Zhang   PetscScalar    *y_data;
307cb94ee6SHong Zhang 
317cb94ee6SHong Zhang   PetscFunctionBegin;
32089b2837SJed Brown   ierr = TSGetRHSJacobian(ts,PETSC_NULL,&Jac,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
337cb94ee6SHong Zhang   /* This allows us to construct preconditioners in-place if we like */
347cb94ee6SHong Zhang   ierr = MatSetUnfactored(Jac);CHKERRQ(ierr);
357cb94ee6SHong Zhang 
367cb94ee6SHong Zhang   /* jok - TRUE means reuse current Jacobian else recompute Jacobian */
377cb94ee6SHong Zhang   if (jok) {
387cb94ee6SHong Zhang     ierr     = MatCopy(cvode->pmat,Jac,str);CHKERRQ(ierr);
397cb94ee6SHong Zhang     *jcurPtr = FALSE;
407cb94ee6SHong Zhang   } else {
417cb94ee6SHong Zhang     /* make PETSc vector yy point to SUNDIALS vector y */
427cb94ee6SHong Zhang     y_data = (PetscScalar *) N_VGetArrayPointer(y);
437cb94ee6SHong Zhang     ierr   = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
444ac7836bSHong Zhang 
457cb94ee6SHong Zhang     /* compute the Jacobian */
467cb94ee6SHong Zhang     ierr = TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);CHKERRQ(ierr);
477cb94ee6SHong Zhang     ierr = VecResetArray(yy); CHKERRQ(ierr);
484ac7836bSHong Zhang 
497cb94ee6SHong Zhang     /* copy the Jacobian matrix */
507cb94ee6SHong Zhang     if (!cvode->pmat) {
517cb94ee6SHong Zhang       ierr = MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);CHKERRQ(ierr);
527cb94ee6SHong Zhang       ierr = PetscLogObjectParent(ts,cvode->pmat);CHKERRQ(ierr);
532c823083SHong Zhang     } else {
547cb94ee6SHong Zhang       ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr);
557cb94ee6SHong Zhang     }
567cb94ee6SHong Zhang     *jcurPtr = TRUE;
577cb94ee6SHong Zhang   }
587cb94ee6SHong Zhang 
597cb94ee6SHong Zhang   /* construct I-gamma*Jac  */
607cb94ee6SHong Zhang   gm   = -_gamma;
617cb94ee6SHong Zhang   ierr = MatScale(Jac,gm);CHKERRQ(ierr);
627cb94ee6SHong Zhang   ierr = MatShift(Jac,one);CHKERRQ(ierr);
637cb94ee6SHong Zhang 
64dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
657cb94ee6SHong Zhang   ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr);
667cb94ee6SHong Zhang   PetscFunctionReturn(0);
677cb94ee6SHong Zhang }
687cb94ee6SHong Zhang 
697cb94ee6SHong Zhang /*
707cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
717cb94ee6SHong Zhang */
727cb94ee6SHong Zhang #undef __FUNCT__
737cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials"
744ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
754ac7836bSHong Zhang                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
767cb94ee6SHong Zhang {
777cb94ee6SHong Zhang   TS              ts = (TS) P_data;
787cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
79dcbc6d53SSean Farley   PC              pc;
807cb94ee6SHong Zhang   Vec             rr = cvode->w1,zz = cvode->w2;
817cb94ee6SHong Zhang   PetscErrorCode  ierr;
827cb94ee6SHong Zhang   PetscScalar     *r_data,*z_data;
837cb94ee6SHong Zhang 
847cb94ee6SHong Zhang   PetscFunctionBegin;
857cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
867cb94ee6SHong Zhang   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
877cb94ee6SHong Zhang   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
887cb94ee6SHong Zhang   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
897cb94ee6SHong Zhang   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
904ac7836bSHong Zhang 
917cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
92dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
937cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
947cb94ee6SHong Zhang   ierr = VecResetArray(rr); CHKERRQ(ierr);
957cb94ee6SHong Zhang   ierr = VecResetArray(zz); CHKERRQ(ierr);
967cb94ee6SHong Zhang   PetscFunctionReturn(0);
977cb94ee6SHong Zhang }
987cb94ee6SHong Zhang 
997cb94ee6SHong Zhang /*
1007cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
1017cb94ee6SHong Zhang */
1027cb94ee6SHong Zhang #undef __FUNCT__
1037cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials"
1044ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
1057cb94ee6SHong Zhang {
1067cb94ee6SHong Zhang   TS              ts = (TS) ctx;
1077adad957SLisandro Dalcin   MPI_Comm        comm = ((PetscObject)ts)->comm;
1087cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
1097cb94ee6SHong Zhang   Vec             yy = cvode->w1,yyd = cvode->w2;
1107cb94ee6SHong Zhang   PetscScalar     *y_data,*ydot_data;
1117cb94ee6SHong Zhang   PetscErrorCode  ierr;
1127cb94ee6SHong Zhang 
1137cb94ee6SHong Zhang   PetscFunctionBegin;
1145ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
1157cb94ee6SHong Zhang   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
1167cb94ee6SHong Zhang   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
1174ac7836bSHong Zhang   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
1184ac7836bSHong Zhang   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
1194ac7836bSHong Zhang 
1207cb94ee6SHong Zhang   /* now compute the right hand side function */
1217cb94ee6SHong Zhang   ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr);
1227cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
1237cb94ee6SHong Zhang   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
1244ac7836bSHong Zhang   PetscFunctionReturn(0);
1257cb94ee6SHong Zhang }
1267cb94ee6SHong Zhang 
1277cb94ee6SHong Zhang /*
128214bc6a2SJed Brown        TSStep_Sundials - Calls Sundials to integrate the ODE.
1297cb94ee6SHong Zhang */
1307cb94ee6SHong Zhang #undef __FUNCT__
131214bc6a2SJed Brown #define __FUNCT__ "TSStep_Sundials"
132193ac0bcSJed Brown PetscErrorCode TSStep_Sundials(TS ts)
1337cb94ee6SHong Zhang {
1347cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1357cb94ee6SHong Zhang   Vec            sol = ts->vec_sol;
1367cb94ee6SHong Zhang   PetscErrorCode ierr;
137186e87acSLisandro Dalcin   PetscInt       i,flag;
1389f94935aSHong Zhang   long int       its,nsteps;
1397cb94ee6SHong Zhang   realtype       t,tout;
1407cb94ee6SHong Zhang   PetscScalar    *y_data;
1417cb94ee6SHong Zhang   void           *mem;
1427cb94ee6SHong Zhang 
1437cb94ee6SHong Zhang   PetscFunctionBegin;
14416016685SHong Zhang   mem  = cvode->mem;
1457cb94ee6SHong Zhang   tout = ts->max_time;
1467cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1477cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data,cvode->y);
1487cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
149186e87acSLisandro Dalcin 
150186e87acSLisandro Dalcin   for (i = 0; i < ts->max_steps; i++) {
1517cb94ee6SHong Zhang     if (ts->ptime >= ts->max_time) break;
1523f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
153186e87acSLisandro Dalcin 
1542bfc04deSHong Zhang     if (cvode->monitorstep) {
15528adb3f7SHong Zhang       flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
1562bfc04deSHong Zhang     } else {
1572bfc04deSHong Zhang       flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
1582bfc04deSHong Zhang     }
1599f94935aSHong Zhang 
1609f94935aSHong Zhang     if (flag){ /* display error message */
1619f94935aSHong Zhang       switch (flag){
1629f94935aSHong Zhang       case CV_ILL_INPUT:
1639f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
1649f94935aSHong Zhang         break;
1659f94935aSHong Zhang       case CV_TOO_CLOSE:
1669f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
1679f94935aSHong Zhang         break;
1683c7fefeeSJed Brown       case CV_TOO_MUCH_WORK: {
1699f94935aSHong Zhang         PetscReal      tcur;
1709f94935aSHong Zhang         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
1719f94935aSHong Zhang         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
1729f94935aSHong 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);
1733c7fefeeSJed Brown       } break;
1749f94935aSHong Zhang       case CV_TOO_MUCH_ACC:
1759f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
1769f94935aSHong Zhang         break;
1779f94935aSHong Zhang       case CV_ERR_FAILURE:
1789f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
1799f94935aSHong Zhang         break;
1809f94935aSHong Zhang       case CV_CONV_FAILURE:
1819f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
1829f94935aSHong Zhang         break;
1839f94935aSHong Zhang       case CV_LINIT_FAIL:
1849f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
1859f94935aSHong Zhang         break;
1869f94935aSHong Zhang       case CV_LSETUP_FAIL:
1879f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
1889f94935aSHong Zhang         break;
1899f94935aSHong Zhang       case CV_LSOLVE_FAIL:
1909f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
1919f94935aSHong Zhang         break;
1929f94935aSHong Zhang       case CV_RHSFUNC_FAIL:
1939f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
1949f94935aSHong Zhang         break;
1959f94935aSHong Zhang       case CV_FIRST_RHSFUNC_ERR:
1969f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
1979f94935aSHong Zhang         break;
1989f94935aSHong Zhang       case CV_REPTD_RHSFUNC_ERR:
1999f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
2009f94935aSHong Zhang         break;
2019f94935aSHong Zhang       case CV_UNREC_RHSFUNC_ERR:
2029f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
2039f94935aSHong Zhang         break;
2049f94935aSHong Zhang       case CV_RTFUNC_FAIL:
2059f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
2069f94935aSHong Zhang         break;
2079f94935aSHong Zhang       default:
2089f94935aSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
2099f94935aSHong Zhang       }
2109f94935aSHong Zhang     }
2119f94935aSHong Zhang 
2127cb94ee6SHong Zhang     if (t > ts->max_time && cvode->exact_final_time) {
2137cb94ee6SHong Zhang       /* interpolate to final requested time */
2147cb94ee6SHong Zhang       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
2157cb94ee6SHong Zhang       t = tout;
2167cb94ee6SHong Zhang     }
2177cb94ee6SHong Zhang 
2187cb94ee6SHong Zhang     /* copy the solution from cvode->y to cvode->update and sol */
2197cb94ee6SHong Zhang     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
2207cb94ee6SHong Zhang     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
2217cb94ee6SHong Zhang     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
2227cb94ee6SHong Zhang     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
2237cb94ee6SHong Zhang     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
2244ac7836bSHong Zhang     ierr = CVSpilsGetNumLinIters(mem, &its);
225186e87acSLisandro Dalcin     ts->nonlinear_its = its; ts->linear_its = its;
226186e87acSLisandro Dalcin 
227186e87acSLisandro Dalcin     ts->time_step = t - ts->ptime;
228186e87acSLisandro Dalcin     ts->ptime     = t;
2297cb94ee6SHong Zhang     ts->steps++;
230186e87acSLisandro Dalcin 
2313f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
2327cb94ee6SHong Zhang     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
2337cb94ee6SHong Zhang   }
2349f94935aSHong Zhang   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
2357cb94ee6SHong Zhang   PetscFunctionReturn(0);
2367cb94ee6SHong Zhang }
2377cb94ee6SHong Zhang 
2387cb94ee6SHong Zhang #undef __FUNCT__
239277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials"
240277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts)
241277b19d0SLisandro Dalcin {
242277b19d0SLisandro Dalcin   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
243277b19d0SLisandro Dalcin   PetscErrorCode ierr;
244277b19d0SLisandro Dalcin 
245277b19d0SLisandro Dalcin   PetscFunctionBegin;
246*5c6b4a3dSJed Brown   ierr = MatDestroy(&cvode->pmat);CHKERRQ(ierr);
247*5c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
248*5c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->func);CHKERRQ(ierr);
249*5c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->rhs);CHKERRQ(ierr);
250*5c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
251*5c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
252277b19d0SLisandro Dalcin   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);
267335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr);
268335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C","",PETSC_NULL);CHKERRQ(ierr);
269335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
270335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr);
271335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
272335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
273335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
274335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr);
275335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr);
276335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C","",PETSC_NULL);CHKERRQ(ierr);
277335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr);
2787cb94ee6SHong Zhang   PetscFunctionReturn(0);
2797cb94ee6SHong Zhang }
2807cb94ee6SHong Zhang 
2817cb94ee6SHong Zhang #undef __FUNCT__
282214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials"
283214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts)
2847cb94ee6SHong Zhang {
2857cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2867cb94ee6SHong Zhang   PetscErrorCode ierr;
28716016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2887cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
28916016685SHong Zhang   void           *mem;
290dcbc6d53SSean Farley   PC             pc;
29116016685SHong Zhang   const PCType   pctype;
292ace3abfcSBarry Smith   PetscBool      pcnone;
29316016685SHong Zhang   Vec            sol = ts->vec_sol;
2947cb94ee6SHong Zhang 
2957cb94ee6SHong Zhang   PetscFunctionBegin;
2969f94935aSHong Zhang 
2977cb94ee6SHong Zhang   /* get the vector size */
2987cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2997cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
3007cb94ee6SHong Zhang 
3017cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
302a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
303e32f2f54SBarry Smith   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
3047cb94ee6SHong Zhang 
30528adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
3067cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
3077cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
3087cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
3097cb94ee6SHong Zhang   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
3107cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
3117cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
3127cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr);
3137cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
3147cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr);
3157cb94ee6SHong Zhang 
3167cb94ee6SHong Zhang   /*
3177cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
3187cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
3197cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
3207cb94ee6SHong Zhang   */
3217adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
3227adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
3237cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
3247cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
32516016685SHong Zhang 
32616016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
32716016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
328e32f2f54SBarry Smith   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
32916016685SHong Zhang   cvode->mem = mem;
33016016685SHong Zhang 
33116016685SHong Zhang   /* Set the pointer to user-defined data */
33216016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
333e32f2f54SBarry Smith   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
33416016685SHong Zhang 
335fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
336fc6b9e64SJed Brown   flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step);
337fc6b9e64SJed Brown   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
338f1cd61daSJed Brown   if (cvode->mindt > 0) {
339f1cd61daSJed Brown     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
3409f94935aSHong Zhang     if (flag){
3419f94935aSHong Zhang       if (flag == CV_MEM_NULL){
3429f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
3439f94935aSHong Zhang       } else if (flag == CV_ILL_INPUT){
3449f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
3459f94935aSHong Zhang       } else {
3469f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
3479f94935aSHong Zhang       }
3489f94935aSHong Zhang     }
349f1cd61daSJed Brown   }
350f1cd61daSJed Brown   if (cvode->maxdt > 0) {
351f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
352f1cd61daSJed Brown     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
353f1cd61daSJed Brown   }
354f1cd61daSJed Brown 
35516016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
35616016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
35716016685SHong Zhang    * the initial dependent variable vector cvode->y */
35816016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
35916016685SHong Zhang   if (flag){
360e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
36116016685SHong Zhang   }
36216016685SHong Zhang 
3639f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
36416016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
36516016685SHong Zhang   if (flag){
366e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
36716016685SHong Zhang   }
36816016685SHong Zhang 
3699f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3709f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
37116016685SHong Zhang   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
37216016685SHong Zhang 
37316016685SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
37416016685SHong Zhang   /* setup the ode integrator with the given preconditioner */
375dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
376dcbc6d53SSean Farley   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
377dcbc6d53SSean Farley   ierr = PetscTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
37816016685SHong Zhang   if (pcnone){
37916016685SHong Zhang     flag  = CVSpgmr(mem,PREC_NONE,0);
380e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
38116016685SHong Zhang   } else {
38216016685SHong Zhang     flag  = CVSpgmr(mem,PREC_LEFT,0);
383e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
38416016685SHong Zhang 
38516016685SHong Zhang     /* Set preconditioner and solve routines Precond and PSolve,
38616016685SHong Zhang      and the pointer to the user-defined block data */
38716016685SHong Zhang     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
388e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
38916016685SHong Zhang   }
39016016685SHong Zhang 
39116016685SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
392e32f2f54SBarry Smith   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
3937cb94ee6SHong Zhang   PetscFunctionReturn(0);
3947cb94ee6SHong Zhang }
3957cb94ee6SHong Zhang 
3966fadb2cdSHong Zhang /* type of CVODE linear multistep method */
397dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
3986fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
399dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
400a04cf4d8SBarry Smith 
4017cb94ee6SHong Zhang #undef __FUNCT__
402214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials"
403214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts)
4047cb94ee6SHong Zhang {
4057cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4067cb94ee6SHong Zhang   PetscErrorCode ierr;
4077cb94ee6SHong Zhang   int            indx;
408ace3abfcSBarry Smith   PetscBool      flag;
4097cb94ee6SHong Zhang 
4107cb94ee6SHong Zhang   PetscFunctionBegin;
4117cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
4126fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
4137cb94ee6SHong Zhang     if (flag) {
4146fadb2cdSHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
4157cb94ee6SHong Zhang     }
4166fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
4177cb94ee6SHong Zhang     if (flag) {
4187cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
4197cb94ee6SHong Zhang     }
4207cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
4217cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
422f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
423f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
4247cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
4257cb94ee6SHong Zhang     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
426acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_sundials_exact_final_time","Interpolate output to stop exactly at the final time","TSSundialsSetExactFinalTime",cvode->exact_final_time,&cvode->exact_final_time,PETSC_NULL);CHKERRQ(ierr);
427acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
4287cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4297cb94ee6SHong Zhang   PetscFunctionReturn(0);
4307cb94ee6SHong Zhang }
4317cb94ee6SHong Zhang 
4327cb94ee6SHong Zhang #undef __FUNCT__
4337cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
4347cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
4357cb94ee6SHong Zhang {
4367cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4377cb94ee6SHong Zhang   PetscErrorCode ierr;
4387cb94ee6SHong Zhang   char           *type;
4397cb94ee6SHong Zhang   char           atype[] = "Adams";
4407cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
441ace3abfcSBarry Smith   PetscBool      iascii,isstring;
4422c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
4432c823083SHong Zhang   PetscInt       qlast,qcur;
4445d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
4457cb94ee6SHong Zhang 
4467cb94ee6SHong Zhang   PetscFunctionBegin;
4477cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
4487cb94ee6SHong Zhang   else                                     {type = btype;}
4497cb94ee6SHong Zhang 
4502692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4512692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
4527cb94ee6SHong Zhang   if (iascii) {
4537cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
4547cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
4557cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
4567cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
4577cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
4587cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4597cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4607cb94ee6SHong Zhang     } else {
4617cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4627cb94ee6SHong Zhang     }
463f1cd61daSJed Brown     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
464f1cd61daSJed Brown     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
4652c823083SHong Zhang 
4665d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4675d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
4685d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
4692c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
4702c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
4712c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
4722c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
4732c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
4742c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
4752c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
4762c823083SHong Zhang 
4772c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
4782c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
4792c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
4802c823083SHong Zhang 
4812c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
4822c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
4832c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
4842c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
4852c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4862c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
4872c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
4882c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
4892c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4902c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
4912c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4922c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
4937cb94ee6SHong Zhang   } else if (isstring) {
4947cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
4957cb94ee6SHong Zhang   } else {
496e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
4977cb94ee6SHong Zhang   }
4987cb94ee6SHong Zhang   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
4997cb94ee6SHong Zhang   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
5007cb94ee6SHong Zhang   PetscFunctionReturn(0);
5017cb94ee6SHong Zhang }
5027cb94ee6SHong Zhang 
5037cb94ee6SHong Zhang 
5047cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
5057cb94ee6SHong Zhang EXTERN_C_BEGIN
5067cb94ee6SHong Zhang #undef __FUNCT__
5077cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
5087087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
5097cb94ee6SHong Zhang {
5107cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5117cb94ee6SHong Zhang 
5127cb94ee6SHong Zhang   PetscFunctionBegin;
5137cb94ee6SHong Zhang   cvode->cvode_type = type;
5147cb94ee6SHong Zhang   PetscFunctionReturn(0);
5157cb94ee6SHong Zhang }
5167cb94ee6SHong Zhang EXTERN_C_END
5177cb94ee6SHong Zhang 
5187cb94ee6SHong Zhang EXTERN_C_BEGIN
5197cb94ee6SHong Zhang #undef __FUNCT__
5207cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
5217087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
5227cb94ee6SHong Zhang {
5237cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5247cb94ee6SHong Zhang 
5257cb94ee6SHong Zhang   PetscFunctionBegin;
5267cb94ee6SHong Zhang   cvode->restart = restart;
5277cb94ee6SHong Zhang   PetscFunctionReturn(0);
5287cb94ee6SHong Zhang }
5297cb94ee6SHong Zhang EXTERN_C_END
5307cb94ee6SHong Zhang 
5317cb94ee6SHong Zhang EXTERN_C_BEGIN
5327cb94ee6SHong Zhang #undef __FUNCT__
5337cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
5347087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
5357cb94ee6SHong Zhang {
5367cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5377cb94ee6SHong Zhang 
5387cb94ee6SHong Zhang   PetscFunctionBegin;
5397cb94ee6SHong Zhang   cvode->linear_tol = tol;
5407cb94ee6SHong Zhang   PetscFunctionReturn(0);
5417cb94ee6SHong Zhang }
5427cb94ee6SHong Zhang EXTERN_C_END
5437cb94ee6SHong Zhang 
5447cb94ee6SHong Zhang EXTERN_C_BEGIN
5457cb94ee6SHong Zhang #undef __FUNCT__
5467cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
5477087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
5487cb94ee6SHong Zhang {
5497cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5507cb94ee6SHong Zhang 
5517cb94ee6SHong Zhang   PetscFunctionBegin;
5527cb94ee6SHong Zhang   cvode->gtype = type;
5537cb94ee6SHong Zhang   PetscFunctionReturn(0);
5547cb94ee6SHong Zhang }
5557cb94ee6SHong Zhang EXTERN_C_END
5567cb94ee6SHong Zhang 
5577cb94ee6SHong Zhang EXTERN_C_BEGIN
5587cb94ee6SHong Zhang #undef __FUNCT__
5597cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
5607087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
5617cb94ee6SHong Zhang {
5627cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5637cb94ee6SHong Zhang 
5647cb94ee6SHong Zhang   PetscFunctionBegin;
5657cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5667cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
5677cb94ee6SHong Zhang   PetscFunctionReturn(0);
5687cb94ee6SHong Zhang }
5697cb94ee6SHong Zhang EXTERN_C_END
5707cb94ee6SHong Zhang 
5717cb94ee6SHong Zhang EXTERN_C_BEGIN
5727cb94ee6SHong Zhang #undef __FUNCT__
573f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
5747087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
575f1cd61daSJed Brown {
576f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
577f1cd61daSJed Brown 
578f1cd61daSJed Brown   PetscFunctionBegin;
579f1cd61daSJed Brown   cvode->mindt = mindt;
580f1cd61daSJed Brown   PetscFunctionReturn(0);
581f1cd61daSJed Brown }
582f1cd61daSJed Brown EXTERN_C_END
583f1cd61daSJed Brown 
584f1cd61daSJed Brown EXTERN_C_BEGIN
585f1cd61daSJed Brown #undef __FUNCT__
586f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
5877087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
588f1cd61daSJed Brown {
589f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
590f1cd61daSJed Brown 
591f1cd61daSJed Brown   PetscFunctionBegin;
592f1cd61daSJed Brown   cvode->maxdt = maxdt;
593f1cd61daSJed Brown   PetscFunctionReturn(0);
594f1cd61daSJed Brown }
595f1cd61daSJed Brown EXTERN_C_END
596f1cd61daSJed Brown 
597f1cd61daSJed Brown EXTERN_C_BEGIN
598f1cd61daSJed Brown #undef __FUNCT__
5997cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
6007087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
6017cb94ee6SHong Zhang {
602dcbc6d53SSean Farley   SNES            snes;
603dcbc6d53SSean Farley   KSP             ksp;
604dcbc6d53SSean Farley   PetscErrorCode  ierr;
6057cb94ee6SHong Zhang 
6067cb94ee6SHong Zhang   PetscFunctionBegin;
607dcbc6d53SSean Farley   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
608dcbc6d53SSean Farley   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
609dcbc6d53SSean Farley   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
6107cb94ee6SHong Zhang   PetscFunctionReturn(0);
6117cb94ee6SHong Zhang }
6127cb94ee6SHong Zhang EXTERN_C_END
6137cb94ee6SHong Zhang 
6147cb94ee6SHong Zhang EXTERN_C_BEGIN
6157cb94ee6SHong Zhang #undef __FUNCT__
6167cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
6177087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
6187cb94ee6SHong Zhang {
6197cb94ee6SHong Zhang   PetscFunctionBegin;
6202c823083SHong Zhang   if (nonlin) *nonlin = ts->nonlinear_its;
6212c823083SHong Zhang   if (lin)    *lin    = ts->linear_its;
6227cb94ee6SHong Zhang   PetscFunctionReturn(0);
6237cb94ee6SHong Zhang }
6247cb94ee6SHong Zhang EXTERN_C_END
6257cb94ee6SHong Zhang 
6267cb94ee6SHong Zhang EXTERN_C_BEGIN
6277cb94ee6SHong Zhang #undef __FUNCT__
6287cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
6297087cfbeSBarry Smith PetscErrorCode  TSSundialsSetExactFinalTime_Sundials(TS ts,PetscBool  s)
6307cb94ee6SHong Zhang {
6317cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6327cb94ee6SHong Zhang 
6337cb94ee6SHong Zhang   PetscFunctionBegin;
6347cb94ee6SHong Zhang   cvode->exact_final_time = s;
6357cb94ee6SHong Zhang   PetscFunctionReturn(0);
6367cb94ee6SHong Zhang }
6377cb94ee6SHong Zhang EXTERN_C_END
6382bfc04deSHong Zhang 
6392bfc04deSHong Zhang EXTERN_C_BEGIN
6402bfc04deSHong Zhang #undef __FUNCT__
6412bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
6427087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
6432bfc04deSHong Zhang {
6442bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6452bfc04deSHong Zhang 
6462bfc04deSHong Zhang   PetscFunctionBegin;
6472bfc04deSHong Zhang   cvode->monitorstep = s;
6482bfc04deSHong Zhang   PetscFunctionReturn(0);
6492bfc04deSHong Zhang }
6502bfc04deSHong Zhang EXTERN_C_END
6517cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
6527cb94ee6SHong Zhang 
6537cb94ee6SHong Zhang #undef __FUNCT__
6547cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
6557cb94ee6SHong Zhang /*@C
6567cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
6577cb94ee6SHong Zhang 
6587cb94ee6SHong Zhang    Not Collective
6597cb94ee6SHong Zhang 
6607cb94ee6SHong Zhang    Input parameters:
6617cb94ee6SHong Zhang .    ts     - the time-step context
6627cb94ee6SHong Zhang 
6637cb94ee6SHong Zhang    Output Parameters:
6647cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
6657cb94ee6SHong Zhang -   lin    - number of linear iterations
6667cb94ee6SHong Zhang 
6677cb94ee6SHong Zhang    Level: advanced
6687cb94ee6SHong Zhang 
6697cb94ee6SHong Zhang    Notes:
6707cb94ee6SHong Zhang     These return the number since the creation of the TS object
6717cb94ee6SHong Zhang 
6727cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
6737cb94ee6SHong Zhang 
6747cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6757cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6767cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
677f1cd61daSJed Brown           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSundialsSetExactFinalTime()
6787cb94ee6SHong Zhang 
6797cb94ee6SHong Zhang @*/
6807087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
6817cb94ee6SHong Zhang {
6824ac538c5SBarry Smith   PetscErrorCode ierr;
6837cb94ee6SHong Zhang 
6847cb94ee6SHong Zhang   PetscFunctionBegin;
6854ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
6867cb94ee6SHong Zhang   PetscFunctionReturn(0);
6877cb94ee6SHong Zhang }
6887cb94ee6SHong Zhang 
6897cb94ee6SHong Zhang #undef __FUNCT__
6907cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
6917cb94ee6SHong Zhang /*@
6927cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
6937cb94ee6SHong Zhang 
6943f9fe445SBarry Smith    Logically Collective on TS
6957cb94ee6SHong Zhang 
6967cb94ee6SHong Zhang    Input parameters:
6977cb94ee6SHong Zhang +    ts     - the time-step context
6987cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
6997cb94ee6SHong Zhang 
7007cb94ee6SHong Zhang    Level: intermediate
7017cb94ee6SHong Zhang 
7027cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
7037cb94ee6SHong Zhang 
7047cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
7057cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7067cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7077cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7087cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7097cb94ee6SHong Zhang @*/
7107087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
7117cb94ee6SHong Zhang {
7124ac538c5SBarry Smith   PetscErrorCode ierr;
7137cb94ee6SHong Zhang 
7147cb94ee6SHong Zhang   PetscFunctionBegin;
7154ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
7167cb94ee6SHong Zhang   PetscFunctionReturn(0);
7177cb94ee6SHong Zhang }
7187cb94ee6SHong Zhang 
7197cb94ee6SHong Zhang #undef __FUNCT__
7207cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart"
7217cb94ee6SHong Zhang /*@
7227cb94ee6SHong Zhang    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
7237cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
7247cb94ee6SHong Zhang        this is ALSO the maximum number of GMRES steps that will be used.
7257cb94ee6SHong Zhang 
7263f9fe445SBarry Smith    Logically Collective on TS
7277cb94ee6SHong Zhang 
7287cb94ee6SHong Zhang    Input parameters:
7297cb94ee6SHong Zhang +    ts      - the time-step context
7307cb94ee6SHong Zhang -    restart - number of direction vectors (the restart size).
7317cb94ee6SHong Zhang 
7327cb94ee6SHong Zhang    Level: advanced
7337cb94ee6SHong Zhang 
7347cb94ee6SHong Zhang .keywords: GMRES, restart
7357cb94ee6SHong Zhang 
7367cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
7377cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7387cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7397cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7407cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7417cb94ee6SHong Zhang 
7427cb94ee6SHong Zhang @*/
7437087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart(TS ts,int restart)
7447cb94ee6SHong Zhang {
7454ac538c5SBarry Smith   PetscErrorCode ierr;
7467cb94ee6SHong Zhang 
7477cb94ee6SHong Zhang   PetscFunctionBegin;
748c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,restart,2);
7494ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr);
7507cb94ee6SHong Zhang   PetscFunctionReturn(0);
7517cb94ee6SHong Zhang }
7527cb94ee6SHong Zhang 
7537cb94ee6SHong Zhang #undef __FUNCT__
7547cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
7557cb94ee6SHong Zhang /*@
7567cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
7577cb94ee6SHong Zhang        system by SUNDIALS.
7587cb94ee6SHong Zhang 
7593f9fe445SBarry Smith    Logically Collective on TS
7607cb94ee6SHong Zhang 
7617cb94ee6SHong Zhang    Input parameters:
7627cb94ee6SHong Zhang +    ts     - the time-step context
7637cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
7647cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
7657cb94ee6SHong Zhang 
7667cb94ee6SHong Zhang    Level: advanced
7677cb94ee6SHong Zhang 
7687cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
7697cb94ee6SHong Zhang 
7707cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7717cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7727cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7737cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7747cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7757cb94ee6SHong Zhang 
7767cb94ee6SHong Zhang @*/
7777087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
7787cb94ee6SHong Zhang {
7794ac538c5SBarry Smith   PetscErrorCode ierr;
7807cb94ee6SHong Zhang 
7817cb94ee6SHong Zhang   PetscFunctionBegin;
782c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,tol,2);
7834ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
7847cb94ee6SHong Zhang   PetscFunctionReturn(0);
7857cb94ee6SHong Zhang }
7867cb94ee6SHong Zhang 
7877cb94ee6SHong Zhang #undef __FUNCT__
7887cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
7897cb94ee6SHong Zhang /*@
7907cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
7917cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
7927cb94ee6SHong Zhang 
7933f9fe445SBarry Smith    Logically Collective on TS
7947cb94ee6SHong Zhang 
7957cb94ee6SHong Zhang    Input parameters:
7967cb94ee6SHong Zhang +    ts  - the time-step context
7977cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
7987cb94ee6SHong Zhang 
7997cb94ee6SHong Zhang    Level: advanced
8007cb94ee6SHong Zhang 
8017cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
8027cb94ee6SHong Zhang 
8037cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8047cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
8057cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8067cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
8077cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
8087cb94ee6SHong Zhang 
8097cb94ee6SHong Zhang @*/
8107087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
8117cb94ee6SHong Zhang {
8124ac538c5SBarry Smith   PetscErrorCode ierr;
8137cb94ee6SHong Zhang 
8147cb94ee6SHong Zhang   PetscFunctionBegin;
8154ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
8167cb94ee6SHong Zhang   PetscFunctionReturn(0);
8177cb94ee6SHong Zhang }
8187cb94ee6SHong Zhang 
8197cb94ee6SHong Zhang #undef __FUNCT__
8207cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
8217cb94ee6SHong Zhang /*@
8227cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
8237cb94ee6SHong Zhang                          Sundials for error control.
8247cb94ee6SHong Zhang 
8253f9fe445SBarry Smith    Logically Collective on TS
8267cb94ee6SHong Zhang 
8277cb94ee6SHong Zhang    Input parameters:
8287cb94ee6SHong Zhang +    ts  - the time-step context
8297cb94ee6SHong Zhang .    aabs - the absolute tolerance
8307cb94ee6SHong Zhang -    rel - the relative tolerance
8317cb94ee6SHong Zhang 
8327cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
8337cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
8347cb94ee6SHong Zhang 
8357cb94ee6SHong Zhang    Level: intermediate
8367cb94ee6SHong Zhang 
8377cb94ee6SHong Zhang .keywords: Sundials, tolerance
8387cb94ee6SHong Zhang 
8397cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8407cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
8417cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8427cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
8437cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
8447cb94ee6SHong Zhang 
8457cb94ee6SHong Zhang @*/
8467087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
8477cb94ee6SHong Zhang {
8484ac538c5SBarry Smith   PetscErrorCode ierr;
8497cb94ee6SHong Zhang 
8507cb94ee6SHong Zhang   PetscFunctionBegin;
8514ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
8527cb94ee6SHong Zhang   PetscFunctionReturn(0);
8537cb94ee6SHong Zhang }
8547cb94ee6SHong Zhang 
8557cb94ee6SHong Zhang #undef __FUNCT__
8567cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
8577cb94ee6SHong Zhang /*@
8587cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
8597cb94ee6SHong Zhang 
8607cb94ee6SHong Zhang    Input Parameter:
8617cb94ee6SHong Zhang .    ts - the time-step context
8627cb94ee6SHong Zhang 
8637cb94ee6SHong Zhang    Output Parameter:
8647cb94ee6SHong Zhang .    pc - the preconditioner context
8657cb94ee6SHong Zhang 
8667cb94ee6SHong Zhang    Level: advanced
8677cb94ee6SHong Zhang 
8687cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8697cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8707cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8717cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
8727cb94ee6SHong Zhang @*/
8737087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
8747cb94ee6SHong Zhang {
8754ac538c5SBarry Smith   PetscErrorCode ierr;
8767cb94ee6SHong Zhang 
8777cb94ee6SHong Zhang   PetscFunctionBegin;
8784ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
8797cb94ee6SHong Zhang   PetscFunctionReturn(0);
8807cb94ee6SHong Zhang }
8817cb94ee6SHong Zhang 
8827cb94ee6SHong Zhang #undef __FUNCT__
8837cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime"
8847cb94ee6SHong Zhang /*@
8857cb94ee6SHong Zhang    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
8867cb94ee6SHong Zhang       exact final time requested by the user or just returns it at the final time
8877cb94ee6SHong Zhang       it computed. (Defaults to true).
8887cb94ee6SHong Zhang 
8897cb94ee6SHong Zhang    Input Parameter:
8907cb94ee6SHong Zhang +   ts - the time-step context
8917cb94ee6SHong Zhang -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
8927cb94ee6SHong Zhang 
8937cb94ee6SHong Zhang    Level: beginner
8947cb94ee6SHong Zhang 
8957cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8967cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8977cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8987cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
8997cb94ee6SHong Zhang @*/
9007087cfbeSBarry Smith PetscErrorCode  TSSundialsSetExactFinalTime(TS ts,PetscBool  ft)
9017cb94ee6SHong Zhang {
9024ac538c5SBarry Smith   PetscErrorCode ierr;
9037cb94ee6SHong Zhang 
9047cb94ee6SHong Zhang   PetscFunctionBegin;
9054ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetExactFinalTime_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
9067cb94ee6SHong Zhang   PetscFunctionReturn(0);
9077cb94ee6SHong Zhang }
9087cb94ee6SHong Zhang 
9092bfc04deSHong Zhang #undef __FUNCT__
910f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep"
911f1cd61daSJed Brown /*@
912f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
913f1cd61daSJed Brown 
914f1cd61daSJed Brown    Input Parameter:
915f1cd61daSJed Brown +   ts - the time-step context
916f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
917f1cd61daSJed Brown 
918fc6b9e64SJed Brown    Note:
919fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
920fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
921fc6b9e64SJed Brown 
922f1cd61daSJed Brown    Level: beginner
923f1cd61daSJed Brown 
924f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
925f1cd61daSJed Brown @*/
9267087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
927f1cd61daSJed Brown {
9284ac538c5SBarry Smith   PetscErrorCode ierr;
929f1cd61daSJed Brown 
930f1cd61daSJed Brown   PetscFunctionBegin;
9314ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
932f1cd61daSJed Brown   PetscFunctionReturn(0);
933f1cd61daSJed Brown }
934f1cd61daSJed Brown 
935f1cd61daSJed Brown #undef __FUNCT__
936f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep"
937f1cd61daSJed Brown /*@
938f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
939f1cd61daSJed Brown 
940f1cd61daSJed Brown    Input Parameter:
941f1cd61daSJed Brown +   ts - the time-step context
942f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
943f1cd61daSJed Brown 
944f1cd61daSJed Brown    Level: beginner
945f1cd61daSJed Brown 
946f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
947f1cd61daSJed Brown @*/
9487087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
949f1cd61daSJed Brown {
9504ac538c5SBarry Smith   PetscErrorCode ierr;
951f1cd61daSJed Brown 
952f1cd61daSJed Brown   PetscFunctionBegin;
9534ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
954f1cd61daSJed Brown   PetscFunctionReturn(0);
955f1cd61daSJed Brown }
956f1cd61daSJed Brown 
957f1cd61daSJed Brown #undef __FUNCT__
9582bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps"
9592bfc04deSHong Zhang /*@
9602bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
9612bfc04deSHong Zhang 
9622bfc04deSHong Zhang    Input Parameter:
9632bfc04deSHong Zhang +   ts - the time-step context
9642bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
9652bfc04deSHong Zhang 
9662bfc04deSHong Zhang    Level: beginner
9672bfc04deSHong Zhang 
9682bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9692bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
9702bfc04deSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9712bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
9722bfc04deSHong Zhang @*/
9737087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
9742bfc04deSHong Zhang {
9754ac538c5SBarry Smith   PetscErrorCode ierr;
9762bfc04deSHong Zhang 
9772bfc04deSHong Zhang   PetscFunctionBegin;
9784ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
9792bfc04deSHong Zhang   PetscFunctionReturn(0);
9802bfc04deSHong Zhang }
9817cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
9827cb94ee6SHong Zhang /*MC
98396f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
9847cb94ee6SHong Zhang 
9857cb94ee6SHong Zhang    Options Database:
9867cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
9877cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
9887cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
9897cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
9907cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
9917cb94ee6SHong Zhang .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
99216016685SHong Zhang .    -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time
99316016685SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
99416016685SHong Zhang 
9957cb94ee6SHong Zhang 
9967cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
9977cb94ee6SHong Zhang            only PETSc PC options
9987cb94ee6SHong Zhang 
9997cb94ee6SHong Zhang     Level: beginner
10007cb94ee6SHong Zhang 
10017cb94ee6SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
10027cb94ee6SHong Zhang            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
10037cb94ee6SHong Zhang 
10047cb94ee6SHong Zhang M*/
10057cb94ee6SHong Zhang EXTERN_C_BEGIN
10067cb94ee6SHong Zhang #undef __FUNCT__
10077cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
10087087cfbeSBarry Smith PetscErrorCode  TSCreate_Sundials(TS ts)
10097cb94ee6SHong Zhang {
10107cb94ee6SHong Zhang   TS_Sundials *cvode;
10117cb94ee6SHong Zhang   PetscErrorCode ierr;
10127cb94ee6SHong Zhang 
10137cb94ee6SHong Zhang   PetscFunctionBegin;
1014277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
101528adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
101628adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
1017214bc6a2SJed Brown   ts->ops->setup          = TSSetUp_Sundials;
1018214bc6a2SJed Brown   ts->ops->step           = TSStep_Sundials;
1019214bc6a2SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
10207cb94ee6SHong Zhang 
102138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
10227cb94ee6SHong Zhang   ts->data          = (void*)cvode;
10236fadb2cdSHong Zhang   cvode->cvode_type = SUNDIALS_BDF;
10246fadb2cdSHong Zhang   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
10257cb94ee6SHong Zhang   cvode->restart    = 5;
10267cb94ee6SHong Zhang   cvode->linear_tol = .05;
10277cb94ee6SHong Zhang 
10282bfc04deSHong Zhang   cvode->exact_final_time = PETSC_TRUE;
10292bfc04deSHong Zhang   cvode->monitorstep      = PETSC_FALSE;
10307cb94ee6SHong Zhang 
1031a07356b8SSatish Balay   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
1032f1cd61daSJed Brown 
1033f1cd61daSJed Brown   cvode->mindt = -1.;
1034f1cd61daSJed Brown   cvode->maxdt = -1.;
1035f1cd61daSJed Brown 
10367cb94ee6SHong Zhang   /* set tolerance for Sundials */
10377cb94ee6SHong Zhang   cvode->reltol = 1e-6;
10382c823083SHong Zhang   cvode->abstol = 1e-6;
10397cb94ee6SHong Zhang 
10407cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
10417cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
10427cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
10437cb94ee6SHong Zhang                     "TSSundialsSetGMRESRestart_Sundials",
10447cb94ee6SHong Zhang                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
10457cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
10467cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
10477cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
10487cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
10497cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
10507cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
10517cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
10527cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
10537cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1054f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1055f1cd61daSJed Brown                     "TSSundialsSetMinTimeStep_Sundials",
1056f1cd61daSJed Brown                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1057f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1058f1cd61daSJed Brown                     "TSSundialsSetMaxTimeStep_Sundials",
1059f1cd61daSJed Brown                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
10607cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
10617cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
10627cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
10637cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
10647cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
10657cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
10667cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
10677cb94ee6SHong Zhang                     "TSSundialsSetExactFinalTime_Sundials",
10687cb94ee6SHong Zhang                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
10692bfc04deSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
10702bfc04deSHong Zhang                     "TSSundialsMonitorInternalSteps_Sundials",
10712bfc04deSHong Zhang                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
10722bfc04deSHong Zhang 
10737cb94ee6SHong Zhang   PetscFunctionReturn(0);
10747cb94ee6SHong Zhang }
10757cb94ee6SHong Zhang EXTERN_C_END
1076