1 #ifndef PETSCDSTYPES_H 2 #define PETSCDSTYPES_H 3 4 #include <petscdmlabel.h> 5 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 .seealso: `DMPlexSNESComputeResidualFEM()`, `DMPlexSNESComputeJacobianFEM()`, `DMPlexSNESComputeBoundaryFEM()` 36 S*/ 37 typedef struct _PetscFormKey { 38 DMLabel label; /* The (label, value) select a subdomain */ 39 PetscInt value; 40 PetscInt field; /* Selects the field for the test function */ 41 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. */ 42 } PetscFormKey; 43 44 /*E 45 PetscWeakFormKind - The kind of weak form. The specific forms are given in the documentation for the integraton functions. 46 47 Values: 48 + OBJECTIVE - Objective form 49 . F0, F1 - Residual forms 50 . G0, G1, G2, G3 - Jacobian forms 51 . GP0, GP1, GP2, GP3 - Jacobian preconditioner matrix forms 52 . GT0, GT1, GT2, GT3 - Dynamic Jacobian matrix forms 53 . BDF0, BDF1 - Boundary Residual forms 54 . BDG0, BDG1, BDG2, BDG3 - Jacobian forms 55 . BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian preconditioner matrix forms 56 - R - Riemann solver 57 58 Level: beginner 59 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_NUM_WF 91 } PetscWeakFormKind; 92 PETSC_EXTERN const char *const PetscWeakFormKinds[]; 93 94 #endif 95