xref: /petsc/src/ts/impls/explicit/rk/rk.h (revision 474dd7730651a168e39c2c0e6a0673d505ddd4a4)
1*474dd773SHong Zhang typedef struct _RKTableau *RKTableau;
2*474dd773SHong Zhang struct _RKTableau {
3*474dd773SHong Zhang   char       *name;
4*474dd773SHong Zhang   PetscInt   order;               /* Classical approximation order of the method i              */
5*474dd773SHong Zhang   PetscInt   s;                   /* Number of stages                                           */
6*474dd773SHong Zhang   PetscInt   p;                   /* Interpolation order                                        */
7*474dd773SHong Zhang   PetscBool  FSAL;                /* flag to indicate if tableau is FSAL                        */
8*474dd773SHong Zhang   PetscReal  *A,*b,*c;            /* Tableau                                                    */
9*474dd773SHong Zhang   PetscReal  *bembed;             /* Embedded formula of order one less (order-1)               */
10*474dd773SHong Zhang   PetscReal  *binterp;            /* Dense output formula                                       */
11*474dd773SHong Zhang   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
12*474dd773SHong Zhang };
13*474dd773SHong Zhang typedef struct _RKTableauLink *RKTableauLink;
14*474dd773SHong Zhang struct _RKTableauLink {
15*474dd773SHong Zhang   struct _RKTableau tab;
16*474dd773SHong Zhang   RKTableauLink     next;
17*474dd773SHong Zhang };
18*474dd773SHong Zhang static RKTableauLink RKTableauList;
19*474dd773SHong Zhang 
20*474dd773SHong Zhang typedef struct {
21*474dd773SHong Zhang   RKTableau    tableau;
22*474dd773SHong Zhang   Vec          X0;
23*474dd773SHong Zhang   Vec          *Y;               /* States computed during the step                                              */
24*474dd773SHong Zhang   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part and contains all components      */
25*474dd773SHong Zhang   Vec          *YdotRHS_fast;    /* Function evaluations for the non-stiff part and contains fast components     */
26*474dd773SHong Zhang   Vec          *YdotRHS_slow;    /* Function evaluations for the non-stiff part and contains slow components     */
27*474dd773SHong Zhang   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage                       */
28*474dd773SHong Zhang   Vec          *VecDeltaMu;      /* Increment of the adjoint sensitivity w.r.t P at stage                        */
29*474dd773SHong Zhang   Vec          VecCostIntegral0; /* backup for roll-backs due to events                                          */
30*474dd773SHong Zhang   PetscScalar  *work;            /* Scalar work                                                                  */
31*474dd773SHong Zhang   PetscInt     slow;             /* flag indicates call slow components solver (0) or fast components solver (1) */
32*474dd773SHong Zhang   PetscReal    stage_time;
33*474dd773SHong Zhang   TSStepStatus status;
34*474dd773SHong Zhang   PetscReal    ptime;
35*474dd773SHong Zhang   PetscReal    time_step;
36*474dd773SHong Zhang   PetscInt     dtratio;          /* ratio between slow time step size and fast step size                         */
37*474dd773SHong Zhang   IS           is_fast,is_slow;
38*474dd773SHong Zhang   TS           subts_fast,subts_slow;
39*474dd773SHong Zhang } TS_RK;
40