xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision c75d03b0f056d483820f5ed49dcd95c2d19090d9)
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 /*
107117ce3f2SJed Brown        TSSolve_Sundials - Calls Sundials to integrate the ODE.
1087cb94ee6SHong Zhang */
1097cb94ee6SHong Zhang #undef __FUNCT__
110117ce3f2SJed Brown #define __FUNCT__ "TSSolve_Sundials"
111117ce3f2SJed Brown PetscErrorCode TSSolve_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;
116186e87acSLisandro Dalcin   PetscInt       i,flag;
1179f94935aSHong Zhang   long int       its,nsteps;
1187cb94ee6SHong Zhang   realtype       t,tout;
1197cb94ee6SHong Zhang   PetscScalar    *y_data;
1207cb94ee6SHong Zhang   void           *mem;
1217cb94ee6SHong Zhang 
1227cb94ee6SHong Zhang   PetscFunctionBegin;
12316016685SHong Zhang   mem  = cvode->mem;
1247cb94ee6SHong Zhang   tout = ts->max_time;
1257cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1267cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data,cvode->y);
1277cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
128186e87acSLisandro Dalcin 
129186e87acSLisandro Dalcin   for (i = 0; i < ts->max_steps; i++) {
1307cb94ee6SHong Zhang     if (ts->ptime >= ts->max_time) break;
1313f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
132186e87acSLisandro Dalcin 
1332bfc04deSHong Zhang     if (cvode->monitorstep) {
13428adb3f7SHong Zhang       flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
1352bfc04deSHong Zhang     } else {
1362bfc04deSHong Zhang       flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
1372bfc04deSHong Zhang     }
1389f94935aSHong Zhang 
1399f94935aSHong Zhang     if (flag){ /* display error message */
1409f94935aSHong Zhang       switch (flag){
1419f94935aSHong Zhang       case CV_ILL_INPUT:
1429f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
1439f94935aSHong Zhang         break;
1449f94935aSHong Zhang       case CV_TOO_CLOSE:
1459f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
1469f94935aSHong Zhang         break;
1473c7fefeeSJed Brown       case CV_TOO_MUCH_WORK: {
1489f94935aSHong Zhang         PetscReal      tcur;
1499f94935aSHong Zhang         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
1509f94935aSHong Zhang         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
1519f94935aSHong 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);
1523c7fefeeSJed Brown       } break;
1539f94935aSHong Zhang       case CV_TOO_MUCH_ACC:
1549f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
1559f94935aSHong Zhang         break;
1569f94935aSHong Zhang       case CV_ERR_FAILURE:
1579f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
1589f94935aSHong Zhang         break;
1599f94935aSHong Zhang       case CV_CONV_FAILURE:
1609f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
1619f94935aSHong Zhang         break;
1629f94935aSHong Zhang       case CV_LINIT_FAIL:
1639f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
1649f94935aSHong Zhang         break;
1659f94935aSHong Zhang       case CV_LSETUP_FAIL:
1669f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
1679f94935aSHong Zhang         break;
1689f94935aSHong Zhang       case CV_LSOLVE_FAIL:
1699f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
1709f94935aSHong Zhang         break;
1719f94935aSHong Zhang       case CV_RHSFUNC_FAIL:
1729f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
1739f94935aSHong Zhang         break;
1749f94935aSHong Zhang       case CV_FIRST_RHSFUNC_ERR:
1759f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
1769f94935aSHong Zhang         break;
1779f94935aSHong Zhang       case CV_REPTD_RHSFUNC_ERR:
1789f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
1799f94935aSHong Zhang         break;
1809f94935aSHong Zhang       case CV_UNREC_RHSFUNC_ERR:
1819f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
1829f94935aSHong Zhang         break;
1839f94935aSHong Zhang       case CV_RTFUNC_FAIL:
1849f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
1859f94935aSHong Zhang         break;
1869f94935aSHong Zhang       default:
1879f94935aSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
1889f94935aSHong Zhang       }
1899f94935aSHong Zhang     }
1909f94935aSHong Zhang 
1917cb94ee6SHong Zhang     if (t > ts->max_time && cvode->exact_final_time) {
1927cb94ee6SHong Zhang       /* interpolate to final requested time */
1937cb94ee6SHong Zhang       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
1947cb94ee6SHong Zhang       t = tout;
1957cb94ee6SHong Zhang     }
1967cb94ee6SHong Zhang 
1977cb94ee6SHong Zhang     /* copy the solution from cvode->y to cvode->update and sol */
1987cb94ee6SHong Zhang     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
1997cb94ee6SHong Zhang     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
2007cb94ee6SHong Zhang     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
2017cb94ee6SHong Zhang     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
2027cb94ee6SHong Zhang     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
2034ac7836bSHong Zhang     ierr = CVSpilsGetNumLinIters(mem, &its);
204186e87acSLisandro Dalcin     ts->nonlinear_its = its; ts->linear_its = its;
205186e87acSLisandro Dalcin 
206186e87acSLisandro Dalcin     ts->time_step = t - ts->ptime;
207186e87acSLisandro Dalcin     ts->ptime     = t;
2087cb94ee6SHong Zhang     ts->steps++;
209186e87acSLisandro Dalcin 
2103f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
2117cb94ee6SHong Zhang     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
2127cb94ee6SHong Zhang   }
2139f94935aSHong Zhang   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
214*c75d03b0SJed Brown   ts->steps = nsteps;
2157cb94ee6SHong Zhang   PetscFunctionReturn(0);
2167cb94ee6SHong Zhang }
2177cb94ee6SHong Zhang 
2187cb94ee6SHong Zhang #undef __FUNCT__
219277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials"
220277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts)
221277b19d0SLisandro Dalcin {
222277b19d0SLisandro Dalcin   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
223277b19d0SLisandro Dalcin   PetscErrorCode ierr;
224277b19d0SLisandro Dalcin 
225277b19d0SLisandro Dalcin   PetscFunctionBegin;
2265c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
2270679a0aeSJed Brown   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
2285c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
2295c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
230277b19d0SLisandro Dalcin   if (cvode->mem)    {CVodeFree(&cvode->mem);}
231277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
232277b19d0SLisandro Dalcin }
233277b19d0SLisandro Dalcin 
234277b19d0SLisandro Dalcin #undef __FUNCT__
2357cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials"
2367cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
2377cb94ee6SHong Zhang {
2387cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2397cb94ee6SHong Zhang   PetscErrorCode ierr;
2407cb94ee6SHong Zhang 
2417cb94ee6SHong Zhang   PetscFunctionBegin;
242277b19d0SLisandro Dalcin   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
243a07356b8SSatish Balay   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
244277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
245335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr);
246335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C","",PETSC_NULL);CHKERRQ(ierr);
247335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
248335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr);
249335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
250335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
251335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
252335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr);
253335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr);
254335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C","",PETSC_NULL);CHKERRQ(ierr);
255335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr);
2567cb94ee6SHong Zhang   PetscFunctionReturn(0);
2577cb94ee6SHong Zhang }
2587cb94ee6SHong Zhang 
2597cb94ee6SHong Zhang #undef __FUNCT__
260214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials"
261214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts)
2627cb94ee6SHong Zhang {
2637cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2647cb94ee6SHong Zhang   PetscErrorCode ierr;
26516016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2667cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
26716016685SHong Zhang   void           *mem;
268dcbc6d53SSean Farley   PC             pc;
26916016685SHong Zhang   const PCType   pctype;
270ace3abfcSBarry Smith   PetscBool      pcnone;
2717cb94ee6SHong Zhang 
2727cb94ee6SHong Zhang   PetscFunctionBegin;
2739f94935aSHong Zhang 
2747cb94ee6SHong Zhang   /* get the vector size */
2757cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2767cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
2777cb94ee6SHong Zhang 
2787cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
279a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
280e32f2f54SBarry Smith   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
2817cb94ee6SHong Zhang 
28228adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
2837cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
2847cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
2857cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
2867cb94ee6SHong Zhang   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
2877cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
2887cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
2890679a0aeSJed Brown   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
2907cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
2910679a0aeSJed Brown   ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr);
2927cb94ee6SHong Zhang 
2937cb94ee6SHong Zhang   /*
2947cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
2957cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
2967cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
2977cb94ee6SHong Zhang   */
2987adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
2997adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
3007cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
3017cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
30216016685SHong Zhang 
30316016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
30416016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
305e32f2f54SBarry Smith   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
30616016685SHong Zhang   cvode->mem = mem;
30716016685SHong Zhang 
30816016685SHong Zhang   /* Set the pointer to user-defined data */
30916016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
310e32f2f54SBarry Smith   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
31116016685SHong Zhang 
312fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
313fc6b9e64SJed Brown   flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step);
314fc6b9e64SJed Brown   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
315f1cd61daSJed Brown   if (cvode->mindt > 0) {
316f1cd61daSJed Brown     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
3179f94935aSHong Zhang     if (flag){
3189f94935aSHong Zhang       if (flag == CV_MEM_NULL){
3199f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
3209f94935aSHong Zhang       } else if (flag == CV_ILL_INPUT){
3219f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
3229f94935aSHong Zhang       } else {
3239f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
3249f94935aSHong Zhang       }
3259f94935aSHong Zhang     }
326f1cd61daSJed Brown   }
327f1cd61daSJed Brown   if (cvode->maxdt > 0) {
328f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
329f1cd61daSJed Brown     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
330f1cd61daSJed Brown   }
331f1cd61daSJed Brown 
33216016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
33316016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
33416016685SHong Zhang    * the initial dependent variable vector cvode->y */
33516016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
33616016685SHong Zhang   if (flag){
337e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
33816016685SHong Zhang   }
33916016685SHong Zhang 
3409f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
34116016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
34216016685SHong Zhang   if (flag){
343e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
34416016685SHong Zhang   }
34516016685SHong Zhang 
3469f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3479f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
34816016685SHong Zhang 
34916016685SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
35016016685SHong Zhang   /* setup the ode integrator with the given preconditioner */
351dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
352dcbc6d53SSean Farley   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
353dcbc6d53SSean Farley   ierr = PetscTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
35416016685SHong Zhang   if (pcnone){
35516016685SHong Zhang     flag  = CVSpgmr(mem,PREC_NONE,0);
356e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
35716016685SHong Zhang   } else {
35816016685SHong Zhang     flag  = CVSpgmr(mem,PREC_LEFT,0);
359e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
36016016685SHong Zhang 
36116016685SHong Zhang     /* Set preconditioner and solve routines Precond and PSolve,
36216016685SHong Zhang      and the pointer to the user-defined block data */
36316016685SHong Zhang     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
364e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
36516016685SHong Zhang   }
36616016685SHong Zhang 
36716016685SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
368e32f2f54SBarry Smith   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
3697cb94ee6SHong Zhang   PetscFunctionReturn(0);
3707cb94ee6SHong Zhang }
3717cb94ee6SHong Zhang 
3726fadb2cdSHong Zhang /* type of CVODE linear multistep method */
373dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
3746fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
375dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
376a04cf4d8SBarry Smith 
3777cb94ee6SHong Zhang #undef __FUNCT__
378214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials"
379214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts)
3807cb94ee6SHong Zhang {
3817cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3827cb94ee6SHong Zhang   PetscErrorCode ierr;
3837cb94ee6SHong Zhang   int            indx;
384ace3abfcSBarry Smith   PetscBool      flag;
3857cb94ee6SHong Zhang 
3867cb94ee6SHong Zhang   PetscFunctionBegin;
3877cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
3886fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
3897cb94ee6SHong Zhang     if (flag) {
3906fadb2cdSHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
3917cb94ee6SHong Zhang     }
3926fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
3937cb94ee6SHong Zhang     if (flag) {
3947cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
3957cb94ee6SHong Zhang     }
3967cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
3977cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
398f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
399f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
4007cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
4017cb94ee6SHong Zhang     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
402acfcf0e5SJed 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);
403acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
4047cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4057cb94ee6SHong Zhang   PetscFunctionReturn(0);
4067cb94ee6SHong Zhang }
4077cb94ee6SHong Zhang 
4087cb94ee6SHong Zhang #undef __FUNCT__
4097cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
4107cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
4117cb94ee6SHong Zhang {
4127cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4137cb94ee6SHong Zhang   PetscErrorCode ierr;
4147cb94ee6SHong Zhang   char           *type;
4157cb94ee6SHong Zhang   char           atype[] = "Adams";
4167cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
417ace3abfcSBarry Smith   PetscBool      iascii,isstring;
4182c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
4192c823083SHong Zhang   PetscInt       qlast,qcur;
4205d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
4217cb94ee6SHong Zhang 
4227cb94ee6SHong Zhang   PetscFunctionBegin;
4237cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
4247cb94ee6SHong Zhang   else                                     {type = btype;}
4257cb94ee6SHong Zhang 
4262692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4272692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
4287cb94ee6SHong Zhang   if (iascii) {
4297cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
4307cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
4317cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
4327cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
4337cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
4347cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4357cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4367cb94ee6SHong Zhang     } else {
4377cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4387cb94ee6SHong Zhang     }
439f1cd61daSJed Brown     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
440f1cd61daSJed Brown     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
4412c823083SHong Zhang 
4425d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4435d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
4445d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
4452c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
4462c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
4472c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
4482c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
4492c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
4502c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
4512c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
4522c823083SHong Zhang 
4532c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
4542c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
4552c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
4562c823083SHong Zhang 
4572c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
4582c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
4592c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
4602c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
4612c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4622c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
4632c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
4642c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
4652c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4662c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
4672c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4682c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
4697cb94ee6SHong Zhang   } else if (isstring) {
4707cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
4717cb94ee6SHong Zhang   } else {
472e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
4737cb94ee6SHong Zhang   }
4747cb94ee6SHong Zhang   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
4757cb94ee6SHong Zhang   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
4767cb94ee6SHong Zhang   PetscFunctionReturn(0);
4777cb94ee6SHong Zhang }
4787cb94ee6SHong Zhang 
4797cb94ee6SHong Zhang 
4807cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
4817cb94ee6SHong Zhang EXTERN_C_BEGIN
4827cb94ee6SHong Zhang #undef __FUNCT__
4837cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
4847087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
4857cb94ee6SHong Zhang {
4867cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4877cb94ee6SHong Zhang 
4887cb94ee6SHong Zhang   PetscFunctionBegin;
4897cb94ee6SHong Zhang   cvode->cvode_type = type;
4907cb94ee6SHong Zhang   PetscFunctionReturn(0);
4917cb94ee6SHong Zhang }
4927cb94ee6SHong Zhang EXTERN_C_END
4937cb94ee6SHong Zhang 
4947cb94ee6SHong Zhang EXTERN_C_BEGIN
4957cb94ee6SHong Zhang #undef __FUNCT__
4967cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
4977087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
4987cb94ee6SHong Zhang {
4997cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5007cb94ee6SHong Zhang 
5017cb94ee6SHong Zhang   PetscFunctionBegin;
5027cb94ee6SHong Zhang   cvode->restart = restart;
5037cb94ee6SHong Zhang   PetscFunctionReturn(0);
5047cb94ee6SHong Zhang }
5057cb94ee6SHong Zhang EXTERN_C_END
5067cb94ee6SHong Zhang 
5077cb94ee6SHong Zhang EXTERN_C_BEGIN
5087cb94ee6SHong Zhang #undef __FUNCT__
5097cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
5107087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
5117cb94ee6SHong Zhang {
5127cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5137cb94ee6SHong Zhang 
5147cb94ee6SHong Zhang   PetscFunctionBegin;
5157cb94ee6SHong Zhang   cvode->linear_tol = tol;
5167cb94ee6SHong Zhang   PetscFunctionReturn(0);
5177cb94ee6SHong Zhang }
5187cb94ee6SHong Zhang EXTERN_C_END
5197cb94ee6SHong Zhang 
5207cb94ee6SHong Zhang EXTERN_C_BEGIN
5217cb94ee6SHong Zhang #undef __FUNCT__
5227cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
5237087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
5247cb94ee6SHong Zhang {
5257cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5267cb94ee6SHong Zhang 
5277cb94ee6SHong Zhang   PetscFunctionBegin;
5287cb94ee6SHong Zhang   cvode->gtype = type;
5297cb94ee6SHong Zhang   PetscFunctionReturn(0);
5307cb94ee6SHong Zhang }
5317cb94ee6SHong Zhang EXTERN_C_END
5327cb94ee6SHong Zhang 
5337cb94ee6SHong Zhang EXTERN_C_BEGIN
5347cb94ee6SHong Zhang #undef __FUNCT__
5357cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
5367087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
5377cb94ee6SHong Zhang {
5387cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5397cb94ee6SHong Zhang 
5407cb94ee6SHong Zhang   PetscFunctionBegin;
5417cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5427cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
5437cb94ee6SHong Zhang   PetscFunctionReturn(0);
5447cb94ee6SHong Zhang }
5457cb94ee6SHong Zhang EXTERN_C_END
5467cb94ee6SHong Zhang 
5477cb94ee6SHong Zhang EXTERN_C_BEGIN
5487cb94ee6SHong Zhang #undef __FUNCT__
549f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
5507087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
551f1cd61daSJed Brown {
552f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
553f1cd61daSJed Brown 
554f1cd61daSJed Brown   PetscFunctionBegin;
555f1cd61daSJed Brown   cvode->mindt = mindt;
556f1cd61daSJed Brown   PetscFunctionReturn(0);
557f1cd61daSJed Brown }
558f1cd61daSJed Brown EXTERN_C_END
559f1cd61daSJed Brown 
560f1cd61daSJed Brown EXTERN_C_BEGIN
561f1cd61daSJed Brown #undef __FUNCT__
562f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
5637087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
564f1cd61daSJed Brown {
565f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
566f1cd61daSJed Brown 
567f1cd61daSJed Brown   PetscFunctionBegin;
568f1cd61daSJed Brown   cvode->maxdt = maxdt;
569f1cd61daSJed Brown   PetscFunctionReturn(0);
570f1cd61daSJed Brown }
571f1cd61daSJed Brown EXTERN_C_END
572f1cd61daSJed Brown 
573f1cd61daSJed Brown EXTERN_C_BEGIN
574f1cd61daSJed Brown #undef __FUNCT__
5757cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
5767087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
5777cb94ee6SHong Zhang {
578dcbc6d53SSean Farley   SNES            snes;
579dcbc6d53SSean Farley   KSP             ksp;
580dcbc6d53SSean Farley   PetscErrorCode  ierr;
5817cb94ee6SHong Zhang 
5827cb94ee6SHong Zhang   PetscFunctionBegin;
583dcbc6d53SSean Farley   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
584dcbc6d53SSean Farley   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
585dcbc6d53SSean Farley   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
5867cb94ee6SHong Zhang   PetscFunctionReturn(0);
5877cb94ee6SHong Zhang }
5887cb94ee6SHong Zhang EXTERN_C_END
5897cb94ee6SHong Zhang 
5907cb94ee6SHong Zhang EXTERN_C_BEGIN
5917cb94ee6SHong Zhang #undef __FUNCT__
5927cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
5937087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
5947cb94ee6SHong Zhang {
5957cb94ee6SHong Zhang   PetscFunctionBegin;
5962c823083SHong Zhang   if (nonlin) *nonlin = ts->nonlinear_its;
5972c823083SHong Zhang   if (lin)    *lin    = ts->linear_its;
5987cb94ee6SHong Zhang   PetscFunctionReturn(0);
5997cb94ee6SHong Zhang }
6007cb94ee6SHong Zhang EXTERN_C_END
6017cb94ee6SHong Zhang 
6027cb94ee6SHong Zhang EXTERN_C_BEGIN
6037cb94ee6SHong Zhang #undef __FUNCT__
6047cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
6057087cfbeSBarry Smith PetscErrorCode  TSSundialsSetExactFinalTime_Sundials(TS ts,PetscBool  s)
6067cb94ee6SHong Zhang {
6077cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6087cb94ee6SHong Zhang 
6097cb94ee6SHong Zhang   PetscFunctionBegin;
6107cb94ee6SHong Zhang   cvode->exact_final_time = s;
6117cb94ee6SHong Zhang   PetscFunctionReturn(0);
6127cb94ee6SHong Zhang }
6137cb94ee6SHong Zhang EXTERN_C_END
6142bfc04deSHong Zhang 
6152bfc04deSHong Zhang EXTERN_C_BEGIN
6162bfc04deSHong Zhang #undef __FUNCT__
6172bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
6187087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
6192bfc04deSHong Zhang {
6202bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6212bfc04deSHong Zhang 
6222bfc04deSHong Zhang   PetscFunctionBegin;
6232bfc04deSHong Zhang   cvode->monitorstep = s;
6242bfc04deSHong Zhang   PetscFunctionReturn(0);
6252bfc04deSHong Zhang }
6262bfc04deSHong Zhang EXTERN_C_END
6277cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
6287cb94ee6SHong Zhang 
6297cb94ee6SHong Zhang #undef __FUNCT__
6307cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
6317cb94ee6SHong Zhang /*@C
6327cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
6337cb94ee6SHong Zhang 
6347cb94ee6SHong Zhang    Not Collective
6357cb94ee6SHong Zhang 
6367cb94ee6SHong Zhang    Input parameters:
6377cb94ee6SHong Zhang .    ts     - the time-step context
6387cb94ee6SHong Zhang 
6397cb94ee6SHong Zhang    Output Parameters:
6407cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
6417cb94ee6SHong Zhang -   lin    - number of linear iterations
6427cb94ee6SHong Zhang 
6437cb94ee6SHong Zhang    Level: advanced
6447cb94ee6SHong Zhang 
6457cb94ee6SHong Zhang    Notes:
6467cb94ee6SHong Zhang     These return the number since the creation of the TS object
6477cb94ee6SHong Zhang 
6487cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
6497cb94ee6SHong Zhang 
6507cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6517cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6527cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
653f1cd61daSJed Brown           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSundialsSetExactFinalTime()
6547cb94ee6SHong Zhang 
6557cb94ee6SHong Zhang @*/
6567087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
6577cb94ee6SHong Zhang {
6584ac538c5SBarry Smith   PetscErrorCode ierr;
6597cb94ee6SHong Zhang 
6607cb94ee6SHong Zhang   PetscFunctionBegin;
6614ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
6627cb94ee6SHong Zhang   PetscFunctionReturn(0);
6637cb94ee6SHong Zhang }
6647cb94ee6SHong Zhang 
6657cb94ee6SHong Zhang #undef __FUNCT__
6667cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
6677cb94ee6SHong Zhang /*@
6687cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
6697cb94ee6SHong Zhang 
6703f9fe445SBarry Smith    Logically Collective on TS
6717cb94ee6SHong Zhang 
6727cb94ee6SHong Zhang    Input parameters:
6737cb94ee6SHong Zhang +    ts     - the time-step context
6747cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
6757cb94ee6SHong Zhang 
6767cb94ee6SHong Zhang    Level: intermediate
6777cb94ee6SHong Zhang 
6787cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
6797cb94ee6SHong Zhang 
6807cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
6817cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6827cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6837cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6847cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6857cb94ee6SHong Zhang @*/
6867087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
6877cb94ee6SHong Zhang {
6884ac538c5SBarry Smith   PetscErrorCode ierr;
6897cb94ee6SHong Zhang 
6907cb94ee6SHong Zhang   PetscFunctionBegin;
6914ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
6927cb94ee6SHong Zhang   PetscFunctionReturn(0);
6937cb94ee6SHong Zhang }
6947cb94ee6SHong Zhang 
6957cb94ee6SHong Zhang #undef __FUNCT__
6967cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart"
6977cb94ee6SHong Zhang /*@
6987cb94ee6SHong Zhang    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
6997cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
7007cb94ee6SHong Zhang        this is ALSO the maximum number of GMRES steps that will be used.
7017cb94ee6SHong Zhang 
7023f9fe445SBarry Smith    Logically Collective on TS
7037cb94ee6SHong Zhang 
7047cb94ee6SHong Zhang    Input parameters:
7057cb94ee6SHong Zhang +    ts      - the time-step context
7067cb94ee6SHong Zhang -    restart - number of direction vectors (the restart size).
7077cb94ee6SHong Zhang 
7087cb94ee6SHong Zhang    Level: advanced
7097cb94ee6SHong Zhang 
7107cb94ee6SHong Zhang .keywords: GMRES, restart
7117cb94ee6SHong Zhang 
7127cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
7137cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7147cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7157cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7167cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7177cb94ee6SHong Zhang 
7187cb94ee6SHong Zhang @*/
7197087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart(TS ts,int restart)
7207cb94ee6SHong Zhang {
7214ac538c5SBarry Smith   PetscErrorCode ierr;
7227cb94ee6SHong Zhang 
7237cb94ee6SHong Zhang   PetscFunctionBegin;
724c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,restart,2);
7254ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr);
7267cb94ee6SHong Zhang   PetscFunctionReturn(0);
7277cb94ee6SHong Zhang }
7287cb94ee6SHong Zhang 
7297cb94ee6SHong Zhang #undef __FUNCT__
7307cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
7317cb94ee6SHong Zhang /*@
7327cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
7337cb94ee6SHong Zhang        system by SUNDIALS.
7347cb94ee6SHong Zhang 
7353f9fe445SBarry Smith    Logically Collective on TS
7367cb94ee6SHong Zhang 
7377cb94ee6SHong Zhang    Input parameters:
7387cb94ee6SHong Zhang +    ts     - the time-step context
7397cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
7407cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
7417cb94ee6SHong Zhang 
7427cb94ee6SHong Zhang    Level: advanced
7437cb94ee6SHong Zhang 
7447cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
7457cb94ee6SHong Zhang 
7467cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7477cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7487cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7497cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7507cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7517cb94ee6SHong Zhang 
7527cb94ee6SHong Zhang @*/
7537087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
7547cb94ee6SHong Zhang {
7554ac538c5SBarry Smith   PetscErrorCode ierr;
7567cb94ee6SHong Zhang 
7577cb94ee6SHong Zhang   PetscFunctionBegin;
758c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,tol,2);
7594ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
7607cb94ee6SHong Zhang   PetscFunctionReturn(0);
7617cb94ee6SHong Zhang }
7627cb94ee6SHong Zhang 
7637cb94ee6SHong Zhang #undef __FUNCT__
7647cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
7657cb94ee6SHong Zhang /*@
7667cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
7677cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
7687cb94ee6SHong Zhang 
7693f9fe445SBarry Smith    Logically Collective on TS
7707cb94ee6SHong Zhang 
7717cb94ee6SHong Zhang    Input parameters:
7727cb94ee6SHong Zhang +    ts  - the time-step context
7737cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
7747cb94ee6SHong Zhang 
7757cb94ee6SHong Zhang    Level: advanced
7767cb94ee6SHong Zhang 
7777cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
7787cb94ee6SHong Zhang 
7797cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7807cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
7817cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7827cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7837cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7847cb94ee6SHong Zhang 
7857cb94ee6SHong Zhang @*/
7867087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
7877cb94ee6SHong Zhang {
7884ac538c5SBarry Smith   PetscErrorCode ierr;
7897cb94ee6SHong Zhang 
7907cb94ee6SHong Zhang   PetscFunctionBegin;
7914ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
7927cb94ee6SHong Zhang   PetscFunctionReturn(0);
7937cb94ee6SHong Zhang }
7947cb94ee6SHong Zhang 
7957cb94ee6SHong Zhang #undef __FUNCT__
7967cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
7977cb94ee6SHong Zhang /*@
7987cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
7997cb94ee6SHong Zhang                          Sundials for error control.
8007cb94ee6SHong Zhang 
8013f9fe445SBarry Smith    Logically Collective on TS
8027cb94ee6SHong Zhang 
8037cb94ee6SHong Zhang    Input parameters:
8047cb94ee6SHong Zhang +    ts  - the time-step context
8057cb94ee6SHong Zhang .    aabs - the absolute tolerance
8067cb94ee6SHong Zhang -    rel - the relative tolerance
8077cb94ee6SHong Zhang 
8087cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
8097cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
8107cb94ee6SHong Zhang 
8117cb94ee6SHong Zhang    Level: intermediate
8127cb94ee6SHong Zhang 
8137cb94ee6SHong Zhang .keywords: Sundials, tolerance
8147cb94ee6SHong Zhang 
8157cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8167cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
8177cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8187cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
8197cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
8207cb94ee6SHong Zhang 
8217cb94ee6SHong Zhang @*/
8227087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
8237cb94ee6SHong Zhang {
8244ac538c5SBarry Smith   PetscErrorCode ierr;
8257cb94ee6SHong Zhang 
8267cb94ee6SHong Zhang   PetscFunctionBegin;
8274ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
8287cb94ee6SHong Zhang   PetscFunctionReturn(0);
8297cb94ee6SHong Zhang }
8307cb94ee6SHong Zhang 
8317cb94ee6SHong Zhang #undef __FUNCT__
8327cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
8337cb94ee6SHong Zhang /*@
8347cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
8357cb94ee6SHong Zhang 
8367cb94ee6SHong Zhang    Input Parameter:
8377cb94ee6SHong Zhang .    ts - the time-step context
8387cb94ee6SHong Zhang 
8397cb94ee6SHong Zhang    Output Parameter:
8407cb94ee6SHong Zhang .    pc - the preconditioner context
8417cb94ee6SHong Zhang 
8427cb94ee6SHong Zhang    Level: advanced
8437cb94ee6SHong Zhang 
8447cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8457cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8467cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8477cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
8487cb94ee6SHong Zhang @*/
8497087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
8507cb94ee6SHong Zhang {
8514ac538c5SBarry Smith   PetscErrorCode ierr;
8527cb94ee6SHong Zhang 
8537cb94ee6SHong Zhang   PetscFunctionBegin;
8544ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
8557cb94ee6SHong Zhang   PetscFunctionReturn(0);
8567cb94ee6SHong Zhang }
8577cb94ee6SHong Zhang 
8587cb94ee6SHong Zhang #undef __FUNCT__
8597cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime"
8607cb94ee6SHong Zhang /*@
8617cb94ee6SHong Zhang    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
8627cb94ee6SHong Zhang       exact final time requested by the user or just returns it at the final time
8637cb94ee6SHong Zhang       it computed. (Defaults to true).
8647cb94ee6SHong Zhang 
8657cb94ee6SHong Zhang    Input Parameter:
8667cb94ee6SHong Zhang +   ts - the time-step context
8677cb94ee6SHong Zhang -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
8687cb94ee6SHong Zhang 
8697cb94ee6SHong Zhang    Level: beginner
8707cb94ee6SHong Zhang 
8717cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8727cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8737cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8747cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
8757cb94ee6SHong Zhang @*/
8767087cfbeSBarry Smith PetscErrorCode  TSSundialsSetExactFinalTime(TS ts,PetscBool  ft)
8777cb94ee6SHong Zhang {
8784ac538c5SBarry Smith   PetscErrorCode ierr;
8797cb94ee6SHong Zhang 
8807cb94ee6SHong Zhang   PetscFunctionBegin;
8814ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetExactFinalTime_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
8827cb94ee6SHong Zhang   PetscFunctionReturn(0);
8837cb94ee6SHong Zhang }
8847cb94ee6SHong Zhang 
8852bfc04deSHong Zhang #undef __FUNCT__
886f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep"
887f1cd61daSJed Brown /*@
888f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
889f1cd61daSJed Brown 
890f1cd61daSJed Brown    Input Parameter:
891f1cd61daSJed Brown +   ts - the time-step context
892f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
893f1cd61daSJed Brown 
894fc6b9e64SJed Brown    Note:
895fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
896fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
897fc6b9e64SJed Brown 
898f1cd61daSJed Brown    Level: beginner
899f1cd61daSJed Brown 
900f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
901f1cd61daSJed Brown @*/
9027087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
903f1cd61daSJed Brown {
9044ac538c5SBarry Smith   PetscErrorCode ierr;
905f1cd61daSJed Brown 
906f1cd61daSJed Brown   PetscFunctionBegin;
9074ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
908f1cd61daSJed Brown   PetscFunctionReturn(0);
909f1cd61daSJed Brown }
910f1cd61daSJed Brown 
911f1cd61daSJed Brown #undef __FUNCT__
912f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep"
913f1cd61daSJed Brown /*@
914f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
915f1cd61daSJed Brown 
916f1cd61daSJed Brown    Input Parameter:
917f1cd61daSJed Brown +   ts - the time-step context
918f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
919f1cd61daSJed Brown 
920f1cd61daSJed Brown    Level: beginner
921f1cd61daSJed Brown 
922f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
923f1cd61daSJed Brown @*/
9247087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
925f1cd61daSJed Brown {
9264ac538c5SBarry Smith   PetscErrorCode ierr;
927f1cd61daSJed Brown 
928f1cd61daSJed Brown   PetscFunctionBegin;
9294ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
930f1cd61daSJed Brown   PetscFunctionReturn(0);
931f1cd61daSJed Brown }
932f1cd61daSJed Brown 
933f1cd61daSJed Brown #undef __FUNCT__
9342bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps"
9352bfc04deSHong Zhang /*@
9362bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
9372bfc04deSHong Zhang 
9382bfc04deSHong Zhang    Input Parameter:
9392bfc04deSHong Zhang +   ts - the time-step context
9402bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
9412bfc04deSHong Zhang 
9422bfc04deSHong Zhang    Level: beginner
9432bfc04deSHong Zhang 
9442bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9452bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
9462bfc04deSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9472bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
9482bfc04deSHong Zhang @*/
9497087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
9502bfc04deSHong Zhang {
9514ac538c5SBarry Smith   PetscErrorCode ierr;
9522bfc04deSHong Zhang 
9532bfc04deSHong Zhang   PetscFunctionBegin;
9544ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
9552bfc04deSHong Zhang   PetscFunctionReturn(0);
9562bfc04deSHong Zhang }
9577cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
9587cb94ee6SHong Zhang /*MC
95996f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
9607cb94ee6SHong Zhang 
9617cb94ee6SHong Zhang    Options Database:
9627cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
9637cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
9647cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
9657cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
9667cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
9677cb94ee6SHong Zhang .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
96816016685SHong Zhang .    -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time
96916016685SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
97016016685SHong Zhang 
9717cb94ee6SHong Zhang 
9727cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
9737cb94ee6SHong Zhang            only PETSc PC options
9747cb94ee6SHong Zhang 
9757cb94ee6SHong Zhang     Level: beginner
9767cb94ee6SHong Zhang 
9777cb94ee6SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
9787cb94ee6SHong Zhang            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
9797cb94ee6SHong Zhang 
9807cb94ee6SHong Zhang M*/
9817cb94ee6SHong Zhang EXTERN_C_BEGIN
9827cb94ee6SHong Zhang #undef __FUNCT__
9837cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
9847087cfbeSBarry Smith PetscErrorCode  TSCreate_Sundials(TS ts)
9857cb94ee6SHong Zhang {
9867cb94ee6SHong Zhang   TS_Sundials *cvode;
9877cb94ee6SHong Zhang   PetscErrorCode ierr;
9887cb94ee6SHong Zhang 
9897cb94ee6SHong Zhang   PetscFunctionBegin;
990277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
99128adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
99228adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
993214bc6a2SJed Brown   ts->ops->setup          = TSSetUp_Sundials;
994117ce3f2SJed Brown   ts->ops->solve          = TSSolve_Sundials;
995214bc6a2SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
9967cb94ee6SHong Zhang 
99738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
9987cb94ee6SHong Zhang   ts->data          = (void*)cvode;
9996fadb2cdSHong Zhang   cvode->cvode_type = SUNDIALS_BDF;
10006fadb2cdSHong Zhang   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
10017cb94ee6SHong Zhang   cvode->restart    = 5;
10027cb94ee6SHong Zhang   cvode->linear_tol = .05;
10037cb94ee6SHong Zhang 
10042bfc04deSHong Zhang   cvode->exact_final_time = PETSC_TRUE;
10052bfc04deSHong Zhang   cvode->monitorstep      = PETSC_FALSE;
10067cb94ee6SHong Zhang 
1007a07356b8SSatish Balay   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
1008f1cd61daSJed Brown 
1009f1cd61daSJed Brown   cvode->mindt = -1.;
1010f1cd61daSJed Brown   cvode->maxdt = -1.;
1011f1cd61daSJed Brown 
10127cb94ee6SHong Zhang   /* set tolerance for Sundials */
10137cb94ee6SHong Zhang   cvode->reltol = 1e-6;
10142c823083SHong Zhang   cvode->abstol = 1e-6;
10157cb94ee6SHong Zhang 
10167cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
10177cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
10187cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
10197cb94ee6SHong Zhang                     "TSSundialsSetGMRESRestart_Sundials",
10207cb94ee6SHong Zhang                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
10217cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
10227cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
10237cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
10247cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
10257cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
10267cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
10277cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
10287cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
10297cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1030f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1031f1cd61daSJed Brown                     "TSSundialsSetMinTimeStep_Sundials",
1032f1cd61daSJed Brown                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1033f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1034f1cd61daSJed Brown                     "TSSundialsSetMaxTimeStep_Sundials",
1035f1cd61daSJed Brown                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
10367cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
10377cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
10387cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
10397cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
10407cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
10417cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
10427cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
10437cb94ee6SHong Zhang                     "TSSundialsSetExactFinalTime_Sundials",
10447cb94ee6SHong Zhang                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
10452bfc04deSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
10462bfc04deSHong Zhang                     "TSSundialsMonitorInternalSteps_Sundials",
10472bfc04deSHong Zhang                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
10482bfc04deSHong Zhang 
10497cb94ee6SHong Zhang   PetscFunctionReturn(0);
10507cb94ee6SHong Zhang }
10517cb94ee6SHong Zhang EXTERN_C_END
1052