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