xref: /petsc/include/petscdstypes.h (revision df4cd43f92eaa320656440c40edb1046daee8f75)
1 #ifndef PETSCDSTYPES_H
2 #define PETSCDSTYPES_H
3 
4 #include <petscdmlabel.h>
5 
6 /* SUBMANSEC = DT */
7 
8 /*S
9   PetscDS - PETSc object that manages a discrete system, which is a set of discretizations + continuum equations from a `PetscWeakForm`
10 
11   Level: intermediate
12 
13 .seealso: `PetscDSCreate()`, `PetscDSSetType()`, `PetscDSType`, `PetscWeakForm`, `PetscFECreate()`, `PetscFVCreate()`
14 S*/
15 typedef struct _p_PetscDS *PetscDS;
16 
17 /*S
18   PetscWeakForm - PETSc object that manages a sets of pointwise functions defining a system of equations
19 
20   Level: intermediate
21 
22 .seealso: `PetscWeakFormCreate()`, `PetscDS`, `PetscFECreate()`, `PetscFVCreate()`
23 S*/
24 typedef struct _p_PetscWeakForm *PetscWeakForm;
25 
26 /*S
27   PetscFormKey - This key indicates how to use a set of pointwise functions defining part of a system of equations
28 
29   The subdomain on which to integrate is specified by (label, value), the test function field by (field), and the
30   piece of the equation by (part). For example, LHS = 0 and RHS = 1 in IMEX methods. More pieces can be present for
31   operator splitting methods.
32 
33   Level: intermediate
34 
35 .seealso: `DMPlexSNESComputeResidualFEM()`, `DMPlexSNESComputeJacobianFEM()`, `DMPlexSNESComputeBoundaryFEM()`
36 S*/
37 typedef struct _PetscFormKey {
38   DMLabel  label; /* The (label, value) select a subdomain */
39   PetscInt value;
40   PetscInt field; /* Selects the field for the test function */
41   PetscInt part;  /* Selects the equation part. For example, LHS = 0 and RHS = 1 in IMEX methods. More pieces can be present for operator splitting methods. */
42 } PetscFormKey;
43 
44 /*E
45   PetscWeakFormKind - The kind of weak form. The specific forms are given in the documentation for the integraton functions.
46 
47   Supported kinds include:
48 $ OBJECTIVE                  - Objective form
49 $ F0, F1                     - Residual forms
50 $ G0, G1, G2, G3             - Jacobian forms
51 $ GP0, GP1, GP2, GP3         - Jacobian preconditioner matrix forms
52 $ GT0, GT1, GT2, GT3         - Dynamic Jacobian matrix forms
53 $ BDF0, BDF1                 - Boundary Residual forms
54 $ BDG0, BDG1, BDG2, BDG3     - Jacobian forms
55 $ BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian preconditioner matrix forms
56 $ R                          - Riemann solver
57 
58   Level: beginner
59 
60 .seealso: `PetscFEIntegrateResidual()`, `PetscFEIntegrateJacobian()`, `PetscFEIntegrateBdResidual()`, `PetscFEIntegrateBdJacobian()`, `PetscFVIntegrateRHSFunction()`, `PetscWeakFormSetIndexResidual()`, `PetscWeakFormClearIndex()`
61 E*/
62 typedef enum {
63   PETSC_WF_OBJECTIVE,
64   PETSC_WF_F0,
65   PETSC_WF_F1,
66   PETSC_WF_G0,
67   PETSC_WF_G1,
68   PETSC_WF_G2,
69   PETSC_WF_G3,
70   PETSC_WF_GP0,
71   PETSC_WF_GP1,
72   PETSC_WF_GP2,
73   PETSC_WF_GP3,
74   PETSC_WF_GT0,
75   PETSC_WF_GT1,
76   PETSC_WF_GT2,
77   PETSC_WF_GT3,
78   PETSC_WF_BDF0,
79   PETSC_WF_BDF1,
80   PETSC_WF_BDG0,
81   PETSC_WF_BDG1,
82   PETSC_WF_BDG2,
83   PETSC_WF_BDG3,
84   PETSC_WF_BDGP0,
85   PETSC_WF_BDGP1,
86   PETSC_WF_BDGP2,
87   PETSC_WF_BDGP3,
88   PETSC_WF_R,
89   PETSC_NUM_WF
90 } PetscWeakFormKind;
91 PETSC_EXTERN const char *const PetscWeakFormKinds[];
92 
93 #endif
94