xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision d155d4104c12796cc6e93a7e2b3cd74e97fcda36)
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;
237cb94ee6SHong Zhang   PC             pc = cvode->pc;
247cb94ee6SHong Zhang   PetscErrorCode ierr;
257cb94ee6SHong Zhang   Mat            Jac = ts->B;
267cb94ee6SHong Zhang   Vec            yy = cvode->w1;
277cb94ee6SHong Zhang   PetscScalar    one = 1.0,gm;
287cb94ee6SHong Zhang   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
297cb94ee6SHong Zhang   PetscScalar    *y_data;
307cb94ee6SHong Zhang 
317cb94ee6SHong Zhang   PetscFunctionBegin;
327cb94ee6SHong Zhang   /* This allows us to construct preconditioners in-place if we like */
337cb94ee6SHong Zhang   ierr = MatSetUnfactored(Jac);CHKERRQ(ierr);
347cb94ee6SHong Zhang 
357cb94ee6SHong Zhang   /* jok - TRUE means reuse current Jacobian else recompute Jacobian */
367cb94ee6SHong Zhang   if (jok) {
377cb94ee6SHong Zhang     ierr     = MatCopy(cvode->pmat,Jac,str);CHKERRQ(ierr);
387cb94ee6SHong Zhang     *jcurPtr = FALSE;
397cb94ee6SHong Zhang   } else {
407cb94ee6SHong Zhang     /* make PETSc vector yy point to SUNDIALS vector y */
417cb94ee6SHong Zhang     y_data = (PetscScalar *) N_VGetArrayPointer(y);
427cb94ee6SHong Zhang     ierr   = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
434ac7836bSHong Zhang 
447cb94ee6SHong Zhang     /* compute the Jacobian */
457cb94ee6SHong Zhang     ierr = TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);CHKERRQ(ierr);
467cb94ee6SHong Zhang     ierr = VecResetArray(yy); CHKERRQ(ierr);
474ac7836bSHong Zhang 
487cb94ee6SHong Zhang     /* copy the Jacobian matrix */
497cb94ee6SHong Zhang     if (!cvode->pmat) {
507cb94ee6SHong Zhang       ierr = MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);CHKERRQ(ierr);
517cb94ee6SHong Zhang       ierr = PetscLogObjectParent(ts,cvode->pmat);CHKERRQ(ierr);
522c823083SHong Zhang     } else {
537cb94ee6SHong Zhang       ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr);
547cb94ee6SHong Zhang     }
557cb94ee6SHong Zhang     *jcurPtr = TRUE;
567cb94ee6SHong Zhang   }
577cb94ee6SHong Zhang 
587cb94ee6SHong Zhang   /* construct I-gamma*Jac  */
597cb94ee6SHong Zhang   gm   = -_gamma;
607cb94ee6SHong Zhang   ierr = MatScale(Jac,gm);CHKERRQ(ierr);
617cb94ee6SHong Zhang   ierr = MatShift(Jac,one);CHKERRQ(ierr);
627cb94ee6SHong Zhang 
637cb94ee6SHong Zhang   ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr);
647cb94ee6SHong Zhang   PetscFunctionReturn(0);
657cb94ee6SHong Zhang }
667cb94ee6SHong Zhang 
677cb94ee6SHong Zhang /*
687cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
697cb94ee6SHong Zhang */
707cb94ee6SHong Zhang #undef __FUNCT__
717cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials"
724ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
734ac7836bSHong Zhang                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
747cb94ee6SHong Zhang {
757cb94ee6SHong Zhang   TS              ts = (TS) P_data;
767cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
777cb94ee6SHong Zhang   PC              pc = cvode->pc;
787cb94ee6SHong Zhang   Vec             rr = cvode->w1,zz = cvode->w2;
797cb94ee6SHong Zhang   PetscErrorCode  ierr;
807cb94ee6SHong Zhang   PetscScalar     *r_data,*z_data;
817cb94ee6SHong Zhang 
827cb94ee6SHong Zhang   PetscFunctionBegin;
837cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
847cb94ee6SHong Zhang   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
857cb94ee6SHong Zhang   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
867cb94ee6SHong Zhang   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
877cb94ee6SHong Zhang   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
884ac7836bSHong Zhang 
897cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
907cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
917cb94ee6SHong Zhang   ierr = VecResetArray(rr); CHKERRQ(ierr);
927cb94ee6SHong Zhang   ierr = VecResetArray(zz); CHKERRQ(ierr);
937cb94ee6SHong Zhang   PetscFunctionReturn(0);
947cb94ee6SHong Zhang }
957cb94ee6SHong Zhang 
967cb94ee6SHong Zhang /*
977cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
987cb94ee6SHong Zhang */
997cb94ee6SHong Zhang #undef __FUNCT__
1007cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials"
1014ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
1027cb94ee6SHong Zhang {
1037cb94ee6SHong Zhang   TS              ts = (TS) ctx;
1047adad957SLisandro Dalcin   MPI_Comm        comm = ((PetscObject)ts)->comm;
1057cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
1067cb94ee6SHong Zhang   Vec             yy = cvode->w1,yyd = cvode->w2;
1077cb94ee6SHong Zhang   PetscScalar     *y_data,*ydot_data;
1087cb94ee6SHong Zhang   PetscErrorCode  ierr;
1097cb94ee6SHong Zhang 
1107cb94ee6SHong Zhang   PetscFunctionBegin;
1115ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
1127cb94ee6SHong Zhang   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
1137cb94ee6SHong Zhang   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
1144ac7836bSHong Zhang   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
1154ac7836bSHong Zhang   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
1164ac7836bSHong Zhang 
1177cb94ee6SHong Zhang   /* now compute the right hand side function */
1187cb94ee6SHong Zhang   ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr);
1197cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
1207cb94ee6SHong Zhang   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
1214ac7836bSHong Zhang   PetscFunctionReturn(0);
1227cb94ee6SHong Zhang }
1237cb94ee6SHong Zhang 
1247cb94ee6SHong Zhang /*
1257cb94ee6SHong Zhang        TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE.
1267cb94ee6SHong Zhang */
1277cb94ee6SHong Zhang #undef __FUNCT__
1287cb94ee6SHong Zhang #define __FUNCT__ "TSStep_Sundials_Nonlinear"
1297cb94ee6SHong Zhang PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
1307cb94ee6SHong Zhang {
1317cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1327cb94ee6SHong Zhang   Vec            sol = ts->vec_sol;
1337cb94ee6SHong Zhang   PetscErrorCode ierr;
134186e87acSLisandro Dalcin   PetscInt       i,flag;
1359f94935aSHong Zhang   long int       its,nsteps;
1367cb94ee6SHong Zhang   realtype       t,tout;
1377cb94ee6SHong Zhang   PetscScalar    *y_data;
1387cb94ee6SHong Zhang   void           *mem;
1397cb94ee6SHong Zhang 
1407cb94ee6SHong Zhang   PetscFunctionBegin;
14116016685SHong Zhang   mem  = cvode->mem;
1427cb94ee6SHong Zhang   tout = ts->max_time;
1437cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1447cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data,cvode->y);
1457cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
146186e87acSLisandro Dalcin 
147186e87acSLisandro Dalcin   for (i = 0; i < ts->max_steps; i++) {
1487cb94ee6SHong Zhang     if (ts->ptime >= ts->max_time) break;
1493f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
150186e87acSLisandro Dalcin 
1512bfc04deSHong Zhang     if (cvode->monitorstep) {
15228adb3f7SHong Zhang       flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
1532bfc04deSHong Zhang     } else {
1542bfc04deSHong Zhang       flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
1552bfc04deSHong Zhang     }
1569f94935aSHong Zhang 
1579f94935aSHong Zhang     if (flag){ /* display error message */
1589f94935aSHong Zhang       switch (flag){
1599f94935aSHong Zhang       case CV_ILL_INPUT:
1609f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
1619f94935aSHong Zhang         break;
1629f94935aSHong Zhang       case CV_TOO_CLOSE:
1639f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
1649f94935aSHong Zhang         break;
1653c7fefeeSJed Brown       case CV_TOO_MUCH_WORK: {
1669f94935aSHong Zhang         PetscReal      tcur;
1679f94935aSHong Zhang         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
1689f94935aSHong Zhang         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
1699f94935aSHong 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);
1703c7fefeeSJed Brown       } break;
1719f94935aSHong Zhang       case CV_TOO_MUCH_ACC:
1729f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
1739f94935aSHong Zhang         break;
1749f94935aSHong Zhang       case CV_ERR_FAILURE:
1759f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
1769f94935aSHong Zhang         break;
1779f94935aSHong Zhang       case CV_CONV_FAILURE:
1789f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
1799f94935aSHong Zhang         break;
1809f94935aSHong Zhang       case CV_LINIT_FAIL:
1819f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
1829f94935aSHong Zhang         break;
1839f94935aSHong Zhang       case CV_LSETUP_FAIL:
1849f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
1859f94935aSHong Zhang         break;
1869f94935aSHong Zhang       case CV_LSOLVE_FAIL:
1879f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
1889f94935aSHong Zhang         break;
1899f94935aSHong Zhang       case CV_RHSFUNC_FAIL:
1909f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
1919f94935aSHong Zhang         break;
1929f94935aSHong Zhang       case CV_FIRST_RHSFUNC_ERR:
1939f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
1949f94935aSHong Zhang         break;
1959f94935aSHong Zhang       case CV_REPTD_RHSFUNC_ERR:
1969f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
1979f94935aSHong Zhang         break;
1989f94935aSHong Zhang       case CV_UNREC_RHSFUNC_ERR:
1999f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
2009f94935aSHong Zhang         break;
2019f94935aSHong Zhang       case CV_RTFUNC_FAIL:
2029f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
2039f94935aSHong Zhang         break;
2049f94935aSHong Zhang       default:
2059f94935aSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
2069f94935aSHong Zhang       }
2079f94935aSHong Zhang     }
2089f94935aSHong Zhang 
2097cb94ee6SHong Zhang     if (t > ts->max_time && cvode->exact_final_time) {
2107cb94ee6SHong Zhang       /* interpolate to final requested time */
2117cb94ee6SHong Zhang       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
2127cb94ee6SHong Zhang       t = tout;
2137cb94ee6SHong Zhang     }
2147cb94ee6SHong Zhang 
2157cb94ee6SHong Zhang     /* copy the solution from cvode->y to cvode->update and sol */
2167cb94ee6SHong Zhang     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
2177cb94ee6SHong Zhang     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
2187cb94ee6SHong Zhang     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
2197cb94ee6SHong Zhang     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
2207cb94ee6SHong Zhang     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
2214ac7836bSHong Zhang     ierr = CVSpilsGetNumLinIters(mem, &its);
222186e87acSLisandro Dalcin     ts->nonlinear_its = its; ts->linear_its = its;
223186e87acSLisandro Dalcin 
224186e87acSLisandro Dalcin     ts->time_step = t - ts->ptime;
225186e87acSLisandro Dalcin     ts->ptime     = t;
2267cb94ee6SHong Zhang     ts->steps++;
227186e87acSLisandro Dalcin 
2283f2090d5SJed Brown     ierr = TSPostStep(ts);CHKERRQ(ierr);
2297cb94ee6SHong Zhang     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
2307cb94ee6SHong Zhang   }
2319f94935aSHong Zhang   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
232186e87acSLisandro Dalcin 
2339f94935aSHong Zhang   *steps = nsteps;
2347cb94ee6SHong Zhang   *time  = t;
2357cb94ee6SHong Zhang   PetscFunctionReturn(0);
2367cb94ee6SHong Zhang }
2377cb94ee6SHong Zhang 
2387cb94ee6SHong Zhang #undef __FUNCT__
239277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials"
240277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts)
241277b19d0SLisandro Dalcin {
242277b19d0SLisandro Dalcin   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
243277b19d0SLisandro Dalcin   PetscErrorCode ierr;
244277b19d0SLisandro Dalcin 
245277b19d0SLisandro Dalcin   PetscFunctionBegin;
246277b19d0SLisandro Dalcin   if (cvode->pc)     {ierr = PCReset(cvode->pc);CHKERRQ(ierr);}
247*d155d410SJed Brown   ierr = MatDestroy(&cvode->pmat);
248*d155d410SJed Brown   ierr = VecDestroy(&cvode->update);
249*d155d410SJed Brown   ierr = VecDestroy(&cvode->func);
250*d155d410SJed Brown   ierr = VecDestroy(&cvode->rhs);
251*d155d410SJed Brown   ierr = VecDestroy(&cvode->w1);
252*d155d410SJed Brown   ierr = VecDestroy(&cvode->w2);
253277b19d0SLisandro Dalcin   if (cvode->mem)    {CVodeFree(&cvode->mem);}
254277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
255277b19d0SLisandro Dalcin }
256277b19d0SLisandro Dalcin 
257277b19d0SLisandro Dalcin #undef __FUNCT__
2587cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials"
2597cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
2607cb94ee6SHong Zhang {
2617cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2627cb94ee6SHong Zhang   PetscErrorCode ierr;
2637cb94ee6SHong Zhang 
2647cb94ee6SHong Zhang   PetscFunctionBegin;
265277b19d0SLisandro Dalcin   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
266*d155d410SJed Brown   ierr = PCDestroy(&cvode->pc);CHKERRQ(ierr);
267a07356b8SSatish Balay   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
268277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
2697cb94ee6SHong Zhang   PetscFunctionReturn(0);
2707cb94ee6SHong Zhang }
2717cb94ee6SHong Zhang 
2727cb94ee6SHong Zhang #undef __FUNCT__
2737cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear"
2747cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
2757cb94ee6SHong Zhang {
2767cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2777cb94ee6SHong Zhang   PetscErrorCode ierr;
27816016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2797cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
28016016685SHong Zhang   void           *mem;
28116016685SHong Zhang   const PCType   pctype;
282ace3abfcSBarry Smith   PetscBool      pcnone;
28316016685SHong Zhang   Vec            sol = ts->vec_sol;
2847cb94ee6SHong Zhang 
2857cb94ee6SHong Zhang   PetscFunctionBegin;
2867cb94ee6SHong Zhang   ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr);
2879f94935aSHong Zhang 
2887cb94ee6SHong Zhang   /* get the vector size */
2897cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2907cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
2917cb94ee6SHong Zhang 
2927cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
293a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
294e32f2f54SBarry Smith   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
2957cb94ee6SHong Zhang 
29628adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
2977cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
2987cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
2997cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
3007cb94ee6SHong Zhang   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
3017cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
3027cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
3037cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr);
3047cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
3057cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr);
3067cb94ee6SHong Zhang 
3077cb94ee6SHong Zhang   /*
3087cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
3097cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
3107cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
3117cb94ee6SHong Zhang   */
3127adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
3137adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
3147cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
3157cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
31616016685SHong Zhang 
31716016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
31816016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
319e32f2f54SBarry Smith   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
32016016685SHong Zhang   cvode->mem = mem;
32116016685SHong Zhang 
32216016685SHong Zhang   /* Set the pointer to user-defined data */
32316016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
324e32f2f54SBarry Smith   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
32516016685SHong Zhang 
326fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
327fc6b9e64SJed Brown   flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step);
328fc6b9e64SJed Brown   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
329f1cd61daSJed Brown   if (cvode->mindt > 0) {
330f1cd61daSJed Brown     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
3319f94935aSHong Zhang     if (flag){
3329f94935aSHong Zhang       if (flag == CV_MEM_NULL){
3339f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
3349f94935aSHong Zhang       } else if (flag == CV_ILL_INPUT){
3359f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
3369f94935aSHong Zhang       } else {
3379f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
3389f94935aSHong Zhang       }
3399f94935aSHong Zhang     }
340f1cd61daSJed Brown   }
341f1cd61daSJed Brown   if (cvode->maxdt > 0) {
342f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
343f1cd61daSJed Brown     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
344f1cd61daSJed Brown   }
345f1cd61daSJed Brown 
34616016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
34716016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
34816016685SHong Zhang    * the initial dependent variable vector cvode->y */
34916016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
35016016685SHong Zhang   if (flag){
351e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
35216016685SHong Zhang   }
35316016685SHong Zhang 
3549f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
35516016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
35616016685SHong Zhang   if (flag){
357e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
35816016685SHong Zhang   }
35916016685SHong Zhang 
3609f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3619f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
36216016685SHong Zhang   ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
36316016685SHong Zhang 
36416016685SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
36516016685SHong Zhang   /* setup the ode integrator with the given preconditioner */
36616016685SHong Zhang   ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr);
36716016685SHong Zhang   ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr);
36816016685SHong Zhang   if (pcnone){
36916016685SHong Zhang     flag  = CVSpgmr(mem,PREC_NONE,0);
370e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
37116016685SHong Zhang   } else {
37216016685SHong Zhang     flag  = CVSpgmr(mem,PREC_LEFT,0);
373e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
37416016685SHong Zhang 
37516016685SHong Zhang     /* Set preconditioner and solve routines Precond and PSolve,
37616016685SHong Zhang      and the pointer to the user-defined block data */
37716016685SHong Zhang     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
378e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
37916016685SHong Zhang   }
38016016685SHong Zhang 
38116016685SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
382e32f2f54SBarry Smith   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
3837cb94ee6SHong Zhang   PetscFunctionReturn(0);
3847cb94ee6SHong Zhang }
3857cb94ee6SHong Zhang 
3866fadb2cdSHong Zhang /* type of CVODE linear multistep method */
387dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
3886fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
389dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
390a04cf4d8SBarry Smith 
3917cb94ee6SHong Zhang #undef __FUNCT__
3927cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear"
3937cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
3947cb94ee6SHong Zhang {
3957cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3967cb94ee6SHong Zhang   PetscErrorCode ierr;
3977cb94ee6SHong Zhang   int            indx;
398ace3abfcSBarry Smith   PetscBool      flag;
3997cb94ee6SHong Zhang 
4007cb94ee6SHong Zhang   PetscFunctionBegin;
4017cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
4026fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
4037cb94ee6SHong Zhang     if (flag) {
4046fadb2cdSHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
4057cb94ee6SHong Zhang     }
4066fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
4077cb94ee6SHong Zhang     if (flag) {
4087cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
4097cb94ee6SHong Zhang     }
4107cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
4117cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
412f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
413f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
4147cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
4157cb94ee6SHong Zhang     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
416acfcf0e5SJed 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);
417acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
4187cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4197cb94ee6SHong Zhang   PetscFunctionReturn(0);
4207cb94ee6SHong Zhang }
4217cb94ee6SHong Zhang 
4227cb94ee6SHong Zhang #undef __FUNCT__
4237cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
4247cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
4257cb94ee6SHong Zhang {
4267cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4277cb94ee6SHong Zhang   PetscErrorCode ierr;
4287cb94ee6SHong Zhang   char           *type;
4297cb94ee6SHong Zhang   char           atype[] = "Adams";
4307cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
431ace3abfcSBarry Smith   PetscBool      iascii,isstring;
4322c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
4332c823083SHong Zhang   PetscInt       qlast,qcur;
4345d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
4357cb94ee6SHong Zhang 
4367cb94ee6SHong Zhang   PetscFunctionBegin;
4377cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
4387cb94ee6SHong Zhang   else                                     {type = btype;}
4397cb94ee6SHong Zhang 
4402692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4412692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
4427cb94ee6SHong Zhang   if (iascii) {
4437cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
4447cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
4457cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
4467cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
4477cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
4487cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4497cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4507cb94ee6SHong Zhang     } else {
4517cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4527cb94ee6SHong Zhang     }
453f1cd61daSJed Brown     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
454f1cd61daSJed Brown     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
4552c823083SHong Zhang 
4565d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4575d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
4585d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
4592c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
4602c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
4612c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
4622c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
4632c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
4642c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
4652c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
4662c823083SHong Zhang 
4672c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
4682c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
4692c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
4702c823083SHong Zhang 
4712c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
4722c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
4732c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
4742c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
4752c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4762c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
4772c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
4782c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
4792c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4802c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
4812c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4822c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
4837cb94ee6SHong Zhang   } else if (isstring) {
4847cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
4857cb94ee6SHong Zhang   } else {
486e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
4877cb94ee6SHong Zhang   }
4887cb94ee6SHong Zhang   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
4897cb94ee6SHong Zhang   ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr);
4907cb94ee6SHong Zhang   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
4917cb94ee6SHong Zhang   PetscFunctionReturn(0);
4927cb94ee6SHong Zhang }
4937cb94ee6SHong Zhang 
4947cb94ee6SHong Zhang 
4957cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
4967cb94ee6SHong Zhang EXTERN_C_BEGIN
4977cb94ee6SHong Zhang #undef __FUNCT__
4987cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
4997087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
5007cb94ee6SHong Zhang {
5017cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5027cb94ee6SHong Zhang 
5037cb94ee6SHong Zhang   PetscFunctionBegin;
5047cb94ee6SHong Zhang   cvode->cvode_type = type;
5057cb94ee6SHong Zhang   PetscFunctionReturn(0);
5067cb94ee6SHong Zhang }
5077cb94ee6SHong Zhang EXTERN_C_END
5087cb94ee6SHong Zhang 
5097cb94ee6SHong Zhang EXTERN_C_BEGIN
5107cb94ee6SHong Zhang #undef __FUNCT__
5117cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
5127087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
5137cb94ee6SHong Zhang {
5147cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5157cb94ee6SHong Zhang 
5167cb94ee6SHong Zhang   PetscFunctionBegin;
5177cb94ee6SHong Zhang   cvode->restart = restart;
5187cb94ee6SHong Zhang   PetscFunctionReturn(0);
5197cb94ee6SHong Zhang }
5207cb94ee6SHong Zhang EXTERN_C_END
5217cb94ee6SHong Zhang 
5227cb94ee6SHong Zhang EXTERN_C_BEGIN
5237cb94ee6SHong Zhang #undef __FUNCT__
5247cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
5257087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
5267cb94ee6SHong Zhang {
5277cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5287cb94ee6SHong Zhang 
5297cb94ee6SHong Zhang   PetscFunctionBegin;
5307cb94ee6SHong Zhang   cvode->linear_tol = tol;
5317cb94ee6SHong Zhang   PetscFunctionReturn(0);
5327cb94ee6SHong Zhang }
5337cb94ee6SHong Zhang EXTERN_C_END
5347cb94ee6SHong Zhang 
5357cb94ee6SHong Zhang EXTERN_C_BEGIN
5367cb94ee6SHong Zhang #undef __FUNCT__
5377cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
5387087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
5397cb94ee6SHong Zhang {
5407cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5417cb94ee6SHong Zhang 
5427cb94ee6SHong Zhang   PetscFunctionBegin;
5437cb94ee6SHong Zhang   cvode->gtype = type;
5447cb94ee6SHong Zhang   PetscFunctionReturn(0);
5457cb94ee6SHong Zhang }
5467cb94ee6SHong Zhang EXTERN_C_END
5477cb94ee6SHong Zhang 
5487cb94ee6SHong Zhang EXTERN_C_BEGIN
5497cb94ee6SHong Zhang #undef __FUNCT__
5507cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
5517087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
5527cb94ee6SHong Zhang {
5537cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5547cb94ee6SHong Zhang 
5557cb94ee6SHong Zhang   PetscFunctionBegin;
5567cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5577cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
5587cb94ee6SHong Zhang   PetscFunctionReturn(0);
5597cb94ee6SHong Zhang }
5607cb94ee6SHong Zhang EXTERN_C_END
5617cb94ee6SHong Zhang 
5627cb94ee6SHong Zhang EXTERN_C_BEGIN
5637cb94ee6SHong Zhang #undef __FUNCT__
564f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
5657087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
566f1cd61daSJed Brown {
567f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
568f1cd61daSJed Brown 
569f1cd61daSJed Brown   PetscFunctionBegin;
570f1cd61daSJed Brown   cvode->mindt = mindt;
571f1cd61daSJed Brown   PetscFunctionReturn(0);
572f1cd61daSJed Brown }
573f1cd61daSJed Brown EXTERN_C_END
574f1cd61daSJed Brown 
575f1cd61daSJed Brown EXTERN_C_BEGIN
576f1cd61daSJed Brown #undef __FUNCT__
577f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
5787087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
579f1cd61daSJed Brown {
580f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
581f1cd61daSJed Brown 
582f1cd61daSJed Brown   PetscFunctionBegin;
583f1cd61daSJed Brown   cvode->maxdt = maxdt;
584f1cd61daSJed Brown   PetscFunctionReturn(0);
585f1cd61daSJed Brown }
586f1cd61daSJed Brown EXTERN_C_END
587f1cd61daSJed Brown 
588f1cd61daSJed Brown EXTERN_C_BEGIN
589f1cd61daSJed Brown #undef __FUNCT__
5907cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
5917087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
5927cb94ee6SHong Zhang {
5937cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5947cb94ee6SHong Zhang 
5957cb94ee6SHong Zhang   PetscFunctionBegin;
5967cb94ee6SHong Zhang   *pc = cvode->pc;
5977cb94ee6SHong Zhang   PetscFunctionReturn(0);
5987cb94ee6SHong Zhang }
5997cb94ee6SHong Zhang EXTERN_C_END
6007cb94ee6SHong Zhang 
6017cb94ee6SHong Zhang EXTERN_C_BEGIN
6027cb94ee6SHong Zhang #undef __FUNCT__
6037cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
6047087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
6057cb94ee6SHong Zhang {
6067cb94ee6SHong Zhang   PetscFunctionBegin;
6072c823083SHong Zhang   if (nonlin) *nonlin = ts->nonlinear_its;
6082c823083SHong Zhang   if (lin)    *lin    = ts->linear_its;
6097cb94ee6SHong Zhang   PetscFunctionReturn(0);
6107cb94ee6SHong Zhang }
6117cb94ee6SHong Zhang EXTERN_C_END
6127cb94ee6SHong Zhang 
6137cb94ee6SHong Zhang EXTERN_C_BEGIN
6147cb94ee6SHong Zhang #undef __FUNCT__
6157cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
6167087cfbeSBarry Smith PetscErrorCode  TSSundialsSetExactFinalTime_Sundials(TS ts,PetscBool  s)
6177cb94ee6SHong Zhang {
6187cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6197cb94ee6SHong Zhang 
6207cb94ee6SHong Zhang   PetscFunctionBegin;
6217cb94ee6SHong Zhang   cvode->exact_final_time = s;
6227cb94ee6SHong Zhang   PetscFunctionReturn(0);
6237cb94ee6SHong Zhang }
6247cb94ee6SHong Zhang EXTERN_C_END
6252bfc04deSHong Zhang 
6262bfc04deSHong Zhang EXTERN_C_BEGIN
6272bfc04deSHong Zhang #undef __FUNCT__
6282bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
6297087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
6302bfc04deSHong Zhang {
6312bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6322bfc04deSHong Zhang 
6332bfc04deSHong Zhang   PetscFunctionBegin;
6342bfc04deSHong Zhang   cvode->monitorstep = s;
6352bfc04deSHong Zhang   PetscFunctionReturn(0);
6362bfc04deSHong Zhang }
6372bfc04deSHong Zhang EXTERN_C_END
6387cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
6397cb94ee6SHong Zhang 
6407cb94ee6SHong Zhang #undef __FUNCT__
6417cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
6427cb94ee6SHong Zhang /*@C
6437cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
6447cb94ee6SHong Zhang 
6457cb94ee6SHong Zhang    Not Collective
6467cb94ee6SHong Zhang 
6477cb94ee6SHong Zhang    Input parameters:
6487cb94ee6SHong Zhang .    ts     - the time-step context
6497cb94ee6SHong Zhang 
6507cb94ee6SHong Zhang    Output Parameters:
6517cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
6527cb94ee6SHong Zhang -   lin    - number of linear iterations
6537cb94ee6SHong Zhang 
6547cb94ee6SHong Zhang    Level: advanced
6557cb94ee6SHong Zhang 
6567cb94ee6SHong Zhang    Notes:
6577cb94ee6SHong Zhang     These return the number since the creation of the TS object
6587cb94ee6SHong Zhang 
6597cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
6607cb94ee6SHong Zhang 
6617cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6627cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6637cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
664f1cd61daSJed Brown           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSundialsSetExactFinalTime()
6657cb94ee6SHong Zhang 
6667cb94ee6SHong Zhang @*/
6677087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
6687cb94ee6SHong Zhang {
6694ac538c5SBarry Smith   PetscErrorCode ierr;
6707cb94ee6SHong Zhang 
6717cb94ee6SHong Zhang   PetscFunctionBegin;
6724ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
6737cb94ee6SHong Zhang   PetscFunctionReturn(0);
6747cb94ee6SHong Zhang }
6757cb94ee6SHong Zhang 
6767cb94ee6SHong Zhang #undef __FUNCT__
6777cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
6787cb94ee6SHong Zhang /*@
6797cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
6807cb94ee6SHong Zhang 
6813f9fe445SBarry Smith    Logically Collective on TS
6827cb94ee6SHong Zhang 
6837cb94ee6SHong Zhang    Input parameters:
6847cb94ee6SHong Zhang +    ts     - the time-step context
6857cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
6867cb94ee6SHong Zhang 
6877cb94ee6SHong Zhang    Level: intermediate
6887cb94ee6SHong Zhang 
6897cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
6907cb94ee6SHong Zhang 
6917cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
6927cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6937cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6947cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6957cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6967cb94ee6SHong Zhang @*/
6977087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
6987cb94ee6SHong Zhang {
6994ac538c5SBarry Smith   PetscErrorCode ierr;
7007cb94ee6SHong Zhang 
7017cb94ee6SHong Zhang   PetscFunctionBegin;
7024ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
7037cb94ee6SHong Zhang   PetscFunctionReturn(0);
7047cb94ee6SHong Zhang }
7057cb94ee6SHong Zhang 
7067cb94ee6SHong Zhang #undef __FUNCT__
7077cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart"
7087cb94ee6SHong Zhang /*@
7097cb94ee6SHong Zhang    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
7107cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
7117cb94ee6SHong Zhang        this is ALSO the maximum number of GMRES steps that will be used.
7127cb94ee6SHong Zhang 
7133f9fe445SBarry Smith    Logically Collective on TS
7147cb94ee6SHong Zhang 
7157cb94ee6SHong Zhang    Input parameters:
7167cb94ee6SHong Zhang +    ts      - the time-step context
7177cb94ee6SHong Zhang -    restart - number of direction vectors (the restart size).
7187cb94ee6SHong Zhang 
7197cb94ee6SHong Zhang    Level: advanced
7207cb94ee6SHong Zhang 
7217cb94ee6SHong Zhang .keywords: GMRES, restart
7227cb94ee6SHong Zhang 
7237cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
7247cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7257cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7267cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7277cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7287cb94ee6SHong Zhang 
7297cb94ee6SHong Zhang @*/
7307087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart(TS ts,int restart)
7317cb94ee6SHong Zhang {
7324ac538c5SBarry Smith   PetscErrorCode ierr;
7337cb94ee6SHong Zhang 
7347cb94ee6SHong Zhang   PetscFunctionBegin;
735c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,restart,2);
7364ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr);
7377cb94ee6SHong Zhang   PetscFunctionReturn(0);
7387cb94ee6SHong Zhang }
7397cb94ee6SHong Zhang 
7407cb94ee6SHong Zhang #undef __FUNCT__
7417cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
7427cb94ee6SHong Zhang /*@
7437cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
7447cb94ee6SHong Zhang        system by SUNDIALS.
7457cb94ee6SHong Zhang 
7463f9fe445SBarry Smith    Logically Collective on TS
7477cb94ee6SHong Zhang 
7487cb94ee6SHong Zhang    Input parameters:
7497cb94ee6SHong Zhang +    ts     - the time-step context
7507cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
7517cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
7527cb94ee6SHong Zhang 
7537cb94ee6SHong Zhang    Level: advanced
7547cb94ee6SHong Zhang 
7557cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
7567cb94ee6SHong Zhang 
7577cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7587cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7597cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7607cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7617cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7627cb94ee6SHong Zhang 
7637cb94ee6SHong Zhang @*/
7647087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
7657cb94ee6SHong Zhang {
7664ac538c5SBarry Smith   PetscErrorCode ierr;
7677cb94ee6SHong Zhang 
7687cb94ee6SHong Zhang   PetscFunctionBegin;
769c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,tol,2);
7704ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
7717cb94ee6SHong Zhang   PetscFunctionReturn(0);
7727cb94ee6SHong Zhang }
7737cb94ee6SHong Zhang 
7747cb94ee6SHong Zhang #undef __FUNCT__
7757cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
7767cb94ee6SHong Zhang /*@
7777cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
7787cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
7797cb94ee6SHong Zhang 
7803f9fe445SBarry Smith    Logically Collective on TS
7817cb94ee6SHong Zhang 
7827cb94ee6SHong Zhang    Input parameters:
7837cb94ee6SHong Zhang +    ts  - the time-step context
7847cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
7857cb94ee6SHong Zhang 
7867cb94ee6SHong Zhang    Level: advanced
7877cb94ee6SHong Zhang 
7887cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
7897cb94ee6SHong Zhang 
7907cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7917cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
7927cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7937cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7947cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7957cb94ee6SHong Zhang 
7967cb94ee6SHong Zhang @*/
7977087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
7987cb94ee6SHong Zhang {
7994ac538c5SBarry Smith   PetscErrorCode ierr;
8007cb94ee6SHong Zhang 
8017cb94ee6SHong Zhang   PetscFunctionBegin;
8024ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
8037cb94ee6SHong Zhang   PetscFunctionReturn(0);
8047cb94ee6SHong Zhang }
8057cb94ee6SHong Zhang 
8067cb94ee6SHong Zhang #undef __FUNCT__
8077cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
8087cb94ee6SHong Zhang /*@
8097cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
8107cb94ee6SHong Zhang                          Sundials for error control.
8117cb94ee6SHong Zhang 
8123f9fe445SBarry Smith    Logically Collective on TS
8137cb94ee6SHong Zhang 
8147cb94ee6SHong Zhang    Input parameters:
8157cb94ee6SHong Zhang +    ts  - the time-step context
8167cb94ee6SHong Zhang .    aabs - the absolute tolerance
8177cb94ee6SHong Zhang -    rel - the relative tolerance
8187cb94ee6SHong Zhang 
8197cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
8207cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
8217cb94ee6SHong Zhang 
8227cb94ee6SHong Zhang    Level: intermediate
8237cb94ee6SHong Zhang 
8247cb94ee6SHong Zhang .keywords: Sundials, tolerance
8257cb94ee6SHong Zhang 
8267cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8277cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
8287cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8297cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
8307cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
8317cb94ee6SHong Zhang 
8327cb94ee6SHong Zhang @*/
8337087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
8347cb94ee6SHong Zhang {
8354ac538c5SBarry Smith   PetscErrorCode ierr;
8367cb94ee6SHong Zhang 
8377cb94ee6SHong Zhang   PetscFunctionBegin;
8384ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
8397cb94ee6SHong Zhang   PetscFunctionReturn(0);
8407cb94ee6SHong Zhang }
8417cb94ee6SHong Zhang 
8427cb94ee6SHong Zhang #undef __FUNCT__
8437cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
8447cb94ee6SHong Zhang /*@
8457cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
8467cb94ee6SHong Zhang 
8477cb94ee6SHong Zhang    Input Parameter:
8487cb94ee6SHong Zhang .    ts - the time-step context
8497cb94ee6SHong Zhang 
8507cb94ee6SHong Zhang    Output Parameter:
8517cb94ee6SHong Zhang .    pc - the preconditioner context
8527cb94ee6SHong Zhang 
8537cb94ee6SHong Zhang    Level: advanced
8547cb94ee6SHong Zhang 
8557cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8567cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8577cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8587cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
8597cb94ee6SHong Zhang @*/
8607087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
8617cb94ee6SHong Zhang {
8624ac538c5SBarry Smith   PetscErrorCode ierr;
8637cb94ee6SHong Zhang 
8647cb94ee6SHong Zhang   PetscFunctionBegin;
8654ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
8667cb94ee6SHong Zhang   PetscFunctionReturn(0);
8677cb94ee6SHong Zhang }
8687cb94ee6SHong Zhang 
8697cb94ee6SHong Zhang #undef __FUNCT__
8707cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime"
8717cb94ee6SHong Zhang /*@
8727cb94ee6SHong Zhang    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
8737cb94ee6SHong Zhang       exact final time requested by the user or just returns it at the final time
8747cb94ee6SHong Zhang       it computed. (Defaults to true).
8757cb94ee6SHong Zhang 
8767cb94ee6SHong Zhang    Input Parameter:
8777cb94ee6SHong Zhang +   ts - the time-step context
8787cb94ee6SHong Zhang -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
8797cb94ee6SHong Zhang 
8807cb94ee6SHong Zhang    Level: beginner
8817cb94ee6SHong Zhang 
8827cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8837cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8847cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8857cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
8867cb94ee6SHong Zhang @*/
8877087cfbeSBarry Smith PetscErrorCode  TSSundialsSetExactFinalTime(TS ts,PetscBool  ft)
8887cb94ee6SHong Zhang {
8894ac538c5SBarry Smith   PetscErrorCode ierr;
8907cb94ee6SHong Zhang 
8917cb94ee6SHong Zhang   PetscFunctionBegin;
8924ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetExactFinalTime_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
8937cb94ee6SHong Zhang   PetscFunctionReturn(0);
8947cb94ee6SHong Zhang }
8957cb94ee6SHong Zhang 
8962bfc04deSHong Zhang #undef __FUNCT__
897f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep"
898f1cd61daSJed Brown /*@
899f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
900f1cd61daSJed Brown 
901f1cd61daSJed Brown    Input Parameter:
902f1cd61daSJed Brown +   ts - the time-step context
903f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
904f1cd61daSJed Brown 
905fc6b9e64SJed Brown    Note:
906fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
907fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
908fc6b9e64SJed Brown 
909f1cd61daSJed Brown    Level: beginner
910f1cd61daSJed Brown 
911f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
912f1cd61daSJed Brown @*/
9137087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
914f1cd61daSJed Brown {
9154ac538c5SBarry Smith   PetscErrorCode ierr;
916f1cd61daSJed Brown 
917f1cd61daSJed Brown   PetscFunctionBegin;
9184ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
919f1cd61daSJed Brown   PetscFunctionReturn(0);
920f1cd61daSJed Brown }
921f1cd61daSJed Brown 
922f1cd61daSJed Brown #undef __FUNCT__
923f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep"
924f1cd61daSJed Brown /*@
925f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
926f1cd61daSJed Brown 
927f1cd61daSJed Brown    Input Parameter:
928f1cd61daSJed Brown +   ts - the time-step context
929f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
930f1cd61daSJed Brown 
931f1cd61daSJed Brown    Level: beginner
932f1cd61daSJed Brown 
933f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
934f1cd61daSJed Brown @*/
9357087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
936f1cd61daSJed Brown {
9374ac538c5SBarry Smith   PetscErrorCode ierr;
938f1cd61daSJed Brown 
939f1cd61daSJed Brown   PetscFunctionBegin;
9404ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
941f1cd61daSJed Brown   PetscFunctionReturn(0);
942f1cd61daSJed Brown }
943f1cd61daSJed Brown 
944f1cd61daSJed Brown #undef __FUNCT__
9452bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps"
9462bfc04deSHong Zhang /*@
9472bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
9482bfc04deSHong Zhang 
9492bfc04deSHong Zhang    Input Parameter:
9502bfc04deSHong Zhang +   ts - the time-step context
9512bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
9522bfc04deSHong Zhang 
9532bfc04deSHong Zhang    Level: beginner
9542bfc04deSHong Zhang 
9552bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9562bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
9572bfc04deSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9582bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
9592bfc04deSHong Zhang @*/
9607087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
9612bfc04deSHong Zhang {
9624ac538c5SBarry Smith   PetscErrorCode ierr;
9632bfc04deSHong Zhang 
9642bfc04deSHong Zhang   PetscFunctionBegin;
9654ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
9662bfc04deSHong Zhang   PetscFunctionReturn(0);
9672bfc04deSHong Zhang }
9687cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
9697cb94ee6SHong Zhang /*MC
97096f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
9717cb94ee6SHong Zhang 
9727cb94ee6SHong Zhang    Options Database:
9737cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
9747cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
9757cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
9767cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
9777cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
9787cb94ee6SHong Zhang .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
97916016685SHong Zhang .    -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time
98016016685SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
98116016685SHong Zhang 
9827cb94ee6SHong Zhang 
9837cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
9847cb94ee6SHong Zhang            only PETSc PC options
9857cb94ee6SHong Zhang 
9867cb94ee6SHong Zhang     Level: beginner
9877cb94ee6SHong Zhang 
9887cb94ee6SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
9897cb94ee6SHong Zhang            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
9907cb94ee6SHong Zhang 
9917cb94ee6SHong Zhang M*/
9927cb94ee6SHong Zhang EXTERN_C_BEGIN
9937cb94ee6SHong Zhang #undef __FUNCT__
9947cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
9957087cfbeSBarry Smith PetscErrorCode  TSCreate_Sundials(TS ts)
9967cb94ee6SHong Zhang {
9977cb94ee6SHong Zhang   TS_Sundials *cvode;
9987cb94ee6SHong Zhang   PetscErrorCode ierr;
9997cb94ee6SHong Zhang 
10007cb94ee6SHong Zhang   PetscFunctionBegin;
100117186662SBarry Smith   if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support for nonlinear problems");
1002277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
100328adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
100428adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
10057cb94ee6SHong Zhang   ts->ops->setup          = TSSetUp_Sundials_Nonlinear;
10067cb94ee6SHong Zhang   ts->ops->step           = TSStep_Sundials_Nonlinear;
10077cb94ee6SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear;
10087cb94ee6SHong Zhang 
100938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
10107adad957SLisandro Dalcin   ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr);
10117cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr);
10127cb94ee6SHong Zhang   ts->data          = (void*)cvode;
10136fadb2cdSHong Zhang   cvode->cvode_type = SUNDIALS_BDF;
10146fadb2cdSHong Zhang   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
10157cb94ee6SHong Zhang   cvode->restart    = 5;
10167cb94ee6SHong Zhang   cvode->linear_tol = .05;
10177cb94ee6SHong Zhang 
10182bfc04deSHong Zhang   cvode->exact_final_time = PETSC_TRUE;
10192bfc04deSHong Zhang   cvode->monitorstep      = PETSC_FALSE;
10207cb94ee6SHong Zhang 
1021a07356b8SSatish Balay   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
1022f1cd61daSJed Brown 
1023f1cd61daSJed Brown   cvode->mindt = -1.;
1024f1cd61daSJed Brown   cvode->maxdt = -1.;
1025f1cd61daSJed Brown 
10267cb94ee6SHong Zhang   /* set tolerance for Sundials */
10277cb94ee6SHong Zhang   cvode->reltol = 1e-6;
10282c823083SHong Zhang   cvode->abstol = 1e-6;
10297cb94ee6SHong Zhang 
10307cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
10317cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
10327cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
10337cb94ee6SHong Zhang                     "TSSundialsSetGMRESRestart_Sundials",
10347cb94ee6SHong Zhang                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
10357cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
10367cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
10377cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
10387cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
10397cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
10407cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
10417cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
10427cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
10437cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1044f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1045f1cd61daSJed Brown                     "TSSundialsSetMinTimeStep_Sundials",
1046f1cd61daSJed Brown                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1047f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1048f1cd61daSJed Brown                     "TSSundialsSetMaxTimeStep_Sundials",
1049f1cd61daSJed Brown                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
10507cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
10517cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
10527cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
10537cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
10547cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
10557cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
10567cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
10577cb94ee6SHong Zhang                     "TSSundialsSetExactFinalTime_Sundials",
10587cb94ee6SHong Zhang                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
10592bfc04deSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
10602bfc04deSHong Zhang                     "TSSundialsMonitorInternalSteps_Sundials",
10612bfc04deSHong Zhang                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
10622bfc04deSHong Zhang 
10637cb94ee6SHong Zhang   PetscFunctionReturn(0);
10647cb94ee6SHong Zhang }
10657cb94ee6SHong Zhang EXTERN_C_END
1066