126bd1501SBarry Smith #if !defined(PETSCDSTYPES_H) 226bd1501SBarry Smith #define PETSCDSTYPES_H 32764a2aaSMatthew G. Knepley 46528b96dSMatthew G. Knepley #include <petscdmlabel.h> 56528b96dSMatthew G. Knepley 6ac09b921SBarry Smith /* SUBMANSEC = DT */ 7ac09b921SBarry Smith 82764a2aaSMatthew G. Knepley /*S 987497f52SBarry Smith 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*/ 37*9371c9d4SSatish Balay typedef struct _PetscFormKey { 3806ad1575SMatthew G. Knepley DMLabel label; /* The (label, value) select a subdomain */ 396528b96dSMatthew G. Knepley PetscInt value; 4006ad1575SMatthew G. Knepley PetscInt field; /* Selects the field for the test function */ 4106ad1575SMatthew 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. */ 4206ad1575SMatthew G. Knepley } PetscFormKey; 4306ad1575SMatthew G. Knepley 4406ad1575SMatthew G. Knepley /*E 4506ad1575SMatthew G. Knepley PetscWeakFormKind - The kind of weak form. The specific forms are given in the documentation for the integraton functions. 4606ad1575SMatthew G. Knepley 4706ad1575SMatthew G. Knepley Supported kinds include: 4806ad1575SMatthew G. Knepley $ OBJECTIVE - Objective form 4906ad1575SMatthew G. Knepley $ F0, F1 - Residual forms 5006ad1575SMatthew G. Knepley $ G0, G1, G2, G3 - Jacobian forms 5106ad1575SMatthew G. Knepley $ GP0, GP1, GP2, GP3 - Jacobian preconditioner matrix forms 5206ad1575SMatthew G. Knepley $ GT0, GT1, GT2, GT3 - Dynamic Jacobian matrix forms 5306ad1575SMatthew G. Knepley $ BDF0, BDF1 - Boundary Residual forms 5406ad1575SMatthew G. Knepley $ BDG0, BDG1, BDG2, BDG3 - Jacobian forms 5506ad1575SMatthew G. Knepley $ BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian preconditioner matrix forms 5606ad1575SMatthew G. Knepley $ R - Riemann solver 5706ad1575SMatthew G. Knepley 5806ad1575SMatthew G. Knepley Level: beginner 5906ad1575SMatthew G. Knepley 60db781477SPatrick Sanan .seealso: `PetscFEIntegrateResidual()`, `PetscFEIntegrateJacobian()`, `PetscFEIntegrateBdResidual()`, `PetscFEIntegrateBdJacobian()`, `PetscFVIntegrateRHSFunction()`, `PetscWeakFormSetIndexResidual()`, `PetscWeakFormClearIndex()` 6106ad1575SMatthew G. Knepley E*/ 62*9371c9d4SSatish Balay typedef enum { 63*9371c9d4SSatish Balay PETSC_WF_OBJECTIVE, 64*9371c9d4SSatish Balay PETSC_WF_F0, 65*9371c9d4SSatish Balay PETSC_WF_F1, 66*9371c9d4SSatish Balay PETSC_WF_G0, 67*9371c9d4SSatish Balay PETSC_WF_G1, 68*9371c9d4SSatish Balay PETSC_WF_G2, 69*9371c9d4SSatish Balay PETSC_WF_G3, 70*9371c9d4SSatish Balay PETSC_WF_GP0, 71*9371c9d4SSatish Balay PETSC_WF_GP1, 72*9371c9d4SSatish Balay PETSC_WF_GP2, 73*9371c9d4SSatish Balay PETSC_WF_GP3, 74*9371c9d4SSatish Balay PETSC_WF_GT0, 75*9371c9d4SSatish Balay PETSC_WF_GT1, 76*9371c9d4SSatish Balay PETSC_WF_GT2, 77*9371c9d4SSatish Balay PETSC_WF_GT3, 78*9371c9d4SSatish Balay PETSC_WF_BDF0, 79*9371c9d4SSatish Balay PETSC_WF_BDF1, 80*9371c9d4SSatish Balay PETSC_WF_BDG0, 81*9371c9d4SSatish Balay PETSC_WF_BDG1, 82*9371c9d4SSatish Balay PETSC_WF_BDG2, 83*9371c9d4SSatish Balay PETSC_WF_BDG3, 84*9371c9d4SSatish Balay PETSC_WF_BDGP0, 85*9371c9d4SSatish Balay PETSC_WF_BDGP1, 86*9371c9d4SSatish Balay PETSC_WF_BDGP2, 87*9371c9d4SSatish Balay PETSC_WF_BDGP3, 88*9371c9d4SSatish Balay PETSC_WF_R, 89*9371c9d4SSatish Balay PETSC_NUM_WF 90*9371c9d4SSatish Balay } PetscWeakFormKind; 9106ad1575SMatthew G. Knepley PETSC_EXTERN const char *const PetscWeakFormKinds[]; 926528b96dSMatthew G. Knepley 932764a2aaSMatthew G. Knepley #endif 94