xref: /petsc/src/ts/impls/implicit/sundials/sundials.h (revision c4e80e11bf140eefea2ededde164a42757847390)
17cb94ee6SHong Zhang 
27cb94ee6SHong Zhang /*
37cb94ee6SHong Zhang     Provides a PETSc interface to SUNDIALS. Alan Hindmarsh's parallel ODE
47cb94ee6SHong Zhang    solver developed at LLNL.
57cb94ee6SHong Zhang */
67cb94ee6SHong Zhang 
742e3df96SBarry Smith #if !defined(PETSC_SUNDIALS_H)
826bd1501SBarry Smith #define PETSC_SUNDIALS_H
97cb94ee6SHong Zhang 
10af0996ceSBarry Smith #include <petsc/private/tsimpl.h>       /*I   "petscts.h"   I*/
11af0996ceSBarry Smith #include <petsc/private/pcimpl.h>
12af0996ceSBarry Smith #include <petsc/private/matimpl.h>
137cb94ee6SHong Zhang 
147cb94ee6SHong Zhang /*
157cb94ee6SHong Zhang    Include files specific for SUNDIALS
167cb94ee6SHong Zhang */
177cb94ee6SHong Zhang #if defined(PETSC_HAVE_SUNDIALS)
187cb94ee6SHong Zhang 
197cb94ee6SHong Zhang EXTERN_C_BEGIN
20c6db04a5SJed Brown #include <cvode/cvode.h>                  /* prototypes for CVODE fcts. */
21c6db04a5SJed Brown #include <cvode/cvode_spgmr.h>            /* prototypes and constants for CVSPGMR solver */
22c6db04a5SJed Brown #include <nvector/nvector_parallel.h>     /* definition N_Vector and macro NV_DATA_P  */
237cb94ee6SHong Zhang EXTERN_C_END
247cb94ee6SHong Zhang 
257cb94ee6SHong Zhang typedef struct {
267cb94ee6SHong Zhang   Vec update;           /* work vector where new solution is formed */
270679a0aeSJed Brown   Vec ydot;             /* work vector the time derivative is stored */
287cb94ee6SHong Zhang   Vec w1,w2;            /* work space vectors for function evaluation */
29b4eba00bSSean Farley 
307cb94ee6SHong Zhang   /* PETSc peconditioner objects used by SUNDIALS */
31f61b2b6aSHong Zhang   PetscInt                  cvode_type;   /* the SUNDIALS method, BDF or ADAMS  */
327cb94ee6SHong Zhang   TSSundialsGramSchmidtType gtype;
33f61b2b6aSHong Zhang   PetscReal                 linear_tol;
34f1cd61daSJed Brown   PetscReal                 mindt,maxdt;
35f1cd61daSJed Brown 
367cb94ee6SHong Zhang   /* Variables used by Sundials */
377cb94ee6SHong Zhang   MPI_Comm  comm_sundials;
38*c4e80e11SFlorian   PetscReal reltol;
39*c4e80e11SFlorian   PetscReal abstol;          /* only for using SS flag in SUNDIALS */
407cb94ee6SHong Zhang   N_Vector  y;               /* current solution */
412c823083SHong Zhang   void      *mem;
42ace3abfcSBarry Smith   PetscBool monitorstep;     /* flag for monitor internal steps; itask=V_ONE_STEP or itask=CV_NORMAL*/
43f61b2b6aSHong Zhang   PetscInt  maxl;            /* max dimension of the Krylov subspace to be used */
44*c4e80e11SFlorian   PetscInt  maxord;          /* max order of BDF / Adams method */
457cb94ee6SHong Zhang } TS_Sundials;
467cb94ee6SHong Zhang #endif
477cb94ee6SHong Zhang 
487cb94ee6SHong Zhang #endif
497cb94ee6SHong Zhang 
507cb94ee6SHong Zhang 
517cb94ee6SHong Zhang 
527cb94ee6SHong Zhang 
53