Actual source code: petscdstypes.h

  1: #if !defined(PETSCDSTYPES_H)
  2: #define PETSCDSTYPES_H

  4: #include <petscdmlabel.h>

  6: /*S
  7:   PetscDS - PETSc object that manages a discrete system, which is a set of discretizations + continuum equations from a PetscWeakForm

  9:   Level: intermediate

 11: .seealso: PetscDSCreate(), PetscDSSetType(), PetscDSType, PetscWeakForm, PetscFECreate(), PetscFVCreate()
 12: S*/
 13: typedef struct _p_PetscDS *PetscDS;

 15: /*S
 16:   PetscWeakForm - PETSc object that manages a sets of pointwise functions defining a system of equations

 18:   Level: intermediate

 20: .seealso: PetscWeakFormCreate(), PetscDS, PetscFECreate(), PetscFVCreate()
 21: S*/
 22: typedef struct _p_PetscWeakForm *PetscWeakForm;

 24: /*S
 25:   PetscFormKey - This key indicates how to use a set of pointwise functions defining part of a system of equations

 27:   The subdomain on which to integrate is specified by (label, value), the test function field by (field), and the
 28:   piece of the equation by (part). For example, LHS = 0 and RHS = 1 in IMEX methods. More pieces can be present for
 29:   operator splitting methods.

 31:   Level: intermediate

 33: .seealso: DMPlexSNESComputeResidualFEM(), DMPlexSNESComputeJacobianFEM(), DMPlexSNESComputeBoundaryFEM()
 34: S*/
 35: typedef struct _PetscFormKey
 36: {
 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:   Supported kinds include:
 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

 57:   Level: beginner

 59: .seealso: PetscFEIntegrateResidual(), PetscFEIntegrateJacobian(), PetscFEIntegrateBdResidual(), PetscFEIntegrateBdJacobian(), PetscFVIntegrateRHSFunction(), PetscWeakFormSetIndexResidual(), PetscWeakFormClearIndex()
 60: E*/
 61: 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;
 62: PETSC_EXTERN const char *const PetscWeakFormKinds[];

 64: #endif