xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision cdbf8f939cdfb1c797c4b7f2cbbd00be19935363)
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;
250679a0aeSJed Brown   Mat            J,P;
260679a0aeSJed Brown   Vec            yy = cvode->w1,yydot = cvode->ydot;
270679a0aeSJed Brown   PetscReal      gm = (PetscReal)_gamma;
287cb94ee6SHong Zhang   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
297cb94ee6SHong Zhang   PetscScalar    *y_data;
307cb94ee6SHong Zhang 
317cb94ee6SHong Zhang   PetscFunctionBegin;
320679a0aeSJed Brown   ierr = TSGetIJacobian(ts,&J,&P,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
337cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(y);
347cb94ee6SHong Zhang   ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
350679a0aeSJed Brown   ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
360679a0aeSJed Brown   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
370679a0aeSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr);
387cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRQ(ierr);
390679a0aeSJed Brown   ierr = MatScale(P,gm);CHKERRQ(ierr);  /* turn into I-gm*Jrest, J is not used by Sundials  */
407cb94ee6SHong Zhang   *jcurPtr = TRUE;
41dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
420679a0aeSJed Brown   ierr = PCSetOperators(pc,J,P,str);CHKERRQ(ierr);
437cb94ee6SHong Zhang   PetscFunctionReturn(0);
447cb94ee6SHong Zhang }
457cb94ee6SHong Zhang 
467cb94ee6SHong Zhang /*
477cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
487cb94ee6SHong Zhang */
497cb94ee6SHong Zhang #undef __FUNCT__
507cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials"
514ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
524ac7836bSHong Zhang                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
537cb94ee6SHong Zhang {
547cb94ee6SHong Zhang   TS              ts = (TS) P_data;
557cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
56dcbc6d53SSean Farley   PC              pc;
577cb94ee6SHong Zhang   Vec             rr = cvode->w1,zz = cvode->w2;
587cb94ee6SHong Zhang   PetscErrorCode  ierr;
597cb94ee6SHong Zhang   PetscScalar     *r_data,*z_data;
607cb94ee6SHong Zhang 
617cb94ee6SHong Zhang   PetscFunctionBegin;
627cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
637cb94ee6SHong Zhang   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
647cb94ee6SHong Zhang   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
657cb94ee6SHong Zhang   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
667cb94ee6SHong Zhang   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
674ac7836bSHong Zhang 
687cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
69dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
707cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
717cb94ee6SHong Zhang   ierr = VecResetArray(rr); CHKERRQ(ierr);
727cb94ee6SHong Zhang   ierr = VecResetArray(zz); CHKERRQ(ierr);
737cb94ee6SHong Zhang   PetscFunctionReturn(0);
747cb94ee6SHong Zhang }
757cb94ee6SHong Zhang 
767cb94ee6SHong Zhang /*
777cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
787cb94ee6SHong Zhang */
797cb94ee6SHong Zhang #undef __FUNCT__
807cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials"
814ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
827cb94ee6SHong Zhang {
837cb94ee6SHong Zhang   TS              ts = (TS) ctx;
847adad957SLisandro Dalcin   MPI_Comm        comm = ((PetscObject)ts)->comm;
857cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
860679a0aeSJed Brown   Vec             yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
877cb94ee6SHong Zhang   PetscScalar     *y_data,*ydot_data;
887cb94ee6SHong Zhang   PetscErrorCode  ierr;
897cb94ee6SHong Zhang 
907cb94ee6SHong Zhang   PetscFunctionBegin;
915ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
927cb94ee6SHong Zhang   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
937cb94ee6SHong Zhang   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
944ac7836bSHong Zhang   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
954ac7836bSHong Zhang   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
960679a0aeSJed Brown   ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
974ac7836bSHong Zhang 
987cb94ee6SHong Zhang   /* now compute the right hand side function */
990679a0aeSJed Brown   ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr);
1000679a0aeSJed Brown   ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
1017cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
1027cb94ee6SHong Zhang   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
1034ac7836bSHong Zhang   PetscFunctionReturn(0);
1047cb94ee6SHong Zhang }
1057cb94ee6SHong Zhang 
1067cb94ee6SHong Zhang /*
107b4eba00bSSean Farley        TSStep_Sundials - Calls Sundials to integrate the ODE.
1087cb94ee6SHong Zhang */
1097cb94ee6SHong Zhang #undef __FUNCT__
110b4eba00bSSean Farley #define __FUNCT__ "TSStep_Sundials"
111b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts)
1127cb94ee6SHong Zhang {
1137cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1147cb94ee6SHong Zhang   Vec            sol = ts->vec_sol;
1157cb94ee6SHong Zhang   PetscErrorCode ierr;
116b4eba00bSSean Farley   PetscInt       flag;
1179f94935aSHong Zhang   long int       its,nsteps;
1187cb94ee6SHong Zhang   realtype       t,tout;
1197cb94ee6SHong Zhang   PetscScalar    *y_data;
1207cb94ee6SHong Zhang   void           *mem;
121abe29043SSean Farley   SNES           snes;
122abe29043SSean Farley   Vec            res; /* This, together with snes, will check if the SNES vec_func has been set */
1237cb94ee6SHong Zhang 
1247cb94ee6SHong Zhang   PetscFunctionBegin;
12516016685SHong Zhang   mem  = cvode->mem;
1267cb94ee6SHong Zhang   tout = ts->max_time;
1277cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1287cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data,cvode->y);
1297cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
130186e87acSLisandro Dalcin 
131b4eba00bSSean Farley   /* Should think about moving this outside the loop */
132abe29043SSean Farley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
133abe29043SSean Farley   ierr = SNESGetFunction(snes, &res, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
134abe29043SSean Farley   if (!res) {
135abe29043SSean Farley     ierr = TSSetIFunction(ts, sol, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
136abe29043SSean Farley   }
137abe29043SSean Farley 
1383f2090d5SJed Brown   ierr = TSPreStep(ts);CHKERRQ(ierr);
139186e87acSLisandro Dalcin 
1402bfc04deSHong Zhang   if (cvode->monitorstep) {
14128adb3f7SHong Zhang     flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
1422bfc04deSHong Zhang   } else {
1432bfc04deSHong Zhang     flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
1442bfc04deSHong Zhang   }
1459f94935aSHong Zhang 
1469f94935aSHong Zhang   if (flag){ /* display error message */
1479f94935aSHong Zhang     switch (flag){
1489f94935aSHong Zhang       case CV_ILL_INPUT:
1499f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
1509f94935aSHong Zhang         break;
1519f94935aSHong Zhang       case CV_TOO_CLOSE:
1529f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
1539f94935aSHong Zhang         break;
1543c7fefeeSJed Brown       case CV_TOO_MUCH_WORK: {
1559f94935aSHong Zhang         PetscReal      tcur;
1569f94935aSHong Zhang         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
1579f94935aSHong Zhang         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
1589f94935aSHong 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);
1593c7fefeeSJed Brown       } break;
1609f94935aSHong Zhang       case CV_TOO_MUCH_ACC:
1619f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
1629f94935aSHong Zhang         break;
1639f94935aSHong Zhang       case CV_ERR_FAILURE:
1649f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
1659f94935aSHong Zhang         break;
1669f94935aSHong Zhang       case CV_CONV_FAILURE:
1679f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
1689f94935aSHong Zhang         break;
1699f94935aSHong Zhang       case CV_LINIT_FAIL:
1709f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
1719f94935aSHong Zhang         break;
1729f94935aSHong Zhang       case CV_LSETUP_FAIL:
1739f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
1749f94935aSHong Zhang         break;
1759f94935aSHong Zhang       case CV_LSOLVE_FAIL:
1769f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
1779f94935aSHong Zhang         break;
1789f94935aSHong Zhang       case CV_RHSFUNC_FAIL:
1799f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
1809f94935aSHong Zhang         break;
1819f94935aSHong Zhang       case CV_FIRST_RHSFUNC_ERR:
1829f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
1839f94935aSHong Zhang         break;
1849f94935aSHong Zhang       case CV_REPTD_RHSFUNC_ERR:
1859f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
1869f94935aSHong Zhang         break;
1879f94935aSHong Zhang       case CV_UNREC_RHSFUNC_ERR:
1889f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
1899f94935aSHong Zhang         break;
1909f94935aSHong Zhang       case CV_RTFUNC_FAIL:
1919f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
1929f94935aSHong Zhang         break;
1939f94935aSHong Zhang       default:
1949f94935aSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
1959f94935aSHong Zhang     }
1969f94935aSHong Zhang   }
1979f94935aSHong Zhang 
1987cb94ee6SHong Zhang   /* copy the solution from cvode->y to cvode->update and sol */
1997cb94ee6SHong Zhang   ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
2007cb94ee6SHong Zhang   ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
2017cb94ee6SHong Zhang   ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
2027cb94ee6SHong Zhang   ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
2037cb94ee6SHong Zhang   ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
2044ac7836bSHong Zhang   ierr = CVSpilsGetNumLinIters(mem, &its);
205186e87acSLisandro Dalcin   ts->nonlinear_its = its; ts->linear_its = its;
206186e87acSLisandro Dalcin 
207186e87acSLisandro Dalcin   ts->time_step = t - ts->ptime;
208186e87acSLisandro Dalcin   ts->ptime     = t;
2097cb94ee6SHong Zhang   ts->steps++;
210186e87acSLisandro Dalcin 
2119f94935aSHong Zhang   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
212b4eba00bSSean Farley   if (!cvode->monitorstep) ts->steps = nsteps;
213b4eba00bSSean Farley   PetscFunctionReturn(0);
214b4eba00bSSean Farley }
215b4eba00bSSean Farley 
216b4eba00bSSean Farley #undef __FUNCT__
217b4eba00bSSean Farley #define __FUNCT__ "TSInterpolate_Sundials"
218b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
219b4eba00bSSean Farley {
220b4eba00bSSean Farley   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
221b4eba00bSSean Farley   N_Vector        y;
222b4eba00bSSean Farley   PetscErrorCode  ierr;
223b4eba00bSSean Farley   PetscScalar     *x_data;
224b4eba00bSSean Farley   PetscInt        glosize,locsize;
225b4eba00bSSean Farley 
226b4eba00bSSean Farley   PetscFunctionBegin;
227b4eba00bSSean Farley 
228b4eba00bSSean Farley   /* get the vector size */
229b4eba00bSSean Farley   ierr = VecGetSize(X,&glosize);CHKERRQ(ierr);
230b4eba00bSSean Farley   ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr);
231b4eba00bSSean Farley 
232b4eba00bSSean Farley   /* allocate the memory for N_Vec y */
233b4eba00bSSean Farley   y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
234b4eba00bSSean Farley   if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");
235b4eba00bSSean Farley 
236b4eba00bSSean Farley   ierr = VecGetArray(X,&x_data);CHKERRQ(ierr);
237b4eba00bSSean Farley   N_VSetArrayPointer((realtype *)x_data,y);
238b4eba00bSSean Farley   ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr);
239b4eba00bSSean Farley   ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr);
240b4eba00bSSean Farley 
2417cb94ee6SHong Zhang   PetscFunctionReturn(0);
2427cb94ee6SHong Zhang }
2437cb94ee6SHong Zhang 
2447cb94ee6SHong Zhang #undef __FUNCT__
245277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials"
246277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts)
247277b19d0SLisandro Dalcin {
248277b19d0SLisandro Dalcin   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
249277b19d0SLisandro Dalcin   PetscErrorCode ierr;
250277b19d0SLisandro Dalcin 
251277b19d0SLisandro Dalcin   PetscFunctionBegin;
2525c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
2530679a0aeSJed Brown   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
2545c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
2555c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
256277b19d0SLisandro Dalcin   if (cvode->mem)    {CVodeFree(&cvode->mem);}
257277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
258277b19d0SLisandro Dalcin }
259277b19d0SLisandro Dalcin 
260277b19d0SLisandro Dalcin #undef __FUNCT__
2617cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials"
2627cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
2637cb94ee6SHong Zhang {
2647cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2657cb94ee6SHong Zhang   PetscErrorCode ierr;
2667cb94ee6SHong Zhang 
2677cb94ee6SHong Zhang   PetscFunctionBegin;
268277b19d0SLisandro Dalcin   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
269a07356b8SSatish Balay   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
270277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
271335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr);
272f61b2b6aSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C","",PETSC_NULL);CHKERRQ(ierr);
273335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
274335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr);
275335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
276335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
277335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
278335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr);
279335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr);
280335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr);
2817cb94ee6SHong Zhang   PetscFunctionReturn(0);
2827cb94ee6SHong Zhang }
2837cb94ee6SHong Zhang 
2847cb94ee6SHong Zhang #undef __FUNCT__
285214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials"
286214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts)
2877cb94ee6SHong Zhang {
2887cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2897cb94ee6SHong Zhang   PetscErrorCode ierr;
29016016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2917cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
29216016685SHong Zhang   void           *mem;
293dcbc6d53SSean Farley   PC             pc;
29416016685SHong Zhang   const PCType   pctype;
295ace3abfcSBarry Smith   PetscBool      pcnone;
2967cb94ee6SHong Zhang 
2977cb94ee6SHong Zhang   PetscFunctionBegin;
2987cb94ee6SHong Zhang   /* get the vector size */
2997cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
3007cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
3017cb94ee6SHong Zhang 
3027cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
303a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
304e32f2f54SBarry Smith   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
3057cb94ee6SHong Zhang 
30628adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
3077cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
3087cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
3097cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
3107cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
31142528757SHong Zhang 
3127cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
3130679a0aeSJed Brown   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
3147cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
3150679a0aeSJed Brown   ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr);
3167cb94ee6SHong Zhang 
3177cb94ee6SHong Zhang   /*
3187cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
3197cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
3207cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
3217cb94ee6SHong Zhang   */
3227adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
3237adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
3247cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
3257cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
32616016685SHong Zhang 
32716016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
32816016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
329e32f2f54SBarry Smith   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
33016016685SHong Zhang   cvode->mem = mem;
33116016685SHong Zhang 
33216016685SHong Zhang   /* Set the pointer to user-defined data */
33316016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
334e32f2f54SBarry Smith   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
33516016685SHong Zhang 
336fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
337*cdbf8f93SLisandro Dalcin   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
338fc6b9e64SJed Brown   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
339f1cd61daSJed Brown   if (cvode->mindt > 0) {
340f1cd61daSJed Brown     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
3419f94935aSHong Zhang     if (flag){
3429f94935aSHong Zhang       if (flag == CV_MEM_NULL){
3439f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
3449f94935aSHong Zhang       } else if (flag == CV_ILL_INPUT){
3459f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
3469f94935aSHong Zhang       } else {
3479f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
3489f94935aSHong Zhang       }
3499f94935aSHong Zhang     }
350f1cd61daSJed Brown   }
351f1cd61daSJed Brown   if (cvode->maxdt > 0) {
352f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
353f1cd61daSJed Brown     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
354f1cd61daSJed Brown   }
355f1cd61daSJed Brown 
35616016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
35716016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
35816016685SHong Zhang    * the initial dependent variable vector cvode->y */
35916016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
36016016685SHong Zhang   if (flag){
361e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
36216016685SHong Zhang   }
36316016685SHong Zhang 
3649f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
36516016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
36616016685SHong Zhang   if (flag){
367e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
36816016685SHong Zhang   }
36916016685SHong Zhang 
3709f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3719f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
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 {
382f61b2b6aSHong Zhang     flag  = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
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;
4097bda3a07SJed Brown   PC             pc;
4107cb94ee6SHong Zhang 
4117cb94ee6SHong Zhang   PetscFunctionBegin;
4127cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
4136fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
4147cb94ee6SHong Zhang     if (flag) {
4156fadb2cdSHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
4167cb94ee6SHong Zhang     }
4176fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
4187cb94ee6SHong Zhang     if (flag) {
4197cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
4207cb94ee6SHong Zhang     }
4217cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
4227cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
423f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
424f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
4257cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
426f61b2b6aSHong Zhang     ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);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);
4297bda3a07SJed Brown   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
4307bda3a07SJed Brown   ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
4317cb94ee6SHong Zhang   PetscFunctionReturn(0);
4327cb94ee6SHong Zhang }
4337cb94ee6SHong Zhang 
4347cb94ee6SHong Zhang #undef __FUNCT__
4357cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
4367cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
4377cb94ee6SHong Zhang {
4387cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4397cb94ee6SHong Zhang   PetscErrorCode ierr;
4407cb94ee6SHong Zhang   char           *type;
4417cb94ee6SHong Zhang   char           atype[] = "Adams";
4427cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
443ace3abfcSBarry Smith   PetscBool      iascii,isstring;
4442c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
4452c823083SHong Zhang   PetscInt       qlast,qcur;
4465d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
44742528757SHong Zhang   PC             pc;
4487cb94ee6SHong Zhang 
4497cb94ee6SHong Zhang   PetscFunctionBegin;
4507cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
4517cb94ee6SHong Zhang   else                                     {type = btype;}
4527cb94ee6SHong Zhang 
4532692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4542692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
4557cb94ee6SHong Zhang   if (iascii) {
4567cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
4577cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
4587cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
4597cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
460f61b2b6aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr);
4617cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4627cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4637cb94ee6SHong Zhang     } else {
4647cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4657cb94ee6SHong Zhang     }
466f1cd61daSJed Brown     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
467f1cd61daSJed Brown     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
4682c823083SHong Zhang 
4695d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4705d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
4715d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
4722c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
4732c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
4742c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
4752c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
4762c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
4772c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
4782c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
4792c823083SHong Zhang 
4802c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
4812c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
4822c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
4832c823083SHong Zhang 
4842c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
4852c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
4862c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
4872c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
48842528757SHong Zhang 
48942528757SHong Zhang     ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
49042528757SHong Zhang     ierr = PCView(pc,viewer);CHKERRQ(ierr);
4912c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4922c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
4932c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
4942c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
49542528757SHong Zhang 
4962c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4972c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
4982c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4992c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
5007cb94ee6SHong Zhang   } else if (isstring) {
5017cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
5027cb94ee6SHong Zhang   }
5037cb94ee6SHong Zhang   PetscFunctionReturn(0);
5047cb94ee6SHong Zhang }
5057cb94ee6SHong Zhang 
5067cb94ee6SHong Zhang 
5077cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
5087cb94ee6SHong Zhang EXTERN_C_BEGIN
5097cb94ee6SHong Zhang #undef __FUNCT__
5107cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
5117087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
5127cb94ee6SHong Zhang {
5137cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5147cb94ee6SHong Zhang 
5157cb94ee6SHong Zhang   PetscFunctionBegin;
5167cb94ee6SHong Zhang   cvode->cvode_type = type;
5177cb94ee6SHong Zhang   PetscFunctionReturn(0);
5187cb94ee6SHong Zhang }
5197cb94ee6SHong Zhang EXTERN_C_END
5207cb94ee6SHong Zhang 
5217cb94ee6SHong Zhang EXTERN_C_BEGIN
5227cb94ee6SHong Zhang #undef __FUNCT__
523f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl_Sundials"
524f61b2b6aSHong Zhang PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
5257cb94ee6SHong Zhang {
5267cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5277cb94ee6SHong Zhang 
5287cb94ee6SHong Zhang   PetscFunctionBegin;
529f61b2b6aSHong Zhang   cvode->maxl = maxl;
5307cb94ee6SHong Zhang   PetscFunctionReturn(0);
5317cb94ee6SHong Zhang }
5327cb94ee6SHong Zhang EXTERN_C_END
5337cb94ee6SHong Zhang 
5347cb94ee6SHong Zhang EXTERN_C_BEGIN
5357cb94ee6SHong Zhang #undef __FUNCT__
5367cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
5377087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
5387cb94ee6SHong Zhang {
5397cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5407cb94ee6SHong Zhang 
5417cb94ee6SHong Zhang   PetscFunctionBegin;
5427cb94ee6SHong Zhang   cvode->linear_tol = tol;
5437cb94ee6SHong Zhang   PetscFunctionReturn(0);
5447cb94ee6SHong Zhang }
5457cb94ee6SHong Zhang EXTERN_C_END
5467cb94ee6SHong Zhang 
5477cb94ee6SHong Zhang EXTERN_C_BEGIN
5487cb94ee6SHong Zhang #undef __FUNCT__
5497cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
5507087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
5517cb94ee6SHong Zhang {
5527cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5537cb94ee6SHong Zhang 
5547cb94ee6SHong Zhang   PetscFunctionBegin;
5557cb94ee6SHong Zhang   cvode->gtype = type;
5567cb94ee6SHong Zhang   PetscFunctionReturn(0);
5577cb94ee6SHong Zhang }
5587cb94ee6SHong Zhang EXTERN_C_END
5597cb94ee6SHong Zhang 
5607cb94ee6SHong Zhang EXTERN_C_BEGIN
5617cb94ee6SHong Zhang #undef __FUNCT__
5627cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
5637087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
5647cb94ee6SHong Zhang {
5657cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5667cb94ee6SHong Zhang 
5677cb94ee6SHong Zhang   PetscFunctionBegin;
5687cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5697cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
5707cb94ee6SHong Zhang   PetscFunctionReturn(0);
5717cb94ee6SHong Zhang }
5727cb94ee6SHong Zhang EXTERN_C_END
5737cb94ee6SHong Zhang 
5747cb94ee6SHong Zhang EXTERN_C_BEGIN
5757cb94ee6SHong Zhang #undef __FUNCT__
576f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
5777087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
578f1cd61daSJed Brown {
579f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
580f1cd61daSJed Brown 
581f1cd61daSJed Brown   PetscFunctionBegin;
582f1cd61daSJed Brown   cvode->mindt = mindt;
583f1cd61daSJed Brown   PetscFunctionReturn(0);
584f1cd61daSJed Brown }
585f1cd61daSJed Brown EXTERN_C_END
586f1cd61daSJed Brown 
587f1cd61daSJed Brown EXTERN_C_BEGIN
588f1cd61daSJed Brown #undef __FUNCT__
589f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
5907087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
591f1cd61daSJed Brown {
592f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
593f1cd61daSJed Brown 
594f1cd61daSJed Brown   PetscFunctionBegin;
595f1cd61daSJed Brown   cvode->maxdt = maxdt;
596f1cd61daSJed Brown   PetscFunctionReturn(0);
597f1cd61daSJed Brown }
598f1cd61daSJed Brown EXTERN_C_END
599f1cd61daSJed Brown 
600f1cd61daSJed Brown EXTERN_C_BEGIN
601f1cd61daSJed Brown #undef __FUNCT__
6027cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
6037087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
6047cb94ee6SHong Zhang {
605dcbc6d53SSean Farley   SNES            snes;
606dcbc6d53SSean Farley   KSP             ksp;
607dcbc6d53SSean Farley   PetscErrorCode  ierr;
6087cb94ee6SHong Zhang 
6097cb94ee6SHong Zhang   PetscFunctionBegin;
610dcbc6d53SSean Farley   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
611dcbc6d53SSean Farley   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
612dcbc6d53SSean Farley   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
6137cb94ee6SHong Zhang   PetscFunctionReturn(0);
6147cb94ee6SHong Zhang }
6157cb94ee6SHong Zhang EXTERN_C_END
6167cb94ee6SHong Zhang 
6177cb94ee6SHong Zhang EXTERN_C_BEGIN
6187cb94ee6SHong Zhang #undef __FUNCT__
6197cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
6207087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
6217cb94ee6SHong Zhang {
6227cb94ee6SHong Zhang   PetscFunctionBegin;
6232c823083SHong Zhang   if (nonlin) *nonlin = ts->nonlinear_its;
6242c823083SHong Zhang   if (lin)    *lin    = ts->linear_its;
6257cb94ee6SHong Zhang   PetscFunctionReturn(0);
6267cb94ee6SHong Zhang }
6277cb94ee6SHong Zhang EXTERN_C_END
6287cb94ee6SHong Zhang 
6297cb94ee6SHong Zhang EXTERN_C_BEGIN
6307cb94ee6SHong Zhang #undef __FUNCT__
6312bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
6327087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
6332bfc04deSHong Zhang {
6342bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6352bfc04deSHong Zhang 
6362bfc04deSHong Zhang   PetscFunctionBegin;
6372bfc04deSHong Zhang   cvode->monitorstep = s;
6382bfc04deSHong Zhang   PetscFunctionReturn(0);
6392bfc04deSHong Zhang }
6402bfc04deSHong Zhang EXTERN_C_END
6417cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
6427cb94ee6SHong Zhang 
6437cb94ee6SHong Zhang #undef __FUNCT__
6447cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
6457cb94ee6SHong Zhang /*@C
6467cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
6477cb94ee6SHong Zhang 
6487cb94ee6SHong Zhang    Not Collective
6497cb94ee6SHong Zhang 
6507cb94ee6SHong Zhang    Input parameters:
6517cb94ee6SHong Zhang .    ts     - the time-step context
6527cb94ee6SHong Zhang 
6537cb94ee6SHong Zhang    Output Parameters:
6547cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
6557cb94ee6SHong Zhang -   lin    - number of linear iterations
6567cb94ee6SHong Zhang 
6577cb94ee6SHong Zhang    Level: advanced
6587cb94ee6SHong Zhang 
6597cb94ee6SHong Zhang    Notes:
6607cb94ee6SHong Zhang     These return the number since the creation of the TS object
6617cb94ee6SHong Zhang 
6627cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
6637cb94ee6SHong Zhang 
664f61b2b6aSHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
6657cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
666f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
667a43b19c4SJed Brown           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
6687cb94ee6SHong Zhang 
6697cb94ee6SHong Zhang @*/
6707087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
6717cb94ee6SHong Zhang {
6724ac538c5SBarry Smith   PetscErrorCode ierr;
6737cb94ee6SHong Zhang 
6747cb94ee6SHong Zhang   PetscFunctionBegin;
6754ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
6767cb94ee6SHong Zhang   PetscFunctionReturn(0);
6777cb94ee6SHong Zhang }
6787cb94ee6SHong Zhang 
6797cb94ee6SHong Zhang #undef __FUNCT__
6807cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
6817cb94ee6SHong Zhang /*@
6827cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
6837cb94ee6SHong Zhang 
6843f9fe445SBarry Smith    Logically Collective on TS
6857cb94ee6SHong Zhang 
6867cb94ee6SHong Zhang    Input parameters:
6877cb94ee6SHong Zhang +    ts     - the time-step context
6887cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
6897cb94ee6SHong Zhang 
6907cb94ee6SHong Zhang    Level: intermediate
6917cb94ee6SHong Zhang 
6927cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
6937cb94ee6SHong Zhang 
694f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetMaxl(),
6957cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
696f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
6977cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
698a43b19c4SJed Brown           TSSetExactFinalTime()
6997cb94ee6SHong Zhang @*/
7007087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
7017cb94ee6SHong Zhang {
7024ac538c5SBarry Smith   PetscErrorCode ierr;
7037cb94ee6SHong Zhang 
7047cb94ee6SHong Zhang   PetscFunctionBegin;
7054ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
7067cb94ee6SHong Zhang   PetscFunctionReturn(0);
7077cb94ee6SHong Zhang }
7087cb94ee6SHong Zhang 
7097cb94ee6SHong Zhang #undef __FUNCT__
710f61b2b6aSHong Zhang #define __FUNCT__ "TSSundialsSetMaxl"
7117cb94ee6SHong Zhang /*@
712f61b2b6aSHong Zhang    TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
7137cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
714f61b2b6aSHong Zhang        this is the maximum number of GMRES steps that will be used.
7157cb94ee6SHong Zhang 
7163f9fe445SBarry Smith    Logically Collective on TS
7177cb94ee6SHong Zhang 
7187cb94ee6SHong Zhang    Input parameters:
7197cb94ee6SHong Zhang +    ts      - the time-step context
720f61b2b6aSHong Zhang -    maxl - number of direction vectors (the dimension of Krylov subspace).
7217cb94ee6SHong Zhang 
7227cb94ee6SHong Zhang    Level: advanced
7237cb94ee6SHong Zhang 
724f61b2b6aSHong Zhang .keywords: GMRES
7257cb94ee6SHong Zhang 
7267cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
7277cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
728f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7297cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
730a43b19c4SJed Brown           TSSetExactFinalTime()
7317cb94ee6SHong Zhang 
7327cb94ee6SHong Zhang @*/
733f61b2b6aSHong Zhang PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
7347cb94ee6SHong Zhang {
7354ac538c5SBarry Smith   PetscErrorCode ierr;
7367cb94ee6SHong Zhang 
7377cb94ee6SHong Zhang   PetscFunctionBegin;
738f61b2b6aSHong Zhang   PetscValidLogicalCollectiveInt(ts,maxl,2);
739f61b2b6aSHong Zhang   ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr);
7407cb94ee6SHong Zhang   PetscFunctionReturn(0);
7417cb94ee6SHong Zhang }
7427cb94ee6SHong Zhang 
7437cb94ee6SHong Zhang #undef __FUNCT__
7447cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
7457cb94ee6SHong Zhang /*@
7467cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
7477cb94ee6SHong Zhang        system by SUNDIALS.
7487cb94ee6SHong Zhang 
7493f9fe445SBarry Smith    Logically Collective on TS
7507cb94ee6SHong Zhang 
7517cb94ee6SHong Zhang    Input parameters:
7527cb94ee6SHong Zhang +    ts     - the time-step context
7537cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
7547cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
7557cb94ee6SHong Zhang 
7567cb94ee6SHong Zhang    Level: advanced
7577cb94ee6SHong Zhang 
7587cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
7597cb94ee6SHong Zhang 
760f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
7617cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
762f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7637cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
764a43b19c4SJed Brown           TSSetExactFinalTime()
7657cb94ee6SHong Zhang 
7667cb94ee6SHong Zhang @*/
7677087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
7687cb94ee6SHong Zhang {
7694ac538c5SBarry Smith   PetscErrorCode ierr;
7707cb94ee6SHong Zhang 
7717cb94ee6SHong Zhang   PetscFunctionBegin;
772c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,tol,2);
7734ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
7747cb94ee6SHong Zhang   PetscFunctionReturn(0);
7757cb94ee6SHong Zhang }
7767cb94ee6SHong Zhang 
7777cb94ee6SHong Zhang #undef __FUNCT__
7787cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
7797cb94ee6SHong Zhang /*@
7807cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
7817cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
7827cb94ee6SHong Zhang 
7833f9fe445SBarry Smith    Logically Collective on TS
7847cb94ee6SHong Zhang 
7857cb94ee6SHong Zhang    Input parameters:
7867cb94ee6SHong Zhang +    ts  - the time-step context
7877cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
7887cb94ee6SHong Zhang 
7897cb94ee6SHong Zhang    Level: advanced
7907cb94ee6SHong Zhang 
7917cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
7927cb94ee6SHong Zhang 
793f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
7947cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
795f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
7967cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
797a43b19c4SJed Brown           TSSetExactFinalTime()
7987cb94ee6SHong Zhang 
7997cb94ee6SHong Zhang @*/
8007087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
8017cb94ee6SHong Zhang {
8024ac538c5SBarry Smith   PetscErrorCode ierr;
8037cb94ee6SHong Zhang 
8047cb94ee6SHong Zhang   PetscFunctionBegin;
8054ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
8067cb94ee6SHong Zhang   PetscFunctionReturn(0);
8077cb94ee6SHong Zhang }
8087cb94ee6SHong Zhang 
8097cb94ee6SHong Zhang #undef __FUNCT__
8107cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
8117cb94ee6SHong Zhang /*@
8127cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
8137cb94ee6SHong Zhang                          Sundials for error control.
8147cb94ee6SHong Zhang 
8153f9fe445SBarry Smith    Logically Collective on TS
8167cb94ee6SHong Zhang 
8177cb94ee6SHong Zhang    Input parameters:
8187cb94ee6SHong Zhang +    ts  - the time-step context
8197cb94ee6SHong Zhang .    aabs - the absolute tolerance
8207cb94ee6SHong Zhang -    rel - the relative tolerance
8217cb94ee6SHong Zhang 
8227cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
8237cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
8247cb94ee6SHong Zhang 
8257cb94ee6SHong Zhang    Level: intermediate
8267cb94ee6SHong Zhang 
8277cb94ee6SHong Zhang .keywords: Sundials, tolerance
8287cb94ee6SHong Zhang 
829f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
8307cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
831f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
8327cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
833a43b19c4SJed Brown           TSSetExactFinalTime()
8347cb94ee6SHong Zhang 
8357cb94ee6SHong Zhang @*/
8367087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
8377cb94ee6SHong Zhang {
8384ac538c5SBarry Smith   PetscErrorCode ierr;
8397cb94ee6SHong Zhang 
8407cb94ee6SHong Zhang   PetscFunctionBegin;
8414ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
8427cb94ee6SHong Zhang   PetscFunctionReturn(0);
8437cb94ee6SHong Zhang }
8447cb94ee6SHong Zhang 
8457cb94ee6SHong Zhang #undef __FUNCT__
8467cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
8477cb94ee6SHong Zhang /*@
8487cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
8497cb94ee6SHong Zhang 
8507cb94ee6SHong Zhang    Input Parameter:
8517cb94ee6SHong Zhang .    ts - the time-step context
8527cb94ee6SHong Zhang 
8537cb94ee6SHong Zhang    Output Parameter:
8547cb94ee6SHong Zhang .    pc - the preconditioner context
8557cb94ee6SHong Zhang 
8567cb94ee6SHong Zhang    Level: advanced
8577cb94ee6SHong Zhang 
858f61b2b6aSHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
8597cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
860f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
8617cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
8627cb94ee6SHong Zhang @*/
8637087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
8647cb94ee6SHong Zhang {
8654ac538c5SBarry Smith   PetscErrorCode ierr;
8667cb94ee6SHong Zhang 
8677cb94ee6SHong Zhang   PetscFunctionBegin;
8684ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
8697cb94ee6SHong Zhang   PetscFunctionReturn(0);
8707cb94ee6SHong Zhang }
8717cb94ee6SHong Zhang 
8727cb94ee6SHong Zhang #undef __FUNCT__
873f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep"
874f1cd61daSJed Brown /*@
875f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
876f1cd61daSJed Brown 
877f1cd61daSJed Brown    Input Parameter:
878f1cd61daSJed Brown +   ts - the time-step context
879f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
880f1cd61daSJed Brown 
881fc6b9e64SJed Brown    Note:
882fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
883fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
884fc6b9e64SJed Brown 
885f1cd61daSJed Brown    Level: beginner
886f1cd61daSJed Brown 
887f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
888f1cd61daSJed Brown @*/
8897087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
890f1cd61daSJed Brown {
8914ac538c5SBarry Smith   PetscErrorCode ierr;
892f1cd61daSJed Brown 
893f1cd61daSJed Brown   PetscFunctionBegin;
8944ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
895f1cd61daSJed Brown   PetscFunctionReturn(0);
896f1cd61daSJed Brown }
897f1cd61daSJed Brown 
898f1cd61daSJed Brown #undef __FUNCT__
899f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep"
900f1cd61daSJed Brown /*@
901f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
902f1cd61daSJed Brown 
903f1cd61daSJed Brown    Input Parameter:
904f1cd61daSJed Brown +   ts - the time-step context
905f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
906f1cd61daSJed Brown 
907f1cd61daSJed Brown    Level: beginner
908f1cd61daSJed Brown 
909f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
910f1cd61daSJed Brown @*/
9117087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
912f1cd61daSJed Brown {
9134ac538c5SBarry Smith   PetscErrorCode ierr;
914f1cd61daSJed Brown 
915f1cd61daSJed Brown   PetscFunctionBegin;
9164ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
917f1cd61daSJed Brown   PetscFunctionReturn(0);
918f1cd61daSJed Brown }
919f1cd61daSJed Brown 
920f1cd61daSJed Brown #undef __FUNCT__
9212bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps"
9222bfc04deSHong Zhang /*@
9232bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
9242bfc04deSHong Zhang 
9252bfc04deSHong Zhang    Input Parameter:
9262bfc04deSHong Zhang +   ts - the time-step context
9272bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
9282bfc04deSHong Zhang 
9292bfc04deSHong Zhang    Level: beginner
9302bfc04deSHong Zhang 
931f61b2b6aSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
9322bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
933f61b2b6aSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(),
9342bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
9352bfc04deSHong Zhang @*/
9367087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
9372bfc04deSHong Zhang {
9384ac538c5SBarry Smith   PetscErrorCode ierr;
9392bfc04deSHong Zhang 
9402bfc04deSHong Zhang   PetscFunctionBegin;
9414ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
9422bfc04deSHong Zhang   PetscFunctionReturn(0);
9432bfc04deSHong Zhang }
9447cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
9457cb94ee6SHong Zhang /*MC
94696f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
9477cb94ee6SHong Zhang 
9487cb94ee6SHong Zhang    Options Database:
9497cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
9507cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
9517cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
9527cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
9537cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
954f61b2b6aSHong Zhang .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
95516016685SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
95616016685SHong Zhang 
9577cb94ee6SHong Zhang 
9587cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
9597cb94ee6SHong Zhang            only PETSc PC options
9607cb94ee6SHong Zhang 
9617cb94ee6SHong Zhang     Level: beginner
9627cb94ee6SHong Zhang 
963f61b2b6aSHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
964a43b19c4SJed Brown            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
9657cb94ee6SHong Zhang 
9667cb94ee6SHong Zhang M*/
9677cb94ee6SHong Zhang EXTERN_C_BEGIN
9687cb94ee6SHong Zhang #undef __FUNCT__
9697cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
9707087cfbeSBarry Smith PetscErrorCode  TSCreate_Sundials(TS ts)
9717cb94ee6SHong Zhang {
9727cb94ee6SHong Zhang   TS_Sundials    *cvode;
9737cb94ee6SHong Zhang   PetscErrorCode ierr;
97442528757SHong Zhang   PC             pc;
9757cb94ee6SHong Zhang 
9767cb94ee6SHong Zhang   PetscFunctionBegin;
977277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
97828adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
97928adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
980214bc6a2SJed Brown   ts->ops->setup          = TSSetUp_Sundials;
981b4eba00bSSean Farley   ts->ops->step           = TSStep_Sundials;
982b4eba00bSSean Farley   ts->ops->interpolate    = TSInterpolate_Sundials;
983214bc6a2SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
9847cb94ee6SHong Zhang 
98538f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
9867cb94ee6SHong Zhang   ts->data                = (void*)cvode;
9876fadb2cdSHong Zhang   cvode->cvode_type       = SUNDIALS_BDF;
9886fadb2cdSHong Zhang   cvode->gtype            = SUNDIALS_CLASSICAL_GS;
989f61b2b6aSHong Zhang   cvode->maxl             = 5;
9907cb94ee6SHong Zhang   cvode->linear_tol       = .05;
9917cb94ee6SHong Zhang 
992b4eba00bSSean Farley   cvode->monitorstep      = PETSC_TRUE;
9937cb94ee6SHong Zhang 
994a07356b8SSatish Balay   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
995f1cd61daSJed Brown 
996f1cd61daSJed Brown   cvode->mindt = -1.;
997f1cd61daSJed Brown   cvode->maxdt = -1.;
998f1cd61daSJed Brown 
9997cb94ee6SHong Zhang   /* set tolerance for Sundials */
10007cb94ee6SHong Zhang   cvode->reltol = 1e-6;
10012c823083SHong Zhang   cvode->abstol = 1e-6;
10027cb94ee6SHong Zhang 
100342528757SHong Zhang   /* set PCNONE as default pctype */
100442528757SHong Zhang   ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr);
100542528757SHong Zhang   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
100642528757SHong Zhang 
1007b4eba00bSSean Farley   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1008a43b19c4SJed Brown 
10097cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
10107cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
1011f61b2b6aSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C",
1012f61b2b6aSHong Zhang                     "TSSundialsSetMaxl_Sundials",
1013f61b2b6aSHong Zhang                     TSSundialsSetMaxl_Sundials);CHKERRQ(ierr);
10147cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
10157cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
10167cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
10177cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
10187cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
10197cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
10207cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
10217cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
10227cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1023f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1024f1cd61daSJed Brown                     "TSSundialsSetMinTimeStep_Sundials",
1025f1cd61daSJed Brown                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1026f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1027f1cd61daSJed Brown                     "TSSundialsSetMaxTimeStep_Sundials",
1028f1cd61daSJed Brown                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
10297cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
10307cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
10317cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
10327cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
10337cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
10347cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
10352bfc04deSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
10362bfc04deSHong Zhang                     "TSSundialsMonitorInternalSteps_Sundials",
10372bfc04deSHong Zhang                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
10382bfc04deSHong Zhang 
10397cb94ee6SHong Zhang   PetscFunctionReturn(0);
10407cb94ee6SHong Zhang }
10417cb94ee6SHong Zhang EXTERN_C_END
1042