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:   Note:
 35:   This is a struct, not a `PetscObject`

 37: .seealso: `DMPlexSNESComputeResidualFEM()`, `DMPlexSNESComputeJacobianFEM()`, `DMPlexSNESComputeBoundaryFEM()`
 38: S*/
 39: typedef struct _PetscFormKey {
 40:   DMLabel  label; /* The (label, value) select a subdomain */
 41:   PetscInt value;
 42:   PetscInt field; /* Selects the field for the test function */
 43:   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. */
 44: } PetscFormKey;

 46: /*E
 47:   PetscWeakFormKind - The kind of weak form. The specific forms are given in the documentation for the integraton functions.

 49:   Values:
 50: + OBJECTIVE                  - Objective form
 51: . F0, F1                     - Residual forms
 52: . G0, G1, G2, G3             - Jacobian forms
 53: . GP0, GP1, GP2, GP3         - Jacobian preconditioner matrix forms
 54: . GT0, GT1, GT2, GT3         - Dynamic Jacobian matrix forms
 55: . BDF0, BDF1                 - Boundary Residual forms
 56: . BDG0, BDG1, BDG2, BDG3     - Jacobian forms
 57: . BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian preconditioner matrix forms
 58: . R                          - Riemann solver
 59: - CEED                       - libCEED QFunction

 61:   Level: beginner

 63: .seealso: `PetscWeakForm`, `PetscFEIntegrateResidual()`, `PetscFEIntegrateJacobian()`, `PetscFEIntegrateBdResidual()`, `PetscFEIntegrateBdJacobian()`,
 64:           `PetscFVIntegrateRHSFunction()`, `PetscWeakFormSetIndexResidual()`, `PetscWeakFormClearIndex()`
 65: E*/
 66: typedef enum {
 67:   PETSC_WF_OBJECTIVE,
 68:   PETSC_WF_F0,
 69:   PETSC_WF_F1,
 70:   PETSC_WF_G0,
 71:   PETSC_WF_G1,
 72:   PETSC_WF_G2,
 73:   PETSC_WF_G3,
 74:   PETSC_WF_GP0,
 75:   PETSC_WF_GP1,
 76:   PETSC_WF_GP2,
 77:   PETSC_WF_GP3,
 78:   PETSC_WF_GT0,
 79:   PETSC_WF_GT1,
 80:   PETSC_WF_GT2,
 81:   PETSC_WF_GT3,
 82:   PETSC_WF_BDF0,
 83:   PETSC_WF_BDF1,
 84:   PETSC_WF_BDG0,
 85:   PETSC_WF_BDG1,
 86:   PETSC_WF_BDG2,
 87:   PETSC_WF_BDG3,
 88:   PETSC_WF_BDGP0,
 89:   PETSC_WF_BDGP1,
 90:   PETSC_WF_BDGP2,
 91:   PETSC_WF_BDGP3,
 92:   PETSC_WF_R,
 93:   PETSC_WF_CEED,
 94:   PETSC_NUM_WF
 95: } PetscWeakFormKind;
 96: PETSC_EXTERN const char *const PetscWeakFormKinds[];