1 #pragma once 2 3 #include <petscdmlabel.h> 4 5 /* MANSEC = DM */ 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 Note: 36 This is a struct, not a `PetscObject` 37 38 .seealso: `DMPlexSNESComputeResidualFEM()`, `DMPlexSNESComputeJacobianFEM()`, `DMPlexSNESComputeBoundaryFEM()` 39 S*/ 40 typedef struct { 41 DMLabel label; /* The (label, value) select a subdomain */ 42 PetscInt value; 43 PetscInt field; /* Selects the field for the test function */ 44 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. */ 45 } PetscFormKey; 46 47 /*E 48 PetscWeakFormKind - The kind of weak form. The specific forms are given in the documentation for the integraton functions. 49 50 Values: 51 + OBJECTIVE - Objective form 52 . F0, F1 - Residual forms 53 . G0, G1, G2, G3 - Jacobian forms 54 . GP0, GP1, GP2, GP3 - Jacobian preconditioner matrix forms 55 . GT0, GT1, GT2, GT3 - Dynamic Jacobian matrix forms 56 . BDF0, BDF1 - Boundary Residual forms 57 . BDG0, BDG1, BDG2, BDG3 - Jacobian forms 58 . BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian preconditioner matrix forms 59 . R - Riemann solver 60 - CEED - libCEED QFunction 61 62 Level: beginner 63 64 .seealso: `PetscWeakForm`, `PetscFEIntegrateResidual()`, `PetscFEIntegrateJacobian()`, `PetscFEIntegrateBdResidual()`, `PetscFEIntegrateBdJacobian()`, 65 `PetscFVIntegrateRHSFunction()`, `PetscWeakFormSetIndexResidual()`, `PetscWeakFormClearIndex()` 66 E*/ 67 typedef enum { 68 PETSC_WF_OBJECTIVE, 69 PETSC_WF_F0, 70 PETSC_WF_F1, 71 PETSC_WF_G0, 72 PETSC_WF_G1, 73 PETSC_WF_G2, 74 PETSC_WF_G3, 75 PETSC_WF_GP0, 76 PETSC_WF_GP1, 77 PETSC_WF_GP2, 78 PETSC_WF_GP3, 79 PETSC_WF_GT0, 80 PETSC_WF_GT1, 81 PETSC_WF_GT2, 82 PETSC_WF_GT3, 83 PETSC_WF_BDF0, 84 PETSC_WF_BDF1, 85 PETSC_WF_BDG0, 86 PETSC_WF_BDG1, 87 PETSC_WF_BDG2, 88 PETSC_WF_BDG3, 89 PETSC_WF_BDGP0, 90 PETSC_WF_BDGP1, 91 PETSC_WF_BDGP2, 92 PETSC_WF_BDGP3, 93 PETSC_WF_R, 94 PETSC_WF_CEED, 95 PETSC_NUM_WF 96 } PetscWeakFormKind; 97 PETSC_EXTERN const char *const PetscWeakFormKinds[]; 98