Actual source code: sectionimpl.h

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: #if !defined(PETSCSECTIONIMPL_H)
  2: #define PETSCSECTIONIMPL_H

  4: #include <petscsection.h>
  5: #include <petsc/private/petscimpl.h>
  6: #include <petsc/private/hashmap.h>

  8: typedef struct PetscSectionClosurePermKey {
  9:   PetscInt    depth, size;
 10: } PetscSectionClosurePermKey;
 11: typedef struct {
 12:   PetscInt *perm, *invPerm;
 13: } PetscSectionClosurePermVal;
 14: PETSC_STATIC_INLINE PetscHash_t PetscSectionClosurePermHash(PetscSectionClosurePermKey k) {
 15:   return PetscHashCombine(PetscHashInt(k.depth), PetscHashInt(k.size));
 16: }
 17: PETSC_STATIC_INLINE int PetscSectionClosurePermEqual(PetscSectionClosurePermKey k1, PetscSectionClosurePermKey k2) {
 18:   return k1.depth == k2.depth && k1.size == k2.size;
 19: }
 20: static PetscSectionClosurePermVal PetscSectionClosurePermVal_Empty = {NULL, NULL};
 21: PETSC_HASH_MAP(ClPerm, PetscSectionClosurePermKey, PetscSectionClosurePermVal, PetscSectionClosurePermHash, PetscSectionClosurePermEqual, PetscSectionClosurePermVal_Empty)

 23: struct _p_PetscSection {
 24:   PETSCHEADER(int);
 25:   PetscInt                      pStart, pEnd; /* The chart: all points are contained in [pStart, pEnd) */
 26:   IS                            perm;         /* A permutation of [0, pEnd-pStart) */
 27:   PetscBool                     pointMajor;   /* True if the offsets are point major, otherwise they are fieldMajor */
 28:   PetscInt                     *atlasDof;     /* Describes layout of storage, point --> # of values */
 29:   PetscInt                     *atlasOff;     /* Describes layout of storage, point --> offset into storage */
 30:   PetscInt                      maxDof;       /* Maximum dof on any point */
 31:   PetscSection                  bc;           /* Describes constraints, point --> # local dofs which are constrained */
 32:   PetscInt                     *bcIndices;    /* Local indices for constrained dofs */
 33:   PetscBool                     setup;

 35:   PetscInt                      numFields;    /* The number of fields making up the degrees of freedom */
 36:   char                        **fieldNames;   /* The field names */
 37:   PetscInt                     *numFieldComponents; /* The number of components in each field */
 38:   PetscSection                 *field;        /* A section describing the layout and constraints for each field */
 39:   PetscBool                     useFieldOff;  /* Use the field offsets directly for the global section, rather than the point offset */
 40:   char                        ***compNames;   /* The component names */

 42:   PetscObject                   clObj;        /* Key for the closure (right now we only have one) */
 43:   PetscClPerm                   clHash;       /* Hash of (depth, size) to perm and invPerm */
 44:   PetscSection                  clSection;    /* Section giving the number of points in each closure */
 45:   IS                            clPoints;     /* Points in each closure */
 46:   PetscSectionSym               sym;          /* Symmetries of the data */
 47: };

 49: struct _PetscSectionSymOps {
 50:   PetscErrorCode (*getpoints)(PetscSectionSym,PetscSection,PetscInt,const PetscInt *,const PetscInt **,const PetscScalar **);
 51:   PetscErrorCode (*destroy)(PetscSectionSym);
 52:   PetscErrorCode (*view)(PetscSectionSym,PetscViewer);
 53: };

 55: typedef struct _n_SymWorkLink *SymWorkLink;

 57: struct _n_SymWorkLink
 58: {
 59:   SymWorkLink         next;
 60:   const PetscInt    **perms;
 61:   const PetscScalar **rots;
 62:   PetscInt           numPoints;
 63: };

 65: struct _p_PetscSectionSym {
 66:   PETSCHEADER(struct _PetscSectionSymOps);
 67:   void *data;
 68:   SymWorkLink workin;
 69:   SymWorkLink workout;
 70: };

 72: PETSC_EXTERN PetscErrorCode PetscSectionSetClosurePermutation_Internal(PetscSection, PetscObject, PetscInt, PetscInt, PetscCopyMode, PetscInt *);
 73: PETSC_EXTERN PetscErrorCode PetscSectionGetClosurePermutation_Internal(PetscSection, PetscObject, PetscInt, PetscInt, const PetscInt *[]);
 74: PETSC_EXTERN PetscErrorCode PetscSectionGetClosureInversePermutation_Internal(PetscSection, PetscObject, PetscInt, PetscInt, const PetscInt *[]);
 75: PETSC_EXTERN PetscErrorCode ISIntersect_Caching_Internal(IS, IS, IS *);

 77: #endif