xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision a43b19c4d848258ea28797885a7bfde089a29a6d)
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);
214c75d03b0SJed 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,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr);
2557cb94ee6SHong Zhang   PetscFunctionReturn(0);
2567cb94ee6SHong Zhang }
2577cb94ee6SHong Zhang 
2587cb94ee6SHong Zhang #undef __FUNCT__
259214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials"
260214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts)
2617cb94ee6SHong Zhang {
2627cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2637cb94ee6SHong Zhang   PetscErrorCode ierr;
26416016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2657cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
26616016685SHong Zhang   void           *mem;
267dcbc6d53SSean Farley   PC             pc;
26816016685SHong Zhang   const PCType   pctype;
269ace3abfcSBarry Smith   PetscBool      pcnone;
2707cb94ee6SHong Zhang 
2717cb94ee6SHong Zhang   PetscFunctionBegin;
2729f94935aSHong Zhang 
2737cb94ee6SHong Zhang   /* get the vector size */
2747cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2757cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
2767cb94ee6SHong Zhang 
2777cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
278a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
279e32f2f54SBarry Smith   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
2807cb94ee6SHong Zhang 
28128adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
2827cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
2837cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
2847cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
2857cb94ee6SHong Zhang   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
2867cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
2877cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
2880679a0aeSJed Brown   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
2897cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
2900679a0aeSJed Brown   ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr);
2917cb94ee6SHong Zhang 
2927cb94ee6SHong Zhang   /*
2937cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
2947cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
2957cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
2967cb94ee6SHong Zhang   */
2977adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
2987adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
2997cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
3007cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
30116016685SHong Zhang 
30216016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
30316016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
304e32f2f54SBarry Smith   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
30516016685SHong Zhang   cvode->mem = mem;
30616016685SHong Zhang 
30716016685SHong Zhang   /* Set the pointer to user-defined data */
30816016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
309e32f2f54SBarry Smith   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
31016016685SHong Zhang 
311fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
312fc6b9e64SJed Brown   flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step);
313fc6b9e64SJed Brown   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
314f1cd61daSJed Brown   if (cvode->mindt > 0) {
315f1cd61daSJed Brown     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
3169f94935aSHong Zhang     if (flag){
3179f94935aSHong Zhang       if (flag == CV_MEM_NULL){
3189f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
3199f94935aSHong Zhang       } else if (flag == CV_ILL_INPUT){
3209f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
3219f94935aSHong Zhang       } else {
3229f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
3239f94935aSHong Zhang       }
3249f94935aSHong Zhang     }
325f1cd61daSJed Brown   }
326f1cd61daSJed Brown   if (cvode->maxdt > 0) {
327f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
328f1cd61daSJed Brown     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
329f1cd61daSJed Brown   }
330f1cd61daSJed Brown 
33116016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
33216016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
33316016685SHong Zhang    * the initial dependent variable vector cvode->y */
33416016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
33516016685SHong Zhang   if (flag){
336e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
33716016685SHong Zhang   }
33816016685SHong Zhang 
3399f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
34016016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
34116016685SHong Zhang   if (flag){
342e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
34316016685SHong Zhang   }
34416016685SHong Zhang 
3459f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3469f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
34716016685SHong Zhang 
34816016685SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
34916016685SHong Zhang   /* setup the ode integrator with the given preconditioner */
350dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
351dcbc6d53SSean Farley   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
352dcbc6d53SSean Farley   ierr = PetscTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
35316016685SHong Zhang   if (pcnone){
35416016685SHong Zhang     flag  = CVSpgmr(mem,PREC_NONE,0);
355e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
35616016685SHong Zhang   } else {
35716016685SHong Zhang     flag  = CVSpgmr(mem,PREC_LEFT,0);
358e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
35916016685SHong Zhang 
36016016685SHong Zhang     /* Set preconditioner and solve routines Precond and PSolve,
36116016685SHong Zhang      and the pointer to the user-defined block data */
36216016685SHong Zhang     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
363e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
36416016685SHong Zhang   }
36516016685SHong Zhang 
36616016685SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
367e32f2f54SBarry Smith   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
3687cb94ee6SHong Zhang   PetscFunctionReturn(0);
3697cb94ee6SHong Zhang }
3707cb94ee6SHong Zhang 
3716fadb2cdSHong Zhang /* type of CVODE linear multistep method */
372dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
3736fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
374dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
375a04cf4d8SBarry Smith 
3767cb94ee6SHong Zhang #undef __FUNCT__
377214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials"
378214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts)
3797cb94ee6SHong Zhang {
3807cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3817cb94ee6SHong Zhang   PetscErrorCode ierr;
3827cb94ee6SHong Zhang   int            indx;
383ace3abfcSBarry Smith   PetscBool      flag;
3847cb94ee6SHong Zhang 
3857cb94ee6SHong Zhang   PetscFunctionBegin;
3867cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
3876fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
3887cb94ee6SHong Zhang     if (flag) {
3896fadb2cdSHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
3907cb94ee6SHong Zhang     }
3916fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
3927cb94ee6SHong Zhang     if (flag) {
3937cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
3947cb94ee6SHong Zhang     }
3957cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
3967cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
397f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
398f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
3997cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
4007cb94ee6SHong Zhang     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
401acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
4027cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4037cb94ee6SHong Zhang   PetscFunctionReturn(0);
4047cb94ee6SHong Zhang }
4057cb94ee6SHong Zhang 
4067cb94ee6SHong Zhang #undef __FUNCT__
4077cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
4087cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
4097cb94ee6SHong Zhang {
4107cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4117cb94ee6SHong Zhang   PetscErrorCode ierr;
4127cb94ee6SHong Zhang   char           *type;
4137cb94ee6SHong Zhang   char           atype[] = "Adams";
4147cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
415ace3abfcSBarry Smith   PetscBool      iascii,isstring;
4162c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
4172c823083SHong Zhang   PetscInt       qlast,qcur;
4185d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
4197cb94ee6SHong Zhang 
4207cb94ee6SHong Zhang   PetscFunctionBegin;
4217cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
4227cb94ee6SHong Zhang   else                                     {type = btype;}
4237cb94ee6SHong Zhang 
4242692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4252692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
4267cb94ee6SHong Zhang   if (iascii) {
4277cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
4287cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
4297cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
4307cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
4317cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
4327cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4337cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4347cb94ee6SHong Zhang     } else {
4357cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4367cb94ee6SHong Zhang     }
437f1cd61daSJed Brown     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
438f1cd61daSJed Brown     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
4392c823083SHong Zhang 
4405d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4415d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
4425d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
4432c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
4442c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
4452c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
4462c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
4472c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
4482c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
4492c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
4502c823083SHong Zhang 
4512c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
4522c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
4532c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
4542c823083SHong Zhang 
4552c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
4562c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
4572c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
4582c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
4592c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4602c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
4612c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
4622c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
4632c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4642c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
4652c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4662c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
4677cb94ee6SHong Zhang   } else if (isstring) {
4687cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
4697cb94ee6SHong Zhang   }
4707cb94ee6SHong Zhang   PetscFunctionReturn(0);
4717cb94ee6SHong Zhang }
4727cb94ee6SHong Zhang 
4737cb94ee6SHong Zhang 
4747cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
4757cb94ee6SHong Zhang EXTERN_C_BEGIN
4767cb94ee6SHong Zhang #undef __FUNCT__
4777cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
4787087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
4797cb94ee6SHong Zhang {
4807cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4817cb94ee6SHong Zhang 
4827cb94ee6SHong Zhang   PetscFunctionBegin;
4837cb94ee6SHong Zhang   cvode->cvode_type = type;
4847cb94ee6SHong Zhang   PetscFunctionReturn(0);
4857cb94ee6SHong Zhang }
4867cb94ee6SHong Zhang EXTERN_C_END
4877cb94ee6SHong Zhang 
4887cb94ee6SHong Zhang EXTERN_C_BEGIN
4897cb94ee6SHong Zhang #undef __FUNCT__
4907cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
4917087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
4927cb94ee6SHong Zhang {
4937cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4947cb94ee6SHong Zhang 
4957cb94ee6SHong Zhang   PetscFunctionBegin;
4967cb94ee6SHong Zhang   cvode->restart = restart;
4977cb94ee6SHong Zhang   PetscFunctionReturn(0);
4987cb94ee6SHong Zhang }
4997cb94ee6SHong Zhang EXTERN_C_END
5007cb94ee6SHong Zhang 
5017cb94ee6SHong Zhang EXTERN_C_BEGIN
5027cb94ee6SHong Zhang #undef __FUNCT__
5037cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
5047087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
5057cb94ee6SHong Zhang {
5067cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5077cb94ee6SHong Zhang 
5087cb94ee6SHong Zhang   PetscFunctionBegin;
5097cb94ee6SHong Zhang   cvode->linear_tol = tol;
5107cb94ee6SHong Zhang   PetscFunctionReturn(0);
5117cb94ee6SHong Zhang }
5127cb94ee6SHong Zhang EXTERN_C_END
5137cb94ee6SHong Zhang 
5147cb94ee6SHong Zhang EXTERN_C_BEGIN
5157cb94ee6SHong Zhang #undef __FUNCT__
5167cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
5177087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
5187cb94ee6SHong Zhang {
5197cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5207cb94ee6SHong Zhang 
5217cb94ee6SHong Zhang   PetscFunctionBegin;
5227cb94ee6SHong Zhang   cvode->gtype = type;
5237cb94ee6SHong Zhang   PetscFunctionReturn(0);
5247cb94ee6SHong Zhang }
5257cb94ee6SHong Zhang EXTERN_C_END
5267cb94ee6SHong Zhang 
5277cb94ee6SHong Zhang EXTERN_C_BEGIN
5287cb94ee6SHong Zhang #undef __FUNCT__
5297cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
5307087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
5317cb94ee6SHong Zhang {
5327cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5337cb94ee6SHong Zhang 
5347cb94ee6SHong Zhang   PetscFunctionBegin;
5357cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5367cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
5377cb94ee6SHong Zhang   PetscFunctionReturn(0);
5387cb94ee6SHong Zhang }
5397cb94ee6SHong Zhang EXTERN_C_END
5407cb94ee6SHong Zhang 
5417cb94ee6SHong Zhang EXTERN_C_BEGIN
5427cb94ee6SHong Zhang #undef __FUNCT__
543f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
5447087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
545f1cd61daSJed Brown {
546f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
547f1cd61daSJed Brown 
548f1cd61daSJed Brown   PetscFunctionBegin;
549f1cd61daSJed Brown   cvode->mindt = mindt;
550f1cd61daSJed Brown   PetscFunctionReturn(0);
551f1cd61daSJed Brown }
552f1cd61daSJed Brown EXTERN_C_END
553f1cd61daSJed Brown 
554f1cd61daSJed Brown EXTERN_C_BEGIN
555f1cd61daSJed Brown #undef __FUNCT__
556f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
5577087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
558f1cd61daSJed Brown {
559f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
560f1cd61daSJed Brown 
561f1cd61daSJed Brown   PetscFunctionBegin;
562f1cd61daSJed Brown   cvode->maxdt = maxdt;
563f1cd61daSJed Brown   PetscFunctionReturn(0);
564f1cd61daSJed Brown }
565f1cd61daSJed Brown EXTERN_C_END
566f1cd61daSJed Brown 
567f1cd61daSJed Brown EXTERN_C_BEGIN
568f1cd61daSJed Brown #undef __FUNCT__
5697cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
5707087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
5717cb94ee6SHong Zhang {
572dcbc6d53SSean Farley   SNES            snes;
573dcbc6d53SSean Farley   KSP             ksp;
574dcbc6d53SSean Farley   PetscErrorCode  ierr;
5757cb94ee6SHong Zhang 
5767cb94ee6SHong Zhang   PetscFunctionBegin;
577dcbc6d53SSean Farley   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
578dcbc6d53SSean Farley   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
579dcbc6d53SSean Farley   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
5807cb94ee6SHong Zhang   PetscFunctionReturn(0);
5817cb94ee6SHong Zhang }
5827cb94ee6SHong Zhang EXTERN_C_END
5837cb94ee6SHong Zhang 
5847cb94ee6SHong Zhang EXTERN_C_BEGIN
5857cb94ee6SHong Zhang #undef __FUNCT__
5867cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
5877087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
5887cb94ee6SHong Zhang {
5897cb94ee6SHong Zhang   PetscFunctionBegin;
5902c823083SHong Zhang   if (nonlin) *nonlin = ts->nonlinear_its;
5912c823083SHong Zhang   if (lin)    *lin    = ts->linear_its;
5927cb94ee6SHong Zhang   PetscFunctionReturn(0);
5937cb94ee6SHong Zhang }
5947cb94ee6SHong Zhang EXTERN_C_END
5957cb94ee6SHong Zhang 
5967cb94ee6SHong Zhang EXTERN_C_BEGIN
5977cb94ee6SHong Zhang #undef __FUNCT__
5982bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
5997087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
6002bfc04deSHong Zhang {
6012bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6022bfc04deSHong Zhang 
6032bfc04deSHong Zhang   PetscFunctionBegin;
6042bfc04deSHong Zhang   cvode->monitorstep = s;
6052bfc04deSHong Zhang   PetscFunctionReturn(0);
6062bfc04deSHong Zhang }
6072bfc04deSHong Zhang EXTERN_C_END
6087cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
6097cb94ee6SHong Zhang 
6107cb94ee6SHong Zhang #undef __FUNCT__
6117cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
6127cb94ee6SHong Zhang /*@C
6137cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
6147cb94ee6SHong Zhang 
6157cb94ee6SHong Zhang    Not Collective
6167cb94ee6SHong Zhang 
6177cb94ee6SHong Zhang    Input parameters:
6187cb94ee6SHong Zhang .    ts     - the time-step context
6197cb94ee6SHong Zhang 
6207cb94ee6SHong Zhang    Output Parameters:
6217cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
6227cb94ee6SHong Zhang -   lin    - number of linear iterations
6237cb94ee6SHong Zhang 
6247cb94ee6SHong Zhang    Level: advanced
6257cb94ee6SHong Zhang 
6267cb94ee6SHong Zhang    Notes:
6277cb94ee6SHong Zhang     These return the number since the creation of the TS object
6287cb94ee6SHong Zhang 
6297cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
6307cb94ee6SHong Zhang 
6317cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6327cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6337cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
634*a43b19c4SJed Brown           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
6357cb94ee6SHong Zhang 
6367cb94ee6SHong Zhang @*/
6377087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
6387cb94ee6SHong Zhang {
6394ac538c5SBarry Smith   PetscErrorCode ierr;
6407cb94ee6SHong Zhang 
6417cb94ee6SHong Zhang   PetscFunctionBegin;
6424ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
6437cb94ee6SHong Zhang   PetscFunctionReturn(0);
6447cb94ee6SHong Zhang }
6457cb94ee6SHong Zhang 
6467cb94ee6SHong Zhang #undef __FUNCT__
6477cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
6487cb94ee6SHong Zhang /*@
6497cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
6507cb94ee6SHong Zhang 
6513f9fe445SBarry Smith    Logically Collective on TS
6527cb94ee6SHong Zhang 
6537cb94ee6SHong Zhang    Input parameters:
6547cb94ee6SHong Zhang +    ts     - the time-step context
6557cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
6567cb94ee6SHong Zhang 
6577cb94ee6SHong Zhang    Level: intermediate
6587cb94ee6SHong Zhang 
6597cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
6607cb94ee6SHong Zhang 
6617cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
6627cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6637cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6647cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
665*a43b19c4SJed Brown           TSSetExactFinalTime()
6667cb94ee6SHong Zhang @*/
6677087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
6687cb94ee6SHong Zhang {
6694ac538c5SBarry Smith   PetscErrorCode ierr;
6707cb94ee6SHong Zhang 
6717cb94ee6SHong Zhang   PetscFunctionBegin;
6724ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
6737cb94ee6SHong Zhang   PetscFunctionReturn(0);
6747cb94ee6SHong Zhang }
6757cb94ee6SHong Zhang 
6767cb94ee6SHong Zhang #undef __FUNCT__
6777cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart"
6787cb94ee6SHong Zhang /*@
6797cb94ee6SHong Zhang    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
6807cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
6817cb94ee6SHong Zhang        this is ALSO the maximum number of GMRES steps that will be used.
6827cb94ee6SHong Zhang 
6833f9fe445SBarry Smith    Logically Collective on TS
6847cb94ee6SHong Zhang 
6857cb94ee6SHong Zhang    Input parameters:
6867cb94ee6SHong Zhang +    ts      - the time-step context
6877cb94ee6SHong Zhang -    restart - number of direction vectors (the restart size).
6887cb94ee6SHong Zhang 
6897cb94ee6SHong Zhang    Level: advanced
6907cb94ee6SHong Zhang 
6917cb94ee6SHong Zhang .keywords: GMRES, restart
6927cb94ee6SHong Zhang 
6937cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
6947cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6957cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6967cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
697*a43b19c4SJed Brown           TSSetExactFinalTime()
6987cb94ee6SHong Zhang 
6997cb94ee6SHong Zhang @*/
7007087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart(TS ts,int restart)
7017cb94ee6SHong Zhang {
7024ac538c5SBarry Smith   PetscErrorCode ierr;
7037cb94ee6SHong Zhang 
7047cb94ee6SHong Zhang   PetscFunctionBegin;
705c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,restart,2);
7064ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr);
7077cb94ee6SHong Zhang   PetscFunctionReturn(0);
7087cb94ee6SHong Zhang }
7097cb94ee6SHong Zhang 
7107cb94ee6SHong Zhang #undef __FUNCT__
7117cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
7127cb94ee6SHong Zhang /*@
7137cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
7147cb94ee6SHong Zhang        system by SUNDIALS.
7157cb94ee6SHong Zhang 
7163f9fe445SBarry Smith    Logically Collective on TS
7177cb94ee6SHong Zhang 
7187cb94ee6SHong Zhang    Input parameters:
7197cb94ee6SHong Zhang +    ts     - the time-step context
7207cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
7217cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
7227cb94ee6SHong Zhang 
7237cb94ee6SHong Zhang    Level: advanced
7247cb94ee6SHong Zhang 
7257cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
7267cb94ee6SHong Zhang 
7277cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7287cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7297cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7307cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
731*a43b19c4SJed Brown           TSSetExactFinalTime()
7327cb94ee6SHong Zhang 
7337cb94ee6SHong Zhang @*/
7347087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
7357cb94ee6SHong Zhang {
7364ac538c5SBarry Smith   PetscErrorCode ierr;
7377cb94ee6SHong Zhang 
7387cb94ee6SHong Zhang   PetscFunctionBegin;
739c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,tol,2);
7404ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
7417cb94ee6SHong Zhang   PetscFunctionReturn(0);
7427cb94ee6SHong Zhang }
7437cb94ee6SHong Zhang 
7447cb94ee6SHong Zhang #undef __FUNCT__
7457cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
7467cb94ee6SHong Zhang /*@
7477cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
7487cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
7497cb94ee6SHong Zhang 
7503f9fe445SBarry Smith    Logically Collective on TS
7517cb94ee6SHong Zhang 
7527cb94ee6SHong Zhang    Input parameters:
7537cb94ee6SHong Zhang +    ts  - the time-step context
7547cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
7557cb94ee6SHong Zhang 
7567cb94ee6SHong Zhang    Level: advanced
7577cb94ee6SHong Zhang 
7587cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
7597cb94ee6SHong Zhang 
7607cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7617cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
7627cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7637cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
764*a43b19c4SJed Brown           TSSetExactFinalTime()
7657cb94ee6SHong Zhang 
7667cb94ee6SHong Zhang @*/
7677087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
7687cb94ee6SHong Zhang {
7694ac538c5SBarry Smith   PetscErrorCode ierr;
7707cb94ee6SHong Zhang 
7717cb94ee6SHong Zhang   PetscFunctionBegin;
7724ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
7737cb94ee6SHong Zhang   PetscFunctionReturn(0);
7747cb94ee6SHong Zhang }
7757cb94ee6SHong Zhang 
7767cb94ee6SHong Zhang #undef __FUNCT__
7777cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
7787cb94ee6SHong Zhang /*@
7797cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
7807cb94ee6SHong Zhang                          Sundials for error control.
7817cb94ee6SHong Zhang 
7823f9fe445SBarry Smith    Logically Collective on TS
7837cb94ee6SHong Zhang 
7847cb94ee6SHong Zhang    Input parameters:
7857cb94ee6SHong Zhang +    ts  - the time-step context
7867cb94ee6SHong Zhang .    aabs - the absolute tolerance
7877cb94ee6SHong Zhang -    rel - the relative tolerance
7887cb94ee6SHong Zhang 
7897cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
7907cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
7917cb94ee6SHong Zhang 
7927cb94ee6SHong Zhang    Level: intermediate
7937cb94ee6SHong Zhang 
7947cb94ee6SHong Zhang .keywords: Sundials, tolerance
7957cb94ee6SHong Zhang 
7967cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7977cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
7987cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7997cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
800*a43b19c4SJed Brown           TSSetExactFinalTime()
8017cb94ee6SHong Zhang 
8027cb94ee6SHong Zhang @*/
8037087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
8047cb94ee6SHong Zhang {
8054ac538c5SBarry Smith   PetscErrorCode ierr;
8067cb94ee6SHong Zhang 
8077cb94ee6SHong Zhang   PetscFunctionBegin;
8084ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
8097cb94ee6SHong Zhang   PetscFunctionReturn(0);
8107cb94ee6SHong Zhang }
8117cb94ee6SHong Zhang 
8127cb94ee6SHong Zhang #undef __FUNCT__
8137cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
8147cb94ee6SHong Zhang /*@
8157cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
8167cb94ee6SHong Zhang 
8177cb94ee6SHong Zhang    Input Parameter:
8187cb94ee6SHong Zhang .    ts - the time-step context
8197cb94ee6SHong Zhang 
8207cb94ee6SHong Zhang    Output Parameter:
8217cb94ee6SHong Zhang .    pc - the preconditioner context
8227cb94ee6SHong Zhang 
8237cb94ee6SHong Zhang    Level: advanced
8247cb94ee6SHong Zhang 
8257cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8267cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8277cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8287cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
8297cb94ee6SHong Zhang @*/
8307087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
8317cb94ee6SHong Zhang {
8324ac538c5SBarry Smith   PetscErrorCode ierr;
8337cb94ee6SHong Zhang 
8347cb94ee6SHong Zhang   PetscFunctionBegin;
8354ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
8367cb94ee6SHong Zhang   PetscFunctionReturn(0);
8377cb94ee6SHong Zhang }
8387cb94ee6SHong Zhang 
8397cb94ee6SHong Zhang #undef __FUNCT__
840f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep"
841f1cd61daSJed Brown /*@
842f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
843f1cd61daSJed Brown 
844f1cd61daSJed Brown    Input Parameter:
845f1cd61daSJed Brown +   ts - the time-step context
846f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
847f1cd61daSJed Brown 
848fc6b9e64SJed Brown    Note:
849fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
850fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
851fc6b9e64SJed Brown 
852f1cd61daSJed Brown    Level: beginner
853f1cd61daSJed Brown 
854f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
855f1cd61daSJed Brown @*/
8567087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
857f1cd61daSJed Brown {
8584ac538c5SBarry Smith   PetscErrorCode ierr;
859f1cd61daSJed Brown 
860f1cd61daSJed Brown   PetscFunctionBegin;
8614ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
862f1cd61daSJed Brown   PetscFunctionReturn(0);
863f1cd61daSJed Brown }
864f1cd61daSJed Brown 
865f1cd61daSJed Brown #undef __FUNCT__
866f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep"
867f1cd61daSJed Brown /*@
868f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
869f1cd61daSJed Brown 
870f1cd61daSJed Brown    Input Parameter:
871f1cd61daSJed Brown +   ts - the time-step context
872f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
873f1cd61daSJed Brown 
874f1cd61daSJed Brown    Level: beginner
875f1cd61daSJed Brown 
876f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
877f1cd61daSJed Brown @*/
8787087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
879f1cd61daSJed Brown {
8804ac538c5SBarry Smith   PetscErrorCode ierr;
881f1cd61daSJed Brown 
882f1cd61daSJed Brown   PetscFunctionBegin;
8834ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
884f1cd61daSJed Brown   PetscFunctionReturn(0);
885f1cd61daSJed Brown }
886f1cd61daSJed Brown 
887f1cd61daSJed Brown #undef __FUNCT__
8882bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps"
8892bfc04deSHong Zhang /*@
8902bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
8912bfc04deSHong Zhang 
8922bfc04deSHong Zhang    Input Parameter:
8932bfc04deSHong Zhang +   ts - the time-step context
8942bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
8952bfc04deSHong Zhang 
8962bfc04deSHong Zhang    Level: beginner
8972bfc04deSHong Zhang 
8982bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8992bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
9002bfc04deSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9012bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
9022bfc04deSHong Zhang @*/
9037087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
9042bfc04deSHong Zhang {
9054ac538c5SBarry Smith   PetscErrorCode ierr;
9062bfc04deSHong Zhang 
9072bfc04deSHong Zhang   PetscFunctionBegin;
9084ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
9092bfc04deSHong Zhang   PetscFunctionReturn(0);
9102bfc04deSHong Zhang }
9117cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
9127cb94ee6SHong Zhang /*MC
91396f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
9147cb94ee6SHong Zhang 
9157cb94ee6SHong Zhang    Options Database:
9167cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
9177cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
9187cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
9197cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
9207cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
9217cb94ee6SHong Zhang .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
92216016685SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
92316016685SHong Zhang 
9247cb94ee6SHong Zhang 
9257cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
9267cb94ee6SHong Zhang            only PETSc PC options
9277cb94ee6SHong Zhang 
9287cb94ee6SHong Zhang     Level: beginner
9297cb94ee6SHong Zhang 
9307cb94ee6SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
931*a43b19c4SJed Brown            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
9327cb94ee6SHong Zhang 
9337cb94ee6SHong Zhang M*/
9347cb94ee6SHong Zhang EXTERN_C_BEGIN
9357cb94ee6SHong Zhang #undef __FUNCT__
9367cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
9377087cfbeSBarry Smith PetscErrorCode  TSCreate_Sundials(TS ts)
9387cb94ee6SHong Zhang {
9397cb94ee6SHong Zhang   TS_Sundials *cvode;
9407cb94ee6SHong Zhang   PetscErrorCode ierr;
9417cb94ee6SHong Zhang 
9427cb94ee6SHong Zhang   PetscFunctionBegin;
943277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
94428adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
94528adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
946214bc6a2SJed Brown   ts->ops->setup          = TSSetUp_Sundials;
947117ce3f2SJed Brown   ts->ops->solve          = TSSolve_Sundials;
948214bc6a2SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
9497cb94ee6SHong Zhang 
95038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
9517cb94ee6SHong Zhang   ts->data          = (void*)cvode;
9526fadb2cdSHong Zhang   cvode->cvode_type = SUNDIALS_BDF;
9536fadb2cdSHong Zhang   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
9547cb94ee6SHong Zhang   cvode->restart    = 5;
9557cb94ee6SHong Zhang   cvode->linear_tol = .05;
9567cb94ee6SHong Zhang 
9572bfc04deSHong Zhang   cvode->monitorstep   = PETSC_FALSE;
9587cb94ee6SHong Zhang 
959a07356b8SSatish Balay   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
960f1cd61daSJed Brown 
961f1cd61daSJed Brown   cvode->mindt = -1.;
962f1cd61daSJed Brown   cvode->maxdt = -1.;
963f1cd61daSJed Brown 
9647cb94ee6SHong Zhang   /* set tolerance for Sundials */
9657cb94ee6SHong Zhang   cvode->reltol = 1e-6;
9662c823083SHong Zhang   cvode->abstol = 1e-6;
9677cb94ee6SHong Zhang 
968*a43b19c4SJed Brown   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_TRUE;
969*a43b19c4SJed Brown 
9707cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
9717cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
9727cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
9737cb94ee6SHong Zhang                     "TSSundialsSetGMRESRestart_Sundials",
9747cb94ee6SHong Zhang                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
9757cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
9767cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
9777cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
9787cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
9797cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
9807cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
9817cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
9827cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
9837cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
984f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
985f1cd61daSJed Brown                     "TSSundialsSetMinTimeStep_Sundials",
986f1cd61daSJed Brown                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
987f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
988f1cd61daSJed Brown                     "TSSundialsSetMaxTimeStep_Sundials",
989f1cd61daSJed Brown                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
9907cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
9917cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
9927cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
9937cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
9947cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
9957cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
9962bfc04deSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
9972bfc04deSHong Zhang                     "TSSundialsMonitorInternalSteps_Sundials",
9982bfc04deSHong Zhang                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
9992bfc04deSHong Zhang 
10007cb94ee6SHong Zhang   PetscFunctionReturn(0);
10017cb94ee6SHong Zhang }
10027cb94ee6SHong Zhang EXTERN_C_END
1003