xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 5ff03cf70564ed4f01f7bb4af843d0d7b4679755)
17cb94ee6SHong Zhang #define PETSCTS_DLL
27cb94ee6SHong Zhang 
37cb94ee6SHong Zhang /*
47cb94ee6SHong Zhang     Provides a PETSc interface to SUNDIALS. Alan Hindmarsh's parallel ODE
57cb94ee6SHong Zhang     solver.
67cb94ee6SHong Zhang     The interface to PVODE (old version of CVODE) was originally contributed
77cb94ee6SHong Zhang     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.
87cb94ee6SHong Zhang */
97cb94ee6SHong Zhang 
107cb94ee6SHong Zhang #include "src/ts/impls/implicit/sundials/sundials.h"  /*I "petscts.h" I*/
117cb94ee6SHong Zhang 
127cb94ee6SHong Zhang /*
137cb94ee6SHong Zhang       TSPrecond_Sundials - function that we provide to SUNDIALS to
147cb94ee6SHong Zhang                         evaluate the preconditioner.
157cb94ee6SHong Zhang */
167cb94ee6SHong Zhang #undef __FUNCT__
177cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials"
187cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
197cb94ee6SHong Zhang                     booleantype jok,booleantype *jcurPtr,
207cb94ee6SHong Zhang                     realtype _gamma,void *P_data,
217cb94ee6SHong Zhang                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
227cb94ee6SHong Zhang {
237cb94ee6SHong Zhang   TS             ts = (TS) P_data;
247cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
257cb94ee6SHong Zhang   PC             pc = cvode->pc;
267cb94ee6SHong Zhang   PetscErrorCode ierr;
277cb94ee6SHong Zhang   Mat            Jac = ts->B;
287cb94ee6SHong Zhang   Vec            yy = cvode->w1;
297cb94ee6SHong Zhang   PetscScalar    one = 1.0,gm;
307cb94ee6SHong Zhang   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
317cb94ee6SHong Zhang   PetscScalar    *y_data;
327cb94ee6SHong Zhang 
337cb94ee6SHong Zhang   PetscFunctionBegin;
347cb94ee6SHong Zhang   /* This allows us to construct preconditioners in-place if we like */
357cb94ee6SHong Zhang   ierr = MatSetUnfactored(Jac);CHKERRQ(ierr);
367cb94ee6SHong Zhang 
377cb94ee6SHong Zhang   /* jok - TRUE means reuse current Jacobian else recompute Jacobian */
387cb94ee6SHong Zhang   if (jok) {
397cb94ee6SHong Zhang     ierr     = MatCopy(cvode->pmat,Jac,str);CHKERRQ(ierr);
407cb94ee6SHong Zhang     *jcurPtr = FALSE;
417cb94ee6SHong Zhang   } else {
427cb94ee6SHong Zhang     /* make PETSc vector yy point to SUNDIALS vector y */
437cb94ee6SHong Zhang     y_data = (PetscScalar *) N_VGetArrayPointer(y);
447cb94ee6SHong Zhang     ierr   = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
457cb94ee6SHong Zhang     /* compute the Jacobian */
467cb94ee6SHong Zhang     ierr = TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);CHKERRQ(ierr);
477cb94ee6SHong Zhang     ierr = VecResetArray(yy); CHKERRQ(ierr);
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);
527cb94ee6SHong Zhang     }
537cb94ee6SHong Zhang     else {
547cb94ee6SHong Zhang       ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr);
557cb94ee6SHong Zhang     }
567cb94ee6SHong Zhang     *jcurPtr = TRUE;
577cb94ee6SHong Zhang   }
587cb94ee6SHong Zhang 
597cb94ee6SHong Zhang   /* construct I-gamma*Jac  */
607cb94ee6SHong Zhang   gm   = -_gamma;
617cb94ee6SHong Zhang   ierr = MatScale(Jac,gm);CHKERRQ(ierr);
627cb94ee6SHong Zhang   ierr = MatShift(Jac,one);CHKERRQ(ierr);
637cb94ee6SHong Zhang 
647cb94ee6SHong Zhang   ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr);
657cb94ee6SHong Zhang   PetscFunctionReturn(0);
667cb94ee6SHong Zhang }
677cb94ee6SHong Zhang 
687cb94ee6SHong Zhang /*
697cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
707cb94ee6SHong Zhang */
717cb94ee6SHong Zhang #undef __FUNCT__
727cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials"
737cb94ee6SHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,
747cb94ee6SHong Zhang                                  N_Vector r,N_Vector z,
757cb94ee6SHong Zhang                                  realtype _gamma,realtype delta,
767cb94ee6SHong Zhang                                  int lr,void *P_data,N_Vector vtemp)
777cb94ee6SHong Zhang {
787cb94ee6SHong Zhang   TS              ts = (TS) P_data;
797cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
807cb94ee6SHong Zhang   PC              pc = cvode->pc;
817cb94ee6SHong Zhang   Vec             rr = cvode->w1,zz = cvode->w2;
827cb94ee6SHong Zhang   PetscErrorCode  ierr;
837cb94ee6SHong Zhang   PetscScalar     *r_data,*z_data;
847cb94ee6SHong Zhang 
857cb94ee6SHong Zhang   PetscFunctionBegin;
867cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
877cb94ee6SHong Zhang   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
887cb94ee6SHong Zhang   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
897cb94ee6SHong Zhang   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
907cb94ee6SHong Zhang   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
917cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
927cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
937cb94ee6SHong Zhang   ierr = VecResetArray(rr); CHKERRQ(ierr);
947cb94ee6SHong Zhang   ierr = VecResetArray(zz); CHKERRQ(ierr);
957cb94ee6SHong Zhang   cvode->linear_solves++;
967cb94ee6SHong Zhang   PetscFunctionReturn(0);
977cb94ee6SHong Zhang }
987cb94ee6SHong Zhang 
997cb94ee6SHong Zhang /*
1007cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
1017cb94ee6SHong Zhang */
1027cb94ee6SHong Zhang #undef __FUNCT__
1037cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials"
1047cb94ee6SHong Zhang void TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
1057cb94ee6SHong Zhang {
1067cb94ee6SHong Zhang   TS              ts = (TS) ctx;
1077cb94ee6SHong Zhang   MPI_Comm        comm = ts->comm;
1087cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
1097cb94ee6SHong Zhang   Vec             yy = cvode->w1,yyd = cvode->w2;
1107cb94ee6SHong Zhang   PetscScalar     *y_data,*ydot_data;
1117cb94ee6SHong Zhang   PetscErrorCode  ierr;
1127cb94ee6SHong Zhang 
1137cb94ee6SHong Zhang   PetscFunctionBegin;
114*5ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
1157cb94ee6SHong Zhang   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
1167cb94ee6SHong Zhang   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
1177cb94ee6SHong Zhang   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr)
1187cb94ee6SHong Zhang   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr)
1197cb94ee6SHong Zhang   /* now compute the right hand side function */
1207cb94ee6SHong Zhang   ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr);
1217cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
1227cb94ee6SHong Zhang   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
1237cb94ee6SHong Zhang   PetscFunctionReturnVoid();
1247cb94ee6SHong Zhang }
1257cb94ee6SHong Zhang 
1267cb94ee6SHong Zhang /*
1277cb94ee6SHong Zhang        TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE.
1287cb94ee6SHong Zhang */
1297cb94ee6SHong Zhang #undef __FUNCT__
1307cb94ee6SHong Zhang #define __FUNCT__ "TSStep_Sundials_Nonlinear"
1317cb94ee6SHong Zhang /*
1327cb94ee6SHong Zhang     TSStep_Sundials_Nonlinear -
1337cb94ee6SHong Zhang 
1347cb94ee6SHong Zhang    steps - number of time steps
1357cb94ee6SHong Zhang    time - time that integrater is  terminated.
1367cb94ee6SHong Zhang */
1377cb94ee6SHong Zhang PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
1387cb94ee6SHong Zhang {
1397cb94ee6SHong Zhang   TS_Sundials  *cvode = (TS_Sundials*)ts->data;
1407cb94ee6SHong Zhang   Vec          sol = ts->vec_sol;
1417cb94ee6SHong Zhang   PetscErrorCode ierr;
1427cb94ee6SHong Zhang   int          i,max_steps = ts->max_steps,flag;
1437cb94ee6SHong Zhang   long int     its;
1447cb94ee6SHong Zhang   realtype     t,tout;
1457cb94ee6SHong Zhang   PetscScalar  *y_data;
1467cb94ee6SHong Zhang   void         *mem;
1477cb94ee6SHong Zhang 
1487cb94ee6SHong Zhang   PetscFunctionBegin;
1497cb94ee6SHong Zhang   /*
1507cb94ee6SHong Zhang      Call CVodeCreate to create the solver memory:
1517cb94ee6SHong Zhang      CV_ADAMS   specifies the Adams Method
1527cb94ee6SHong Zhang      CV_FUNCTIONAL  specifies functional iteration
1537cb94ee6SHong Zhang      A pointer to the integrator memory is returned and stored in cvode_mem.
1547cb94ee6SHong Zhang   */
1557cb94ee6SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
1567cb94ee6SHong Zhang   if (!mem) SETERRQ(1,"CVodeCreate() fails");
1577cb94ee6SHong Zhang   flag = CVodeSetFdata(mem,ts);
1587cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVodeSetFdata() fails");
1597cb94ee6SHong Zhang 
1607cb94ee6SHong Zhang   /*
1617cb94ee6SHong Zhang      Call CVodeMalloc to initialize the integrator memory:
1627cb94ee6SHong Zhang      mem is the pointer to the integrator memory returned by CVodeCreate
1637cb94ee6SHong Zhang      f       is the user's right hand side function in y'=f(t,y)
1647cb94ee6SHong Zhang      T0      is the initial time
1657cb94ee6SHong Zhang      u       is the initial dependent variable vector
1667cb94ee6SHong Zhang      CV_SS   specifies scalar relative and absolute tolerances
1677cb94ee6SHong Zhang      reltol  is the relative tolerance
1687cb94ee6SHong Zhang      &abstol is a pointer to the scalar absolute tolerance
1697cb94ee6SHong Zhang   */
1707cb94ee6SHong Zhang   flag = CVodeMalloc(mem,TSFunction_Sundials,ts->ptime,cvode->y,CV_SS,cvode->reltol,&cvode->abstol);
1717cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVodeMalloc() fails");
1727cb94ee6SHong Zhang 
1737cb94ee6SHong Zhang   /* initialize the number of steps */
1747cb94ee6SHong Zhang   *steps = -ts->steps;
1757cb94ee6SHong Zhang   ierr   = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1767cb94ee6SHong Zhang 
1777cb94ee6SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver. */
1787cb94ee6SHong Zhang   /* setup the ode integrator with the given preconditioner */
1797cb94ee6SHong Zhang   /* flag  = CVSpgmr(mem,PREC_LEFT,cvode->restart); */
1807cb94ee6SHong Zhang   flag  = CVSpgmr(mem,PREC_LEFT,0);
1817cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVSpgmr() fails");
1827cb94ee6SHong Zhang   flag = CVSpgmrSetGSType(mem,MODIFIED_GS);
1837cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVSpgmrSetGSType() fails");
1847cb94ee6SHong Zhang 
1857cb94ee6SHong Zhang   /* Set preconditioner setup and solve routines Precond and PSolve,
1867cb94ee6SHong Zhang      and the pointer to the user-defined block data */
1877cb94ee6SHong Zhang   flag = CVSpgmrSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials,ts);
1887cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVSpgmrSetPreconditioner() fails");
1897cb94ee6SHong Zhang 
1907cb94ee6SHong Zhang   tout = ts->max_time;
1917cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1927cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data,cvode->y);
1937cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
1947cb94ee6SHong Zhang   for (i = 0; i < max_steps; i++) {
1957cb94ee6SHong Zhang     if (ts->ptime >= ts->max_time) break;
1967cb94ee6SHong Zhang     ierr = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);CHKERRQ(ierr);
1977cb94ee6SHong Zhang     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
1987cb94ee6SHong Zhang     cvode->nonlinear_solves += its;
1997cb94ee6SHong Zhang 
2007cb94ee6SHong Zhang     if (t > ts->max_time && cvode->exact_final_time) {
2017cb94ee6SHong Zhang       /* interpolate to final requested time */
2027cb94ee6SHong Zhang       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
2037cb94ee6SHong Zhang       t = tout;
2047cb94ee6SHong Zhang     }
2057cb94ee6SHong Zhang     ts->time_step = t - ts->ptime;
2067cb94ee6SHong Zhang     ts->ptime     = t;
2077cb94ee6SHong Zhang 
2087cb94ee6SHong Zhang     /* copy the solution from cvode->y to cvode->update and sol */
2097cb94ee6SHong Zhang     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
2107cb94ee6SHong Zhang     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
2117cb94ee6SHong Zhang     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
2127cb94ee6SHong Zhang     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
2137cb94ee6SHong Zhang     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
2147cb94ee6SHong Zhang     ts->nonlinear_its = its;
2157cb94ee6SHong Zhang     ierr = CVSpgmrGetNumLinIters(mem, &its);
2167cb94ee6SHong Zhang     ts->linear_its = its;
2177cb94ee6SHong Zhang     ts->steps++;
2187cb94ee6SHong Zhang     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
2197cb94ee6SHong Zhang   }
2207cb94ee6SHong Zhang   CVodeFree(mem);
2217cb94ee6SHong Zhang   *steps += ts->steps;
2227cb94ee6SHong Zhang   *time   = t;
2237cb94ee6SHong Zhang   PetscFunctionReturn(0);
2247cb94ee6SHong Zhang }
2257cb94ee6SHong Zhang 
2267cb94ee6SHong Zhang #undef __FUNCT__
2277cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials"
2287cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
2297cb94ee6SHong Zhang {
2307cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2317cb94ee6SHong Zhang   PetscErrorCode ierr;
2327cb94ee6SHong Zhang 
2337cb94ee6SHong Zhang   PetscFunctionBegin;
2347cb94ee6SHong Zhang   if (cvode->pmat)   {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);}
2357cb94ee6SHong Zhang   if (cvode->pc)     {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);}
2367cb94ee6SHong Zhang   if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);}
2377cb94ee6SHong Zhang   if (cvode->func)   {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);}
2387cb94ee6SHong Zhang   if (cvode->rhs)    {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);}
2397cb94ee6SHong Zhang   if (cvode->w1)     {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);}
2407cb94ee6SHong Zhang   if (cvode->w2)     {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);}
2417cb94ee6SHong Zhang   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
2427cb94ee6SHong Zhang   ierr = PetscFree(cvode);CHKERRQ(ierr);
2437cb94ee6SHong Zhang   PetscFunctionReturn(0);
2447cb94ee6SHong Zhang }
2457cb94ee6SHong Zhang 
2467cb94ee6SHong Zhang #undef __FUNCT__
2477cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear"
2487cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
2497cb94ee6SHong Zhang {
2507cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2517cb94ee6SHong Zhang   PetscErrorCode ierr;
2527cb94ee6SHong Zhang   int            glosize,locsize,i;
2537cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
2547cb94ee6SHong Zhang 
2557cb94ee6SHong Zhang   PetscFunctionBegin;
2567cb94ee6SHong Zhang   ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr);
2577cb94ee6SHong Zhang   /* get the vector size */
2587cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2597cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
2607cb94ee6SHong Zhang 
2617cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
2627cb94ee6SHong Zhang   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
2637cb94ee6SHong Zhang   if (!cvode->y) SETERRQ(1,"cvode->y is not allocated");
2647cb94ee6SHong Zhang 
2657cb94ee6SHong Zhang   /* initialize N_Vec y */
2667cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
2677cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
2687cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
2697cb94ee6SHong Zhang   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
2707cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
2717cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
2727cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr);
2737cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
2747cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr);
2757cb94ee6SHong Zhang 
2767cb94ee6SHong Zhang   /*
2777cb94ee6SHong Zhang       Create work vectors for the TSPSolve_Sundials() routine. Note these are
2787cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
2797cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
2807cb94ee6SHong Zhang   */
2817cb94ee6SHong Zhang   ierr = VecCreateMPIWithArray(ts->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
2827cb94ee6SHong Zhang   ierr = VecCreateMPIWithArray(ts->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
2837cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
2847cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
2857cb94ee6SHong Zhang   PetscFunctionReturn(0);
2867cb94ee6SHong Zhang }
2877cb94ee6SHong Zhang 
2887cb94ee6SHong Zhang #undef __FUNCT__
2897cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear"
2907cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
2917cb94ee6SHong Zhang {
2927cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2937cb94ee6SHong Zhang   PetscErrorCode ierr;
2947cb94ee6SHong Zhang   int            indx;
2957cb94ee6SHong Zhang   const char     *btype[] = {"bdf","adams"},*otype[] = {"modified","unmodified"};
2967cb94ee6SHong Zhang   PetscTruth     flag;
2977cb94ee6SHong Zhang 
2987cb94ee6SHong Zhang   PetscFunctionBegin;
2997cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
3007cb94ee6SHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",btype,2,"bdf",&indx,&flag);CHKERRQ(ierr);
3017cb94ee6SHong Zhang     if (flag) {
3027cb94ee6SHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsType)indx);CHKERRQ(ierr);
3037cb94ee6SHong Zhang     }
3047cb94ee6SHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",otype,2,"unmodified",&indx,&flag);CHKERRQ(ierr);
3057cb94ee6SHong Zhang     if (flag) {
3067cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
3077cb94ee6SHong Zhang     }
3087cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
3097cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
3107cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
3117cb94ee6SHong Zhang     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
3127cb94ee6SHong Zhang     ierr = PetscOptionsName("-ts_sundials_exact_final_time","Allow SUNDIALS to stop near the final time, not exactly on it","TSSundialsSetExactFinalTime",&cvode->exact_final_time);CHKERRQ(ierr);
3137cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3147cb94ee6SHong Zhang   PetscFunctionReturn(0);
3157cb94ee6SHong Zhang }
3167cb94ee6SHong Zhang 
3177cb94ee6SHong Zhang #undef __FUNCT__
3187cb94ee6SHong Zhang #define __FUNCT__ "TSPrintHelp_Sundials"
3197cb94ee6SHong Zhang PetscErrorCode TSPrintHelp_Sundials(TS ts,char *p)
3207cb94ee6SHong Zhang {
3217cb94ee6SHong Zhang   PetscErrorCode ierr;
3227cb94ee6SHong Zhang 
3237cb94ee6SHong Zhang   PetscFunctionBegin;
3247cb94ee6SHong Zhang   ierr = (*PetscHelpPrintf)(ts->comm," Options for TSSUNDIALS integrater:\n");CHKERRQ(ierr);
3257cb94ee6SHong Zhang   ierr = (*PetscHelpPrintf)(ts->comm," -ts_sundials_type <bdf,adams>: integration approach\n",p);CHKERRQ(ierr);
3267cb94ee6SHong Zhang   ierr = (*PetscHelpPrintf)(ts->comm," -ts_sundials_atol aabs: absolute tolerance of ODE solution\n",p);CHKERRQ(ierr);
3277cb94ee6SHong Zhang   ierr = (*PetscHelpPrintf)(ts->comm," -ts_sundials_rtol rel: relative tolerance of ODE solution\n",p);CHKERRQ(ierr);
3287cb94ee6SHong Zhang   ierr = (*PetscHelpPrintf)(ts->comm," -ts_sundials_gramschmidt_type <unmodified,modified>\n");CHKERRQ(ierr);
3297cb94ee6SHong Zhang   ierr = (*PetscHelpPrintf)(ts->comm," -ts_sundials_gmres_restart <restart_size> (also max. GMRES its)\n");CHKERRQ(ierr);
3307cb94ee6SHong Zhang   ierr = (*PetscHelpPrintf)(ts->comm," -ts_sundials_linear_tolerance <tol>\n");CHKERRQ(ierr);
3317cb94ee6SHong Zhang   ierr = (*PetscHelpPrintf)(ts->comm," -ts_sundials_not_exact_final_time\n");CHKERRQ(ierr);
3327cb94ee6SHong Zhang   PetscFunctionReturn(0);
3337cb94ee6SHong Zhang }
3347cb94ee6SHong Zhang 
3357cb94ee6SHong Zhang #undef __FUNCT__
3367cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
3377cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
3387cb94ee6SHong Zhang {
3397cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3407cb94ee6SHong Zhang   PetscErrorCode ierr;
3417cb94ee6SHong Zhang   char           *type;
3427cb94ee6SHong Zhang   char            atype[] = "Adams";
3437cb94ee6SHong Zhang   char            btype[] = "BDF: backward differentiation formula";
3447cb94ee6SHong Zhang   PetscTruth     iascii,isstring;
3457cb94ee6SHong Zhang 
3467cb94ee6SHong Zhang   PetscFunctionBegin;
3477cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
3487cb94ee6SHong Zhang   else                                     {type = btype;}
3497cb94ee6SHong Zhang 
3507cb94ee6SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
3517cb94ee6SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
3527cb94ee6SHong Zhang   if (iascii) {
3537cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
3547cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
3557cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
3567cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
3577cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
3587cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
3597cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
3607cb94ee6SHong Zhang     } else {
3617cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
3627cb94ee6SHong Zhang     }
3637cb94ee6SHong Zhang   } else if (isstring) {
3647cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
3657cb94ee6SHong Zhang   } else {
3667cb94ee6SHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
3677cb94ee6SHong Zhang   }
3687cb94ee6SHong Zhang   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3697cb94ee6SHong Zhang   ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr);
3707cb94ee6SHong Zhang   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3717cb94ee6SHong Zhang   PetscFunctionReturn(0);
3727cb94ee6SHong Zhang }
3737cb94ee6SHong Zhang 
3747cb94ee6SHong Zhang 
3757cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
3767cb94ee6SHong Zhang EXTERN_C_BEGIN
3777cb94ee6SHong Zhang #undef __FUNCT__
3787cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
3797cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsType type)
3807cb94ee6SHong Zhang {
3817cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
3827cb94ee6SHong Zhang 
3837cb94ee6SHong Zhang   PetscFunctionBegin;
3847cb94ee6SHong Zhang   cvode->cvode_type = type;
3857cb94ee6SHong Zhang   PetscFunctionReturn(0);
3867cb94ee6SHong Zhang }
3877cb94ee6SHong Zhang EXTERN_C_END
3887cb94ee6SHong Zhang 
3897cb94ee6SHong Zhang EXTERN_C_BEGIN
3907cb94ee6SHong Zhang #undef __FUNCT__
3917cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
3927cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
3937cb94ee6SHong Zhang {
3947cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
3957cb94ee6SHong Zhang 
3967cb94ee6SHong Zhang   PetscFunctionBegin;
3977cb94ee6SHong Zhang   cvode->restart = restart;
3987cb94ee6SHong Zhang   PetscFunctionReturn(0);
3997cb94ee6SHong Zhang }
4007cb94ee6SHong Zhang EXTERN_C_END
4017cb94ee6SHong Zhang 
4027cb94ee6SHong Zhang EXTERN_C_BEGIN
4037cb94ee6SHong Zhang #undef __FUNCT__
4047cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
4057cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
4067cb94ee6SHong Zhang {
4077cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4087cb94ee6SHong Zhang 
4097cb94ee6SHong Zhang   PetscFunctionBegin;
4107cb94ee6SHong Zhang   cvode->linear_tol = tol;
4117cb94ee6SHong Zhang   PetscFunctionReturn(0);
4127cb94ee6SHong Zhang }
4137cb94ee6SHong Zhang EXTERN_C_END
4147cb94ee6SHong Zhang 
4157cb94ee6SHong Zhang EXTERN_C_BEGIN
4167cb94ee6SHong Zhang #undef __FUNCT__
4177cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
4187cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
4197cb94ee6SHong Zhang {
4207cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4217cb94ee6SHong Zhang 
4227cb94ee6SHong Zhang   PetscFunctionBegin;
4237cb94ee6SHong Zhang   cvode->gtype = type;
4247cb94ee6SHong Zhang   PetscFunctionReturn(0);
4257cb94ee6SHong Zhang }
4267cb94ee6SHong Zhang EXTERN_C_END
4277cb94ee6SHong Zhang 
4287cb94ee6SHong Zhang EXTERN_C_BEGIN
4297cb94ee6SHong Zhang #undef __FUNCT__
4307cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
4317cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
4327cb94ee6SHong Zhang {
4337cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4347cb94ee6SHong Zhang 
4357cb94ee6SHong Zhang   PetscFunctionBegin;
4367cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
4377cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
4387cb94ee6SHong Zhang   PetscFunctionReturn(0);
4397cb94ee6SHong Zhang }
4407cb94ee6SHong Zhang EXTERN_C_END
4417cb94ee6SHong Zhang 
4427cb94ee6SHong Zhang EXTERN_C_BEGIN
4437cb94ee6SHong Zhang #undef __FUNCT__
4447cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
4457cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc)
4467cb94ee6SHong Zhang {
4477cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4487cb94ee6SHong Zhang 
4497cb94ee6SHong Zhang   PetscFunctionBegin;
4507cb94ee6SHong Zhang   *pc = cvode->pc;
4517cb94ee6SHong Zhang   PetscFunctionReturn(0);
4527cb94ee6SHong Zhang }
4537cb94ee6SHong Zhang EXTERN_C_END
4547cb94ee6SHong Zhang 
4557cb94ee6SHong Zhang EXTERN_C_BEGIN
4567cb94ee6SHong Zhang #undef __FUNCT__
4577cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
4587cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
4597cb94ee6SHong Zhang {
4607cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4617cb94ee6SHong Zhang 
4627cb94ee6SHong Zhang   PetscFunctionBegin;
4637cb94ee6SHong Zhang   if (nonlin) *nonlin = cvode->nonlinear_solves;
4647cb94ee6SHong Zhang   if (lin)    *lin    = cvode->linear_solves;
4657cb94ee6SHong Zhang   PetscFunctionReturn(0);
4667cb94ee6SHong Zhang }
4677cb94ee6SHong Zhang EXTERN_C_END
4687cb94ee6SHong Zhang 
4697cb94ee6SHong Zhang EXTERN_C_BEGIN
4707cb94ee6SHong Zhang #undef __FUNCT__
4717cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
4727cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
4737cb94ee6SHong Zhang {
4747cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4757cb94ee6SHong Zhang 
4767cb94ee6SHong Zhang   PetscFunctionBegin;
4777cb94ee6SHong Zhang   cvode->exact_final_time = s;
4787cb94ee6SHong Zhang   PetscFunctionReturn(0);
4797cb94ee6SHong Zhang }
4807cb94ee6SHong Zhang EXTERN_C_END
4817cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
4827cb94ee6SHong Zhang 
4837cb94ee6SHong Zhang #undef __FUNCT__
4847cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
4857cb94ee6SHong Zhang /*@C
4867cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
4877cb94ee6SHong Zhang 
4887cb94ee6SHong Zhang    Not Collective
4897cb94ee6SHong Zhang 
4907cb94ee6SHong Zhang    Input parameters:
4917cb94ee6SHong Zhang .    ts     - the time-step context
4927cb94ee6SHong Zhang 
4937cb94ee6SHong Zhang    Output Parameters:
4947cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
4957cb94ee6SHong Zhang -   lin    - number of linear iterations
4967cb94ee6SHong Zhang 
4977cb94ee6SHong Zhang    Level: advanced
4987cb94ee6SHong Zhang 
4997cb94ee6SHong Zhang    Notes:
5007cb94ee6SHong Zhang     These return the number since the creation of the TS object
5017cb94ee6SHong Zhang 
5027cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
5037cb94ee6SHong Zhang 
5047cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5057cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
5067cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5077cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
5087cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
5097cb94ee6SHong Zhang 
5107cb94ee6SHong Zhang @*/
5117cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
5127cb94ee6SHong Zhang {
5137cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,int*,int*);
5147cb94ee6SHong Zhang 
5157cb94ee6SHong Zhang   PetscFunctionBegin;
5167cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
5177cb94ee6SHong Zhang   if (f) {
5187cb94ee6SHong Zhang     ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr);
5197cb94ee6SHong Zhang   }
5207cb94ee6SHong Zhang   PetscFunctionReturn(0);
5217cb94ee6SHong Zhang }
5227cb94ee6SHong Zhang 
5237cb94ee6SHong Zhang #undef __FUNCT__
5247cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
5257cb94ee6SHong Zhang /*@
5267cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
5277cb94ee6SHong Zhang 
5287cb94ee6SHong Zhang    Collective on TS
5297cb94ee6SHong Zhang 
5307cb94ee6SHong Zhang    Input parameters:
5317cb94ee6SHong Zhang +    ts     - the time-step context
5327cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
5337cb94ee6SHong Zhang 
5347cb94ee6SHong Zhang    Level: intermediate
5357cb94ee6SHong Zhang 
5367cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
5377cb94ee6SHong Zhang 
5387cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
5397cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
5407cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5417cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
5427cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
5437cb94ee6SHong Zhang @*/
5447cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsType type)
5457cb94ee6SHong Zhang {
5467cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,TSSundialsType);
5477cb94ee6SHong Zhang 
5487cb94ee6SHong Zhang   PetscFunctionBegin;
5497cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
5507cb94ee6SHong Zhang   if (f) {
5517cb94ee6SHong Zhang     ierr = (*f)(ts,type);CHKERRQ(ierr);
5527cb94ee6SHong Zhang   }
5537cb94ee6SHong Zhang   PetscFunctionReturn(0);
5547cb94ee6SHong Zhang }
5557cb94ee6SHong Zhang 
5567cb94ee6SHong Zhang #undef __FUNCT__
5577cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart"
5587cb94ee6SHong Zhang /*@
5597cb94ee6SHong Zhang    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
5607cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
5617cb94ee6SHong Zhang        this is ALSO the maximum number of GMRES steps that will be used.
5627cb94ee6SHong Zhang 
5637cb94ee6SHong Zhang    Collective on TS
5647cb94ee6SHong Zhang 
5657cb94ee6SHong Zhang    Input parameters:
5667cb94ee6SHong Zhang +    ts      - the time-step context
5677cb94ee6SHong Zhang -    restart - number of direction vectors (the restart size).
5687cb94ee6SHong Zhang 
5697cb94ee6SHong Zhang    Level: advanced
5707cb94ee6SHong Zhang 
5717cb94ee6SHong Zhang .keywords: GMRES, restart
5727cb94ee6SHong Zhang 
5737cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
5747cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
5757cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5767cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
5777cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
5787cb94ee6SHong Zhang 
5797cb94ee6SHong Zhang @*/
5807cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart)
5817cb94ee6SHong Zhang {
5827cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,int);
5837cb94ee6SHong Zhang 
5847cb94ee6SHong Zhang   PetscFunctionBegin;
5857cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr);
5867cb94ee6SHong Zhang   if (f) {
5877cb94ee6SHong Zhang     ierr = (*f)(ts,restart);CHKERRQ(ierr);
5887cb94ee6SHong Zhang   }
5897cb94ee6SHong Zhang   PetscFunctionReturn(0);
5907cb94ee6SHong Zhang }
5917cb94ee6SHong Zhang 
5927cb94ee6SHong Zhang #undef __FUNCT__
5937cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
5947cb94ee6SHong Zhang /*@
5957cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
5967cb94ee6SHong Zhang        system by SUNDIALS.
5977cb94ee6SHong Zhang 
5987cb94ee6SHong Zhang    Collective on TS
5997cb94ee6SHong Zhang 
6007cb94ee6SHong Zhang    Input parameters:
6017cb94ee6SHong Zhang +    ts     - the time-step context
6027cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
6037cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
6047cb94ee6SHong Zhang 
6057cb94ee6SHong Zhang    Level: advanced
6067cb94ee6SHong Zhang 
6077cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
6087cb94ee6SHong Zhang 
6097cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6107cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6117cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6127cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6137cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6147cb94ee6SHong Zhang 
6157cb94ee6SHong Zhang @*/
6167cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol)
6177cb94ee6SHong Zhang {
6187cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,double);
6197cb94ee6SHong Zhang 
6207cb94ee6SHong Zhang   PetscFunctionBegin;
6217cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
6227cb94ee6SHong Zhang   if (f) {
6237cb94ee6SHong Zhang     ierr = (*f)(ts,tol);CHKERRQ(ierr);
6247cb94ee6SHong Zhang   }
6257cb94ee6SHong Zhang   PetscFunctionReturn(0);
6267cb94ee6SHong Zhang }
6277cb94ee6SHong Zhang 
6287cb94ee6SHong Zhang #undef __FUNCT__
6297cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
6307cb94ee6SHong Zhang /*@
6317cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
6327cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
6337cb94ee6SHong Zhang 
6347cb94ee6SHong Zhang    Collective on TS
6357cb94ee6SHong Zhang 
6367cb94ee6SHong Zhang    Input parameters:
6377cb94ee6SHong Zhang +    ts  - the time-step context
6387cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
6397cb94ee6SHong Zhang 
6407cb94ee6SHong Zhang    Level: advanced
6417cb94ee6SHong Zhang 
6427cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
6437cb94ee6SHong Zhang 
6447cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6457cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
6467cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6477cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6487cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6497cb94ee6SHong Zhang 
6507cb94ee6SHong Zhang @*/
6517cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
6527cb94ee6SHong Zhang {
6537cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
6547cb94ee6SHong Zhang 
6557cb94ee6SHong Zhang   PetscFunctionBegin;
6567cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr);
6577cb94ee6SHong Zhang   if (f) {
6587cb94ee6SHong Zhang     ierr = (*f)(ts,type);CHKERRQ(ierr);
6597cb94ee6SHong Zhang   }
6607cb94ee6SHong Zhang   PetscFunctionReturn(0);
6617cb94ee6SHong Zhang }
6627cb94ee6SHong Zhang 
6637cb94ee6SHong Zhang #undef __FUNCT__
6647cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
6657cb94ee6SHong Zhang /*@
6667cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
6677cb94ee6SHong Zhang                          Sundials for error control.
6687cb94ee6SHong Zhang 
6697cb94ee6SHong Zhang    Collective on TS
6707cb94ee6SHong Zhang 
6717cb94ee6SHong Zhang    Input parameters:
6727cb94ee6SHong Zhang +    ts  - the time-step context
6737cb94ee6SHong Zhang .    aabs - the absolute tolerance
6747cb94ee6SHong Zhang -    rel - the relative tolerance
6757cb94ee6SHong Zhang 
6767cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
6777cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
6787cb94ee6SHong Zhang 
6797cb94ee6SHong Zhang    Level: intermediate
6807cb94ee6SHong Zhang 
6817cb94ee6SHong Zhang .keywords: Sundials, tolerance
6827cb94ee6SHong Zhang 
6837cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6847cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
6857cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6867cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6877cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6887cb94ee6SHong Zhang 
6897cb94ee6SHong Zhang @*/
6907cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel)
6917cb94ee6SHong Zhang {
6927cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,double,double);
6937cb94ee6SHong Zhang 
6947cb94ee6SHong Zhang   PetscFunctionBegin;
6957cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
6967cb94ee6SHong Zhang   if (f) {
6977cb94ee6SHong Zhang     ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr);
6987cb94ee6SHong Zhang   }
6997cb94ee6SHong Zhang   PetscFunctionReturn(0);
7007cb94ee6SHong Zhang }
7017cb94ee6SHong Zhang 
7027cb94ee6SHong Zhang #undef __FUNCT__
7037cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
7047cb94ee6SHong Zhang /*@
7057cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
7067cb94ee6SHong Zhang 
7077cb94ee6SHong Zhang    Input Parameter:
7087cb94ee6SHong Zhang .    ts - the time-step context
7097cb94ee6SHong Zhang 
7107cb94ee6SHong Zhang    Output Parameter:
7117cb94ee6SHong Zhang .    pc - the preconditioner context
7127cb94ee6SHong Zhang 
7137cb94ee6SHong Zhang    Level: advanced
7147cb94ee6SHong Zhang 
7157cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7167cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7177cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7187cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
7197cb94ee6SHong Zhang @*/
7207cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc)
7217cb94ee6SHong Zhang {
7227cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,PC *);
7237cb94ee6SHong Zhang 
7247cb94ee6SHong Zhang   PetscFunctionBegin;
7257cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
7267cb94ee6SHong Zhang   if (f) {
7277cb94ee6SHong Zhang     ierr = (*f)(ts,pc);CHKERRQ(ierr);
7287cb94ee6SHong Zhang   } else {
7297cb94ee6SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
7307cb94ee6SHong Zhang   }
7317cb94ee6SHong Zhang 
7327cb94ee6SHong Zhang   PetscFunctionReturn(0);
7337cb94ee6SHong Zhang }
7347cb94ee6SHong Zhang 
7357cb94ee6SHong Zhang #undef __FUNCT__
7367cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime"
7377cb94ee6SHong Zhang /*@
7387cb94ee6SHong Zhang    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
7397cb94ee6SHong Zhang       exact final time requested by the user or just returns it at the final time
7407cb94ee6SHong Zhang       it computed. (Defaults to true).
7417cb94ee6SHong Zhang 
7427cb94ee6SHong Zhang    Input Parameter:
7437cb94ee6SHong Zhang +   ts - the time-step context
7447cb94ee6SHong Zhang -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
7457cb94ee6SHong Zhang 
7467cb94ee6SHong Zhang    Level: beginner
7477cb94ee6SHong Zhang 
7487cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7497cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7507cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7517cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
7527cb94ee6SHong Zhang @*/
7537cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
7547cb94ee6SHong Zhang {
7557cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,PetscTruth);
7567cb94ee6SHong Zhang 
7577cb94ee6SHong Zhang   PetscFunctionBegin;
7587cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr);
7597cb94ee6SHong Zhang   if (f) {
7607cb94ee6SHong Zhang     ierr = (*f)(ts,ft);CHKERRQ(ierr);
7617cb94ee6SHong Zhang   }
7627cb94ee6SHong Zhang   PetscFunctionReturn(0);
7637cb94ee6SHong Zhang }
7647cb94ee6SHong Zhang 
7657cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
7667cb94ee6SHong Zhang /*MC
7677cb94ee6SHong Zhang       TS_Sundials - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
7687cb94ee6SHong Zhang 
7697cb94ee6SHong Zhang    Options Database:
7707cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
7717cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
7727cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
7737cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
7747cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
7757cb94ee6SHong Zhang .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
7767cb94ee6SHong Zhang -    -ts_sundials_not_exact_final_time -Allow SUNDIALS to stop near the final time, not exactly on it
7777cb94ee6SHong Zhang 
7787cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
7797cb94ee6SHong Zhang            only PETSc PC options
7807cb94ee6SHong Zhang 
7817cb94ee6SHong Zhang     Level: beginner
7827cb94ee6SHong Zhang 
7837cb94ee6SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
7847cb94ee6SHong Zhang            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
7857cb94ee6SHong Zhang 
7867cb94ee6SHong Zhang M*/
7877cb94ee6SHong Zhang EXTERN_C_BEGIN
7887cb94ee6SHong Zhang #undef __FUNCT__
7897cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
7907cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts)
7917cb94ee6SHong Zhang {
7927cb94ee6SHong Zhang   TS_Sundials *cvode;
7937cb94ee6SHong Zhang   PetscErrorCode ierr;
7947cb94ee6SHong Zhang 
7957cb94ee6SHong Zhang   PetscFunctionBegin;
7967cb94ee6SHong Zhang   ts->ops->destroy         = TSDestroy_Sundials;
7977cb94ee6SHong Zhang   ts->ops->view            = TSView_Sundials;
7987cb94ee6SHong Zhang 
7997cb94ee6SHong Zhang   if (ts->problem_type != TS_NONLINEAR) {
8007cb94ee6SHong Zhang     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
8017cb94ee6SHong Zhang   }
8027cb94ee6SHong Zhang   ts->ops->setup           = TSSetUp_Sundials_Nonlinear;
8037cb94ee6SHong Zhang   ts->ops->step            = TSStep_Sundials_Nonlinear;
8047cb94ee6SHong Zhang   ts->ops->setfromoptions  = TSSetFromOptions_Sundials_Nonlinear;
8057cb94ee6SHong Zhang 
8067cb94ee6SHong Zhang   ierr = PetscNew(TS_Sundials,&cvode);CHKERRQ(ierr);
8077cb94ee6SHong Zhang   ierr = PCCreate(ts->comm,&cvode->pc);CHKERRQ(ierr);
8087cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr);
8097cb94ee6SHong Zhang   ts->data          = (void*)cvode;
8107cb94ee6SHong Zhang   cvode->cvode_type = CV_BDF;
8117cb94ee6SHong Zhang   cvode->gtype      = SUNDIALS_UNMODIFIED_GS;
8127cb94ee6SHong Zhang   cvode->restart    = 5;
8137cb94ee6SHong Zhang   cvode->linear_tol = .05;
8147cb94ee6SHong Zhang 
8157cb94ee6SHong Zhang   cvode->exact_final_time = PETSC_FALSE;
8167cb94ee6SHong Zhang 
8177cb94ee6SHong Zhang   ierr = MPI_Comm_dup(ts->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
8187cb94ee6SHong Zhang   /* set tolerance for Sundials */
8197cb94ee6SHong Zhang   cvode->abstol = 1e-6;
8207cb94ee6SHong Zhang   cvode->reltol = 1e-6;
8217cb94ee6SHong Zhang 
8227cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
8237cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
8247cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
8257cb94ee6SHong Zhang                     "TSSundialsSetGMRESRestart_Sundials",
8267cb94ee6SHong Zhang                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
8277cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
8287cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
8297cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
8307cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
8317cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
8327cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
8337cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
8347cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
8357cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
8367cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
8377cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
8387cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
8397cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
8407cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
8417cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
8427cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
8437cb94ee6SHong Zhang                     "TSSundialsSetExactFinalTime_Sundials",
8447cb94ee6SHong Zhang                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
8457cb94ee6SHong Zhang   PetscFunctionReturn(0);
8467cb94ee6SHong Zhang }
8477cb94ee6SHong Zhang EXTERN_C_END
8487cb94ee6SHong Zhang 
8497cb94ee6SHong Zhang 
8507cb94ee6SHong Zhang 
8517cb94ee6SHong Zhang 
8527cb94ee6SHong Zhang 
8537cb94ee6SHong Zhang 
8547cb94ee6SHong Zhang 
8557cb94ee6SHong Zhang 
8567cb94ee6SHong Zhang 
8577cb94ee6SHong Zhang 
858