Actual source code: arkimex.h
1: typedef struct _ARKTableau *ARKTableau;
2: struct _ARKTableau {
3: char *name;
4: PetscBool additive; /* If False, it is a DIRK method */
5: PetscInt order; /* Classical approximation order of the method */
6: PetscInt s; /* Number of stages */
7: PetscBool stiffly_accurate; /* The implicit part is stiffly accurate */
8: PetscBool FSAL_implicit; /* The implicit part is FSAL */
9: PetscBool explicit_first_stage; /* The implicit part has an explicit first stage */
10: PetscInt pinterp; /* Interpolation order */
11: PetscReal *At, *bt, *ct; /* Stiff tableau */
12: PetscReal *A, *b, *c; /* Non-stiff tableau */
13: PetscReal *bembedt, *bembed; /* Embedded formula of order one less (order-1) */
14: PetscReal *binterpt, *binterp; /* Dense output formula */
15: PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */
16: };
17: typedef struct _ARKTableauLink *ARKTableauLink;
18: struct _ARKTableauLink {
19: struct _ARKTableau tab;
20: ARKTableauLink next;
21: };
23: typedef struct {
24: ARKTableau tableau;
25: Vec *Y; /* States computed during the step */
26: Vec *YdotI; /* Time derivatives for the stiff part */
27: Vec *YdotRHS; /* Function evaluations for the non-stiff part */
28: Vec *Y_prev; /* States computed during the previous time step */
29: Vec *YdotI_prev; /* Time derivatives for the stiff part for the previous time step*/
30: Vec *YdotRHS_prev; /* Function evaluations for the non-stiff part for the previous time step*/
31: Vec Ydot0; /* Holds the slope from the previous step in FSAL case */
32: Vec Ydot; /* Work vector holding Ydot during residual evaluation */
33: Vec Z; /* Ydot = shift(Y-Z) */
34: IS alg_is; /* Index set for algebraic variables, needed when restarting with DIRK */
35: PetscScalar *work; /* Scalar work */
36: PetscReal scoeff; /* shift = scoeff/dt */
37: PetscReal stage_time;
38: PetscBool imex;
39: PetscBool extrapolate; /* Extrapolate initial guess from previous time-step stage values */
40: TSStepStatus status;
42: /* context for fast-slow split */
43: Vec Y_snes; /* Work vector for SNES */
44: Vec *YdotI_fast; /* Function evaluations for the fast components in YdotI */
45: Vec *YdotRHS_fast; /* Function evaluations for the fast components in YdotRHS */
46: Vec *YdotRHS_slow; /* Function evaluations for the slow components in YdotRHS */
47: IS is_slow, is_fast;
48: TS subts_slow, subts_fast;
49: PetscBool fastslowsplit;
51: /* context for sensitivity analysis */
52: Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */
53: Vec *VecsSensiTemp; /* Vectors to be multiplied with Jacobian transpose */
54: Vec *VecsSensiPTemp; /* Temporary Vectors to store JacobianP-transpose-vector product */
55: } TS_ARKIMEX;