xref: /petsc/include/petscdstypes.h (revision ac09b9214d23ea9ad238aa607de9fa447fd4e91b)
126bd1501SBarry Smith #if !defined(PETSCDSTYPES_H)
226bd1501SBarry Smith #define PETSCDSTYPES_H
32764a2aaSMatthew G. Knepley 
46528b96dSMatthew G. Knepley #include <petscdmlabel.h>
56528b96dSMatthew G. Knepley 
6*ac09b921SBarry Smith /* SUBMANSEC = DT */
7*ac09b921SBarry Smith 
82764a2aaSMatthew G. Knepley /*S
96528b96dSMatthew G. Knepley   PetscDS - PETSc object that manages a discrete system, which is a set of discretizations + continuum equations from a PetscWeakForm
102764a2aaSMatthew G. Knepley 
112764a2aaSMatthew G. Knepley   Level: intermediate
122764a2aaSMatthew G. Knepley 
13db781477SPatrick Sanan .seealso: `PetscDSCreate()`, `PetscDSSetType()`, `PetscDSType`, `PetscWeakForm`, `PetscFECreate()`, `PetscFVCreate()`
142764a2aaSMatthew G. Knepley S*/
152764a2aaSMatthew G. Knepley typedef struct _p_PetscDS *PetscDS;
162764a2aaSMatthew G. Knepley 
176528b96dSMatthew G. Knepley /*S
186528b96dSMatthew G. Knepley   PetscWeakForm - PETSc object that manages a sets of pointwise functions defining a system of equations
196528b96dSMatthew G. Knepley 
206528b96dSMatthew G. Knepley   Level: intermediate
216528b96dSMatthew G. Knepley 
22db781477SPatrick Sanan .seealso: `PetscWeakFormCreate()`, `PetscDS`, `PetscFECreate()`, `PetscFVCreate()`
236528b96dSMatthew G. Knepley S*/
246528b96dSMatthew G. Knepley typedef struct _p_PetscWeakForm *PetscWeakForm;
256528b96dSMatthew G. Knepley 
2606ad1575SMatthew G. Knepley /*S
2706ad1575SMatthew G. Knepley   PetscFormKey - This key indicates how to use a set of pointwise functions defining part of a system of equations
2806ad1575SMatthew G. Knepley 
29a5b23f4aSJose E. Roman   The subdomain on which to integrate is specified by (label, value), the test function field by (field), and the
3006ad1575SMatthew G. Knepley   piece of the equation by (part). For example, LHS = 0 and RHS = 1 in IMEX methods. More pieces can be present for
3106ad1575SMatthew G. Knepley   operator splitting methods.
3206ad1575SMatthew G. Knepley 
3306ad1575SMatthew G. Knepley   Level: intermediate
3406ad1575SMatthew G. Knepley 
35db781477SPatrick Sanan .seealso: `DMPlexSNESComputeResidualFEM()`, `DMPlexSNESComputeJacobianFEM()`, `DMPlexSNESComputeBoundaryFEM()`
3606ad1575SMatthew G. Knepley S*/
3706ad1575SMatthew G. Knepley typedef struct _PetscFormKey
386528b96dSMatthew G. Knepley {
3906ad1575SMatthew G. Knepley   DMLabel  label; /* The (label, value) select a subdomain */
406528b96dSMatthew G. Knepley   PetscInt value;
4106ad1575SMatthew G. Knepley   PetscInt field; /* Selects the field for the test function */
4206ad1575SMatthew G. Knepley   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. */
4306ad1575SMatthew G. Knepley } PetscFormKey;
4406ad1575SMatthew G. Knepley 
4506ad1575SMatthew G. Knepley /*E
4606ad1575SMatthew G. Knepley   PetscWeakFormKind - The kind of weak form. The specific forms are given in the documentation for the integraton functions.
4706ad1575SMatthew G. Knepley 
4806ad1575SMatthew G. Knepley   Supported kinds include:
4906ad1575SMatthew G. Knepley $ OBJECTIVE                  - Objective form
5006ad1575SMatthew G. Knepley $ F0, F1                     - Residual forms
5106ad1575SMatthew G. Knepley $ G0, G1, G2, G3             - Jacobian forms
5206ad1575SMatthew G. Knepley $ GP0, GP1, GP2, GP3         - Jacobian preconditioner matrix forms
5306ad1575SMatthew G. Knepley $ GT0, GT1, GT2, GT3         - Dynamic Jacobian matrix forms
5406ad1575SMatthew G. Knepley $ BDF0, BDF1                 - Boundary Residual forms
5506ad1575SMatthew G. Knepley $ BDG0, BDG1, BDG2, BDG3     - Jacobian forms
5606ad1575SMatthew G. Knepley $ BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian preconditioner matrix forms
5706ad1575SMatthew G. Knepley $ R                          - Riemann solver
5806ad1575SMatthew G. Knepley 
5906ad1575SMatthew G. Knepley   Level: beginner
6006ad1575SMatthew G. Knepley 
61db781477SPatrick Sanan .seealso: `PetscFEIntegrateResidual()`, `PetscFEIntegrateJacobian()`, `PetscFEIntegrateBdResidual()`, `PetscFEIntegrateBdJacobian()`, `PetscFVIntegrateRHSFunction()`, `PetscWeakFormSetIndexResidual()`, `PetscWeakFormClearIndex()`
6206ad1575SMatthew G. Knepley E*/
6306ad1575SMatthew G. Knepley typedef enum {PETSC_WF_OBJECTIVE, PETSC_WF_F0, PETSC_WF_F1, PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3, PETSC_WF_GP0, PETSC_WF_GP1, PETSC_WF_GP2, PETSC_WF_GP3, PETSC_WF_GT0, PETSC_WF_GT1, PETSC_WF_GT2, PETSC_WF_GT3, PETSC_WF_BDF0, PETSC_WF_BDF1, PETSC_WF_BDG0, PETSC_WF_BDG1, PETSC_WF_BDG2, PETSC_WF_BDG3, PETSC_WF_BDGP0, PETSC_WF_BDGP1, PETSC_WF_BDGP2, PETSC_WF_BDGP3, PETSC_WF_R, PETSC_NUM_WF} PetscWeakFormKind;
6406ad1575SMatthew G. Knepley PETSC_EXTERN const char *const PetscWeakFormKinds[];
656528b96dSMatthew G. Knepley 
662764a2aaSMatthew G. Knepley #endif
67