Actual source code: rk.h
1: #pragma once
2: typedef struct _RKTableau *RKTableau;
3: struct _RKTableau {
4: char *name;
5: PetscInt order; /* Classical approximation order of the method i */
6: PetscInt s; /* Number of stages */
7: PetscInt p; /* Interpolation order */
8: PetscBool FSAL; /* flag to indicate if tableau is FSAL */
9: PetscReal *A, *b, *c; /* Tableau */
10: PetscReal *bembed; /* Embedded formula of order one less (order-1) */
11: PetscReal *binterp; /* Dense output formula */
12: PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */
13: };
14: typedef struct _RKTableauLink *RKTableauLink;
15: struct _RKTableauLink {
16: struct _RKTableau tab;
17: RKTableauLink next;
18: };
20: typedef struct {
21: RKTableau tableau;
22: PetscBool newtableau; /* flag to indicate if tableau has changed */
23: Vec X0;
24: Vec *Y; /* States computed during the step */
25: Vec *YdotRHS; /* Function evaluations for the non-stiff part and contains all components */
26: Vec *YdotRHS_fast; /* Function evaluations for the non-stiff part and contains fast components */
27: Vec *YdotRHS_slow; /* Function evaluations for the non-stiff part and contains slow components */
28: Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */
29: Vec *VecsSensiTemp;
30: Vec VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */
31: Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
32: Vec VecDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
33: Vec *VecsSensi2Temp;
34: PetscScalar *work; /* Scalar work */
35: PetscInt slow; /* flag indicates call slow components solver (0) or fast components solver (1) */
36: PetscReal stage_time;
37: TSStepStatus status;
38: PetscReal ptime;
39: PetscReal time_step;
40: PetscInt dtratio; /* ratio between slow time step size and fast step size */
41: IS is_fast, is_slow;
42: TS subts_fast, subts_slow, subts_current, ts_root;
43: PetscBool use_multirate;
44: Mat MatFwdSensip0;
45: Mat *MatsFwdStageSensip;
46: Mat *MatsFwdSensipTemp;
47: Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */
48: } TS_RK;