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