Actual source code: petscdsimpl.h
1: #pragma once
3: #include <petscds.h>
4: #include <petsc/private/petscimpl.h>
5: #include <petsc/private/hashmap.h>
7: PETSC_EXTERN PetscBool PetscDSRegisterAllCalled;
8: PETSC_EXTERN PetscErrorCode PetscDSRegisterAll(void);
10: typedef struct _n_DSBoundary *DSBoundary;
12: struct _n_DSBoundary {
13: const char *name; /* A unique identifier for the condition */
14: DMLabel label; /* The DMLabel indicating the mesh region over which the condition holds */
15: const char *lname; /* The label name if the label is missing from the DM */
16: PetscInt Nv; /* The number of label values defining the region */
17: PetscInt *values; /* The labels values defining the region */
18: PetscWeakForm wf; /* Holds the pointwise functions defining the form (only for NATURAL conditions) */
19: DMBoundaryConditionType type; /* The type of condition, usually either ESSENTIAL or NATURAL */
20: PetscInt field; /* The field constrained by the condition */
21: PetscInt Nc; /* The number of constrained field components */
22: PetscInt *comps; /* The constrained field components */
23: void (*func)(void); /* Function that provides the boundary values (only for ESSENTIAL conditions) */
24: void (*func_t)(void); /* Function that provides the time derivative of the boundary values (only for ESSENTIAL conditions) */
25: void *ctx; /* The user context for func and func_t */
26: DSBoundary next;
27: };
29: typedef struct {
30: PetscInt start; /* Starting entry of the chunk in an array (in bytes) */
31: PetscInt size; /* Current number of entries of the chunk */
32: PetscInt reserved; /* Number of reserved entries in the chunk */
33: } PetscChunk;
35: typedef struct {
36: size_t size; /* Current number of entries used in array */
37: size_t alloc; /* Number of bytes allocated for array */
38: size_t unitbytes; /* Number of bytes per entry */
39: char *array;
40: } PetscChunkBuffer;
42: #define PetscFormKeyHash(key) PetscHashCombine(PetscHashCombine(PetscHashCombine(PetscHashPointer((key).label), PetscHashInt((key).value)), PetscHashInt((key).field)), PetscHashInt((key).part))
44: #define PetscFormKeyEqual(k1, k2) (((k1).label == (k2).label) ? ((k1).value == (k2).value) ? ((k1).field == (k2).field) ? ((k1).part == (k2).part) : 0 : 0 : 0)
46: static PetscChunk _PetscInvalidChunk = {-1, -1, -1};
48: PETSC_HASH_MAP(HMapForm, PetscFormKey, PetscChunk, PetscFormKeyHash, PetscFormKeyEqual, _PetscInvalidChunk)
50: /*
51: We sort lexicographically on the structure.
52: Returns
53: -1: left < right
54: 0: left = right
55: 1: left > right
56: */
57: static inline int Compare_PetscFormKey_Private(const void *left, const void *right, PETSC_UNUSED void *ctx)
58: {
59: PetscFormKey l = *(const PetscFormKey *)left;
60: PetscFormKey r = *(const PetscFormKey *)right;
61: return (l.label < r.label) ? -1 : ((l.label > r.label) ? 1 : ((l.value < r.value) ? -1 : (l.value > r.value) ? 1 : ((l.field < r.field) ? -1 : (l.field > r.field) ? 1 : ((l.part < r.part) ? -1 : (l.part > r.part)))));
62: }
64: typedef struct _PetscWeakFormOps *PetscWeakFormOps;
65: struct _PetscWeakFormOps {
66: PetscErrorCode (*setfromoptions)(PetscWeakForm);
67: PetscErrorCode (*setup)(PetscWeakForm);
68: PetscErrorCode (*view)(PetscWeakForm, PetscViewer);
69: PetscErrorCode (*destroy)(PetscWeakForm);
70: };
72: struct _p_PetscWeakForm {
73: PETSCHEADER(struct _PetscWeakFormOps);
74: void *data; /* Implementation object */
76: PetscInt Nf; /* The number of fields in the system */
77: PetscChunkBuffer *funcs; /* Buffer holding all function pointers */
78: PetscHMapForm *form; /* Stores function pointers for forms */
79: };
81: typedef struct _PetscDSOps *PetscDSOps;
82: struct _PetscDSOps {
83: PetscErrorCode (*setfromoptions)(PetscDS);
84: PetscErrorCode (*setup)(PetscDS);
85: PetscErrorCode (*view)(PetscDS, PetscViewer);
86: PetscErrorCode (*destroy)(PetscDS);
87: };
89: struct _p_PetscDS {
90: PETSCHEADER(struct _PetscDSOps);
91: void *data; /* Implementation object */
92: PetscDS *subprobs; /* The subspaces for each dimension */
93: PetscBool setup; /* Flag for setup */
94: PetscInt dimEmbed; /* The real space coordinate dimension */
95: PetscInt Nf; /* The number of solution fields */
96: PetscObject *disc; /* The discretization for each solution field (PetscFE, PetscFV, etc.) */
97: PetscBool *cohesive; /* Flag for cohesive discretization */
98: PetscBool isCohesive; /* We are on a cohesive cell, meaning lower dimensional FE used on a 0-volume cell. Normal fields appear on both endcaps, whereas cohesive field only appear once in the middle */
99: /* Quadrature */
100: PetscBool forceQuad; /* Flag to force matching quadratures in discretizations */
101: IS *quadPerm[DM_NUM_POLYTOPES]; /* qP[ct][o]: q point permutation for orientation o of integ domain */
102: /* Equations */
103: DSBoundary boundary; /* Linked list of boundary conditions */
104: PetscBool useJacPre; /* Flag for using the Jacobian preconditioner */
105: PetscBool *implicit; /* Flag for implicit or explicit solve for each field */
106: PetscInt *jetDegree; /* The highest derivative for each field equation, or the k-jet that each discretization needs to tabulate */
107: PetscWeakForm wf; /* The PetscWeakForm holding all pointwise functions */
108: PetscPointFunc *update; /* Direct update of field coefficients */
109: PetscSimplePointFunc *exactSol; /* Exact solutions for each field */
110: void **exactCtx; /* Contexts for the exact solution functions */
111: PetscSimplePointFunc *exactSol_t; /* Time derivative of the exact solutions for each field */
112: void **exactCtx_t; /* Contexts for the time derivative of the exact solution functions */
113: PetscInt numConstants; /* Number of constants passed to point functions */
114: PetscScalar *constants; /* Array of constants passed to point functions */
115: void **ctx; /* User contexts for each field */
116: /* Computed sizes */
117: PetscInt totDim; /* Total system dimension */
118: PetscInt totComp; /* Total field components */
119: PetscInt *Nc; /* Number of components for each field */
120: PetscInt *Nb; /* Number of basis functions for each field */
121: PetscInt *off; /* Offsets for each field */
122: PetscInt *offDer; /* Derivative offsets for each field */
123: PetscInt *offCohesive[3]; /* Offsets for each field on side s of a cohesive cell */
124: PetscInt *offDerCohesive[3]; /* Derivative offsets for each field on side s of a cohesive cell */
125: PetscTabulation *T; /* Basis function and derivative tabulation for each field */
126: PetscTabulation *Tf; /* Basis function and derivative tabulation for each local face and field */
127: /* Work space */
128: PetscScalar *u; /* Field evaluation */
129: PetscScalar *u_t; /* Field time derivative evaluation */
130: PetscScalar *u_x; /* Field gradient evaluation */
131: PetscScalar *basisReal; /* Workspace for pushforward */
132: PetscScalar *basisDerReal; /* Workspace for derivative pushforward */
133: PetscScalar *testReal; /* Workspace for pushforward */
134: PetscScalar *testDerReal; /* Workspace for derivative pushforward */
135: PetscReal *x; /* Workspace for computing real coordinates */
136: PetscScalar *f0, *f1; /* Point evaluations of weak form residual integrands */
137: PetscScalar *g0, *g1, *g2, *g3; /* Point evaluations of weak form Jacobian integrands */
138: };
140: typedef struct {
141: PetscInt dummy; /* */
142: } PetscDS_Basic;
144: PETSC_INTERN PetscErrorCode PetscDSGetDiscType_Internal(PetscDS, PetscInt, PetscDiscType *);