xref: /petsc/src/ts/impls/explicit/rk/rk.h (revision 630f8c86da6f35c90b6972842e3fb3a8d0c014ed)
1474dd773SHong Zhang typedef struct _RKTableau *RKTableau;
2474dd773SHong Zhang struct _RKTableau {
3474dd773SHong Zhang   char      *name;
4474dd773SHong Zhang   PetscInt   order;     /* Classical approximation order of the method i              */
5474dd773SHong Zhang   PetscInt   s;         /* Number of stages                                           */
6474dd773SHong Zhang   PetscInt   p;         /* Interpolation order                                        */
7474dd773SHong Zhang   PetscBool  FSAL;      /* flag to indicate if tableau is FSAL                        */
8474dd773SHong Zhang   PetscReal *A, *b, *c; /* Tableau                                                    */
9474dd773SHong Zhang   PetscReal *bembed;    /* Embedded formula of order one less (order-1)               */
10474dd773SHong Zhang   PetscReal *binterp;   /* Dense output formula                                       */
11474dd773SHong Zhang   PetscReal  ccfl;      /* Placeholder for CFL coefficient relative to forward Euler  */
12474dd773SHong Zhang };
13474dd773SHong Zhang typedef struct _RKTableauLink *RKTableauLink;
14474dd773SHong Zhang struct _RKTableauLink {
15474dd773SHong Zhang   struct _RKTableau tab;
16474dd773SHong Zhang   RKTableauLink     next;
17474dd773SHong Zhang };
18474dd773SHong Zhang 
19474dd773SHong Zhang typedef struct {
20474dd773SHong Zhang   RKTableau    tableau;
21*630f8c86SStefano Zampini   PetscBool    newtableau; /* flag to indicate if tableau has changed */
22474dd773SHong Zhang   Vec          X0;
23474dd773SHong Zhang   Vec         *Y;            /* States computed during the step                                              */
24474dd773SHong Zhang   Vec         *YdotRHS;      /* Function evaluations for the non-stiff part and contains all components      */
25474dd773SHong Zhang   Vec         *YdotRHS_fast; /* Function evaluations for the non-stiff part and contains fast components     */
26474dd773SHong Zhang   Vec         *YdotRHS_slow; /* Function evaluations for the non-stiff part and contains slow components     */
272e7b7f96SHong Zhang   Vec         *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage                       */
282e7b7f96SHong Zhang   Vec         *VecsSensiTemp;
292e7b7f96SHong Zhang   Vec          VecDeltaMu;    /* Increment of the adjoint sensitivity w.r.t P at stage                        */
3013af1a74SHong Zhang   Vec         *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
3113af1a74SHong Zhang   Vec          VecDeltaMu2;   /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
3213af1a74SHong Zhang   Vec         *VecsSensi2Temp;
33474dd773SHong Zhang   PetscScalar *work; /* Scalar work                                                                  */
34474dd773SHong Zhang   PetscInt     slow; /* flag indicates call slow components solver (0) or fast components solver (1) */
35474dd773SHong Zhang   PetscReal    stage_time;
36474dd773SHong Zhang   TSStepStatus status;
37474dd773SHong Zhang   PetscReal    ptime;
38474dd773SHong Zhang   PetscReal    time_step;
39474dd773SHong Zhang   PetscInt     dtratio; /* ratio between slow time step size and fast step size                         */
40474dd773SHong Zhang   IS           is_fast, is_slow;
4163a6f1b4SHong Zhang   TS           subts_fast, subts_slow, subts_current, ts_root;
420fe4d17eSHong Zhang   PetscBool    use_multirate;
43922a638cSHong Zhang   Mat          MatFwdSensip0;
44922a638cSHong Zhang   Mat         *MatsFwdStageSensip;
45922a638cSHong Zhang   Mat         *MatsFwdSensipTemp;
46922a638cSHong Zhang   Vec          VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */
47474dd773SHong Zhang } TS_RK;
48