Actual source code: sectionimpl.h
1: #pragma once
3: #include <petscsection.h>
4: #include <petsc/private/petscimpl.h>
5: #include <petsc/private/hashmap.h>
7: typedef struct PetscSectionClosurePermKey {
8: PetscInt depth, size;
9: } PetscSectionClosurePermKey;
11: typedef struct {
12: PetscInt *perm, *invPerm;
13: } PetscSectionClosurePermVal;
15: static inline PetscHash_t PetscSectionClosurePermHash(PetscSectionClosurePermKey k)
16: {
17: return PetscHashCombine(PetscHashInt(k.depth), PetscHashInt(k.size));
18: }
20: static inline int PetscSectionClosurePermEqual(PetscSectionClosurePermKey k1, PetscSectionClosurePermKey k2)
21: {
22: return k1.depth == k2.depth && k1.size == k2.size;
23: }
25: static PetscSectionClosurePermVal PetscSectionClosurePermVal_Empty = {PETSC_NULLPTR, PETSC_NULLPTR};
26: PETSC_HASH_MAP(ClPerm, PetscSectionClosurePermKey, PetscSectionClosurePermVal, PetscSectionClosurePermHash, PetscSectionClosurePermEqual, PetscSectionClosurePermVal_Empty)
28: struct _p_PetscSection {
29: PETSCHEADER(int);
30: PetscInt pStart, pEnd; /* The chart: all points are contained in [pStart, pEnd) */
31: IS perm; /* A permutation of [0, pEnd-pStart), so perm[i] is the ith permuted point */
32: PetscBT blockStarts; /* If present, bit is high for each point that starts a block */
33: PetscBool pointMajor; /* True if the offsets are point major, otherwise they are fieldMajor */
34: PetscBool includesConstraints; /* True if constrained dofs are included when computing offsets */
35: PetscInt *atlasDof; /* Describes layout of storage, point --> # of values */
36: PetscInt *atlasOff; /* Describes layout of storage, point --> offset into storage */
37: PetscInt maxDof; /* Maximum dof on any point */
38: PetscSection bc; /* Describes constraints, point --> # local dofs which are constrained */
39: PetscInt *bcIndices; /* Local indices for constrained dofs */
40: PetscBool setup;
42: PetscInt numFields; /* The number of fields making up the degrees of freedom */
43: char **fieldNames; /* The field names */
44: PetscInt *numFieldComponents; /* The number of components in each field */
45: PetscSection *field; /* A section describing the layout and constraints for each field */
46: PetscBool useFieldOff; /* Use the field offsets directly for the global section, rather than the point offset */
47: char ***compNames; /* The component names */
49: PetscObject clObj; /* Key for the closure (right now we only have one) */
50: PetscClPerm clHash; /* Hash of (depth, size) to perm and invPerm */
51: PetscSection clSection; /* Section giving the number of points in each closure */
52: IS clPoints; /* Points in each closure */
53: PetscSectionSym sym; /* Symmetries of the data */
54: };
56: struct _PetscSectionSymOps {
57: PetscErrorCode (*getpoints)(PetscSectionSym, PetscSection, PetscInt, const PetscInt *, const PetscInt **, const PetscScalar **);
58: PetscErrorCode (*distribute)(PetscSectionSym, PetscSF, PetscSectionSym *);
59: PetscErrorCode (*copy)(PetscSectionSym, PetscSectionSym);
60: PetscErrorCode (*destroy)(PetscSectionSym);
61: PetscErrorCode (*view)(PetscSectionSym, PetscViewer);
62: };
64: typedef struct _n_SymWorkLink *SymWorkLink;
66: struct _n_SymWorkLink {
67: SymWorkLink next;
68: const PetscInt **perms;
69: const PetscScalar **rots;
70: PetscInt numPoints;
71: };
73: struct _p_PetscSectionSym {
74: PETSCHEADER(struct _PetscSectionSymOps);
75: void *data;
76: SymWorkLink workin;
77: SymWorkLink workout;
78: };
80: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PetscSectionCopy_Internal(PetscSection, PetscSection, PetscBT);
81: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PetscSectionSetClosurePermutation_Internal(PetscSection, PetscObject, PetscInt, PetscInt, PetscCopyMode, PetscInt *);
82: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PetscSectionGetClosureInversePermutation_Internal(PetscSection, PetscObject, PetscInt, PetscInt, const PetscInt *[]);
83: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode ISIntersect_Caching_Internal(IS, IS, IS *);
84: #if defined(PETSC_HAVE_HDF5)
85: PETSC_INTERN PetscErrorCode PetscSectionView_HDF5_Internal(PetscSection, PetscViewer);
86: PETSC_INTERN PetscErrorCode PetscSectionLoad_HDF5_Internal(PetscSection, PetscViewer);
87: #endif
89: static inline PetscErrorCode PetscSectionCheckConstraints_Private(PetscSection s)
90: {
91: PetscFunctionBegin;
92: if (!s->bc) {
93: PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s->bc));
94: PetscCall(PetscSectionSetChart(s->bc, s->pStart, s->pEnd));
95: }
96: PetscFunctionReturn(PETSC_SUCCESS);
97: }
99: /* Call this if you directly modify atlasDof so that maxDof gets recalculated on next PetscSectionGetMaxDof() */
100: static inline PetscErrorCode PetscSectionInvalidateMaxDof_Internal(PetscSection s)
101: {
102: s->maxDof = PETSC_MIN_INT;
103: return PETSC_SUCCESS;
104: }
106: #if defined(PETSC_CLANG_STATIC_ANALYZER)
107: void PetscSectionCheckValidField(PetscInt, PetscInt);
108: void PetscSectionCheckValidFieldComponent(PetscInt, PetscInt);
109: #else
110: #define PetscSectionCheckValid_(description, item, nitems) \
111: do { \
112: PetscCheck(((item) >= 0) && ((item) < (nitems)), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid " description " %" PetscInt_FMT "; not in [0, %" PetscInt_FMT ")", (item), (nitems)); \
113: } while (0)
115: #define PetscSectionCheckValidFieldComponent(comp, nfieldcomp) PetscSectionCheckValid_("section field component", comp, nfieldcomp)
117: #define PetscSectionCheckValidField(field, nfields) PetscSectionCheckValid_("field number", field, nfields)
118: #endif /* PETSC_CLANG_STATIC_ANALYZER */