xref: /petsc/src/ts/impls/implicit/sundials/sundials.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840) !
17cb94ee6SHong Zhang /*
27cb94ee6SHong Zhang     Provides a PETSc interface to SUNDIALS. Alan Hindmarsh's parallel ODE
37cb94ee6SHong Zhang    solver developed at LLNL.
47cb94ee6SHong Zhang */
57cb94ee6SHong Zhang 
6*a4963045SJacob Faibussowitsch #pragma once
77cb94ee6SHong Zhang 
8af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
9af0996ceSBarry Smith #include <petsc/private/pcimpl.h>
10af0996ceSBarry Smith #include <petsc/private/matimpl.h>
117cb94ee6SHong Zhang 
127cb94ee6SHong Zhang /*
137cb94ee6SHong Zhang    Include files specific for SUNDIALS
147cb94ee6SHong Zhang */
15e808b789SPatrick Sanan #if defined(PETSC_HAVE_SUNDIALS2)
167cb94ee6SHong Zhang 
177cb94ee6SHong Zhang EXTERN_C_BEGIN
18c6db04a5SJed Brown   #include <cvode/cvode.h>              /* prototypes for CVODE fcts. */
19c6db04a5SJed Brown   #include <cvode/cvode_spgmr.h>        /* prototypes and constants for CVSPGMR solver */
20918687b8SPatrick Sanan   #include <cvode/cvode_dense.h>        /* prototypes and constants for CVDense solver */
21c6db04a5SJed Brown   #include <nvector/nvector_parallel.h> /* definition N_Vector and macro NV_DATA_P  */
22918687b8SPatrick Sanan   #include <nvector/nvector_serial.h>
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 
30918687b8SPatrick Sanan   /* PETSc preconditioner 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;
38c4e80e11SFlorian   PetscReal reltol;
39c4e80e11SFlorian   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 */
44c4e80e11SFlorian   PetscInt  maxord;      /* max order of BDF / Adams method */
45918687b8SPatrick Sanan   PetscBool use_dense;   /* Use a dense instead of iterative solve within SUNDIALS (serial only) */
467cb94ee6SHong Zhang } TS_Sundials;
477cb94ee6SHong Zhang #endif
48