Actual source code: dmpleximpl.h

petsc-3.9.2 2018-05-20
Report Typos and Errors
  1: #if !defined(_PLEXIMPL_H)
  2: #define _PLEXIMPL_H

  4:  #include <petscmat.h>
  5:  #include <petscdmplex.h>
  6:  #include <petscbt.h>
  7:  #include <petscsf.h>
  8:  #include <petsc/private/dmimpl.h>
  9:  #include <petsc/private/isimpl.h>
 10:  #include <petsc/private/hash.h>

 12: PETSC_EXTERN PetscLogEvent DMPLEX_Interpolate;
 13: PETSC_EXTERN PetscLogEvent PETSCPARTITIONER_Partition;
 14: PETSC_EXTERN PetscLogEvent DMPLEX_Distribute;
 15: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeCones;
 16: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeLabels;
 17: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeSF;
 18: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeOverlap;
 19: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeField;
 20: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeData;
 21: PETSC_EXTERN PetscLogEvent DMPLEX_Migrate;
 22: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolateSF;
 23: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalBegin;
 24: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalEnd;
 25: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalBegin;
 26: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalEnd;
 27: PETSC_EXTERN PetscLogEvent DMPLEX_Stratify;
 28: PETSC_EXTERN PetscLogEvent DMPLEX_Preallocate;
 29: PETSC_EXTERN PetscLogEvent DMPLEX_ResidualFEM;
 30: PETSC_EXTERN PetscLogEvent DMPLEX_JacobianFEM;
 31: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolatorFEM;
 32: PETSC_EXTERN PetscLogEvent DMPLEX_InjectorFEM;
 33: PETSC_EXTERN PetscLogEvent DMPLEX_IntegralFEM;
 34: PETSC_EXTERN PetscLogEvent DMPLEX_CreateGmsh;

 36: PETSC_EXTERN PetscBool      PetscPartitionerRegisterAllCalled;
 37: PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterAll(void);
 38: PETSC_INTERN PetscErrorCode PetscPartitionerSetTypeFromOptions_Internal(PetscPartitioner);

 40: typedef enum {REFINER_NOOP = 0,
 41:               REFINER_SIMPLEX_1D,
 42:               REFINER_SIMPLEX_2D,
 43:               REFINER_HYBRID_SIMPLEX_2D,
 44:               REFINER_SIMPLEX_TO_HEX_2D,
 45:               REFINER_HEX_2D,
 46:               REFINER_HYBRID_HEX_2D,
 47:               REFINER_SIMPLEX_3D,
 48:               REFINER_HYBRID_SIMPLEX_3D,
 49:               REFINER_SIMPLEX_TO_HEX_3D,
 50:               REFINER_HEX_3D,
 51:               REFINER_HYBRID_HEX_3D} CellRefiner;

 53: typedef struct _PetscPartitionerOps *PetscPartitionerOps;
 54: struct _PetscPartitionerOps {
 55:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,PetscPartitioner);
 56:   PetscErrorCode (*setup)(PetscPartitioner);
 57:   PetscErrorCode (*view)(PetscPartitioner,PetscViewer);
 58:   PetscErrorCode (*destroy)(PetscPartitioner);
 59:   PetscErrorCode (*partition)(PetscPartitioner, DM, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscSection, IS *);
 60: };

 62: struct _p_PetscPartitioner {
 63:   PETSCHEADER(struct _PetscPartitionerOps);
 64:   void           *data;             /* Implementation object */
 65:   PetscInt        height;           /* Height of points to partition into non-overlapping subsets */
 66: };

 68: typedef struct {
 69:   PetscInt dummy;
 70: } PetscPartitioner_Chaco;

 72: typedef struct {
 73:   PetscInt dummy;
 74: } PetscPartitioner_ParMetis;

 76: typedef struct {
 77:   PetscSection section;   /* Sizes for each partition */
 78:   IS           partition; /* Points in each partition */
 79:   PetscBool    random;    /* Flag for a random partition */
 80: } PetscPartitioner_Shell;

 82: typedef struct {
 83:   PetscInt dummy;
 84: } PetscPartitioner_Simple;

 86: typedef struct {
 87:   PetscInt dummy;
 88: } PetscPartitioner_Gather;

 90: /* Utility struct to store the contents of a Gmsh file in memory */
 91: typedef struct {
 92:   PetscInt dim;      /* Entity dimension */
 93:   PetscInt id;       /* Element number */
 94:   PetscInt numNodes; /* Size of node array */
 95:   int nodes[8];      /* Node array */
 96:   PetscInt numTags;  /* Size of tag array */
 97:   int tags[4];       /* Tag array */
 98: } GmshElement;

100: /* Utility struct to store the contents of a Fluent file in memory */
101: typedef struct {
102:   int          index;    /* Type of section */
103:   unsigned int zoneID;
104:   unsigned int first;
105:   unsigned int last;
106:   int          type;
107:   int          nd;       /* Either ND or element-type */
108:   void        *data;
109: } FluentSection;

111: struct _PetscGridHash {
112:   PetscInt     dim;
113:   PetscReal    lower[3];    /* The lower-left corner */
114:   PetscReal    upper[3];    /* The upper-right corner */
115:   PetscReal    extent[3];   /* The box size */
116:   PetscReal    h[3];        /* The subbox size */
117:   PetscInt     n[3];        /* The number of subboxes */
118:   PetscSection cellSection; /* Offsets for cells in each subbox*/
119:   IS           cells;       /* List of cells in each subbox */
120:   DMLabel      cellsSparse; /* Sparse storage for cell map */
121: };

123: typedef struct {
124:   PetscInt             refct;

126:   /* Sieve */
127:   PetscSection         coneSection;       /* Layout of cones (inedges for DAG) */
128:   PetscInt             maxConeSize;       /* Cached for fast lookup */
129:   PetscInt            *cones;             /* Cone for each point */
130:   PetscInt            *coneOrientations;  /* Orientation of each cone point, means cone traveral should start on point 'o', and if negative start on -(o+1) and go in reverse */
131:   PetscSection         supportSection;    /* Layout of cones (inedges for DAG) */
132:   PetscInt             maxSupportSize;    /* Cached for fast lookup */
133:   PetscInt            *supports;          /* Cone for each point */
134:   PetscBool            refinementUniform; /* Flag for uniform cell refinement */
135:   PetscReal            refinementLimit;   /* Maximum volume for refined cell */
136:   PetscErrorCode     (*refinementFunc)(const PetscReal [], PetscReal *); /* Function giving the maximum volume for refined cell */
137:   PetscInt             hybridPointMax[8]; /* Allow segregation of some points, each dimension has a divider (used in VTK output and refinement) */

139:   PetscInt            *facesTmp;          /* Work space for faces operation */

141:   /* Hierarchy */
142:   PetscBool            regularRefinement; /* This flag signals that we are a regular refinement of coarseMesh */

144:   /* Generation */
145:   char                *tetgenOpts;
146:   char                *triangleOpts;
147:   PetscPartitioner     partitioner;
148:   PetscBool            remeshBd;

150:   /* Submesh */
151:   DMLabel              subpointMap;       /* Label each original mesh point in the submesh with its depth, subpoint are the implicit numbering */

153:   /* Labels and numbering */
154:   PetscObjectState     depthState;        /* State of depth label, so that we can determine if a user changes it */
155:   IS                   globalVertexNumbers;
156:   IS                   globalCellNumbers;

158:   /* Constraints */
159:   PetscSection         anchorSection;      /* maps constrained points to anchor points */
160:   IS                   anchorIS;           /* anchors indexed by the above section */
161:   PetscErrorCode     (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
162:   PetscErrorCode     (*computeanchormatrix)(DM,PetscSection,PetscSection,Mat);

164:   /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
165:   PetscSection         parentSection;     /* dof == 1 if point has parent */
166:   PetscInt            *parents;           /* point to parent */
167:   PetscInt            *childIDs;          /* point to child ID */
168:   PetscSection         childSection;      /* inverse of parent section */
169:   PetscInt            *children;          /* point to children */
170:   DM                   referenceTree;     /* reference tree to which child ID's refer */
171:   PetscErrorCode      (*getchildsymmetry)(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*);

173:   /* MATIS support */
174:   PetscSection         subdomainSection;

176:   /* Adjacency */
177:   PetscBool            useAnchors;        /* Replace constrained points with their anchors in adjacency lists */
178:   PetscErrorCode      (*useradjacency)(DM,PetscInt,PetscInt*,PetscInt[],void*); /* User callback for adjacency */
179:   void                *useradjacencyctx;  /* User context for callback */

181:   /* Projection */
182:   PetscInt             maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */

184:   /* Output */
185:   PetscInt             vtkCellHeight;            /* The height of cells for output, default is 0 */
186:   PetscReal            scale[NUM_PETSC_UNITS];   /* The scale for each SI unit */

188:   /* Geometry */
189:   PetscReal            minradius;         /* Minimum distance from cell centroid to face */
190:   PetscBool            useHashLocation;   /* Use grid hashing for point location */
191:   PetscGridHash        lbox;              /* Local box for searching */

193:   /* Debugging */
194:   PetscBool            printSetValues;
195:   PetscInt             printFEM;
196:   PetscReal            printTol;
197: } DM_Plex;

199: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM,PetscViewer);
200: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec,PetscViewer);
201: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec,PetscViewer);
202: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec,PetscViewer);
203: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec,PetscViewer);
204: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec,PetscViewer);
205: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec,PetscViewer);
206: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
207: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM,PetscViewer);
208: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject,PetscViewer);
209: #if defined(PETSC_HAVE_HDF5)
210: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
211: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
212: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
213: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
214: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
215: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
216: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
217: #endif

219: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *, DM);
220: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
221: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM []);
222: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
223: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM []);
224: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, DMLabel, DM *);
225: PETSC_INTERN PetscErrorCode DMAdaptMetric_Plex(DM, Vec, DMLabel, DM *);
226: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
227: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
228: PETSC_INTERN PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
229: PETSC_INTERN PetscErrorCode DMProjectFieldLocal_Plex(DM,PetscReal,Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
230: PETSC_INTERN PetscErrorCode DMProjectFieldLabelLocal_Plex(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
231: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
232: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *);
233: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
234: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);

236: PETSC_INTERN PetscErrorCode DMPlexView_HDF5_Internal(DM, PetscViewer);
237: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Internal(DM, PetscViewer);
238: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Internal(Vec, PetscViewer);
239: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec, PetscViewer);
240: PETSC_INTERN PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec, PetscViewer);
241: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec, PetscViewer);
242: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec, PetscViewer);
243: /* TODO Make these INTERN */
244: PETSC_EXTERN PetscErrorCode DMPlexView_ExodusII_Internal(DM, int, PetscInt);
245: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec, int, int);
246: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec, int, int);
247: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec, int, int);
248: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec, int, int);
249: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM,PetscInt,PetscInt,PetscInt*);
250: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM,PetscInt,PetscBool,PetscBool,PetscBool,PetscInt*,PetscInt*[]);
251: PETSC_INTERN PetscErrorCode DMPlexGetFaces_Internal(DM,PetscInt,PetscInt,PetscInt*,PetscInt*,const PetscInt*[]);
252: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM,PetscInt,PetscInt,const PetscInt[], PetscInt*,PetscInt*,const PetscInt*[]);
253: PETSC_INTERN PetscErrorCode DMPlexRestoreFaces_Internal(DM,PetscInt,PetscInt,PetscInt*,PetscInt*,const PetscInt*[]);
254: PETSC_INTERN PetscErrorCode DMPlexRefineUniform_Internal(DM,CellRefiner,DM*);
255: PETSC_INTERN PetscErrorCode DMPlexGetCellRefiner_Internal(DM,CellRefiner*);
256: PETSC_INTERN PetscErrorCode CellRefinerGetAffineTransforms_Internal(CellRefiner, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
257: PETSC_INTERN PetscErrorCode CellRefinerRestoreAffineTransforms_Internal(CellRefiner, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
258: PETSC_INTERN PetscErrorCode CellRefinerInCellTest_Internal(CellRefiner, const PetscReal[], PetscBool *);
259: PETSC_INTERN PetscErrorCode DMPlexInvertCell_Internal(PetscInt, PetscInt, PetscInt[]);
260: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
261: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
262: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM,DM,const char *,DM*);
263: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM,DM,PetscSF,PetscInt *,Mat);
264: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM,DM,PetscSF,PetscInt *,Mat);
265: PETSC_EXTERN PetscErrorCode DMPlexAnchorsModifyMat(DM,PetscSection,PetscInt,PetscInt,const PetscInt[],const PetscInt ***,const PetscScalar[],PetscInt*,PetscInt*,PetscInt*[],PetscScalar*[],PetscInt[],PetscBool);
266: PETSC_EXTERN PetscErrorCode indicesPoint_private(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,PetscInt,PetscInt []);
267: PETSC_EXTERN PetscErrorCode indicesPointFields_private(PetscSection,PetscInt,PetscInt,PetscInt [],PetscBool,PetscInt,PetscInt []);
268: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM,PetscInt,const PetscScalar [],PetscInt,PetscInt *);

270: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
271: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
272: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, DMLabel, DM *);
273: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, DMLabel, DM *);

275: /* invert dihedral symmetry: return a^-1,
276:  * using the representation described in
277:  * DMPlexGetConeOrientation() */
278: PETSC_STATIC_INLINE PetscInt DihedralInvert(PetscInt N, PetscInt a)
279: {
280:   return (a <= 0) ? a : (N - a);
281: }

283: /* invert dihedral symmetry: return b * a,
284:  * using the representation described in
285:  * DMPlexGetConeOrientation() */
286: PETSC_STATIC_INLINE PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
287: {
288:   if (!N) return 0;
289:   return  (a >= 0) ?
290:          ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) :
291:          ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
292: }

294: /* swap dihedral symmetries: return b * a^-1,
295:  * using the representation described in
296:  * DMPlexGetConeOrientation() */
297: PETSC_STATIC_INLINE PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
298: {
299:   return DihedralCompose(N,DihedralInvert(N,a),b);
300: }

302: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, PetscInt, PetscInt, PetscReal, Vec, Vec, PetscReal, Vec, void *);
303: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, PetscInt, PetscInt, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
304: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);

306: PETSC_STATIC_INLINE void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
307: {
308:   const PetscReal invDet = 1.0/detJ;

310:   invJ[0] =  invDet*J[3];
311:   invJ[1] = -invDet*J[1];
312:   invJ[2] = -invDet*J[2];
313:   invJ[3] =  invDet*J[0];
314:   (void)PetscLogFlops(5.0);
315: }

317: PETSC_STATIC_INLINE void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
318: {
319:   const PetscReal invDet = 1.0/detJ;

321:   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
322:   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
323:   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
324:   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
325:   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
326:   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
327:   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
328:   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
329:   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
330:   (void)PetscLogFlops(37.0);
331: }

333: PETSC_STATIC_INLINE void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
334: {
335:   *detJ = J[0]*J[3] - J[1]*J[2];
336:   (void)PetscLogFlops(3.0);
337: }

339: PETSC_STATIC_INLINE void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
340: {
341:   *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
342:            J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
343:            J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
344:   (void)PetscLogFlops(12.0);
345: }

347: PETSC_STATIC_INLINE void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
348: {
349:   *detJ = PetscRealPart(J[0])*PetscRealPart(J[3]) - PetscRealPart(J[1])*PetscRealPart(J[2]);
350:   (void)PetscLogFlops(3.0);
351: }

353: PETSC_STATIC_INLINE void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
354: {
355:   *detJ = (PetscRealPart(J[0*3+0])*(PetscRealPart(J[1*3+1])*PetscRealPart(J[2*3+2]) - PetscRealPart(J[1*3+2])*PetscRealPart(J[2*3+1])) +
356:            PetscRealPart(J[0*3+1])*(PetscRealPart(J[1*3+2])*PetscRealPart(J[2*3+0]) - PetscRealPart(J[1*3+0])*PetscRealPart(J[2*3+2])) +
357:            PetscRealPart(J[0*3+2])*(PetscRealPart(J[1*3+0])*PetscRealPart(J[2*3+1]) - PetscRealPart(J[1*3+1])*PetscRealPart(J[2*3+0])));
358:   (void)PetscLogFlops(12.0);
359: }

361: PETSC_STATIC_INLINE void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w) {PetscInt d; for (d = 0; d < dim; ++d) w[d] = a*x[d] + y[d];}

363: PETSC_STATIC_INLINE PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d])*y[d]; return sum;}

365: PETSC_STATIC_INLINE PetscReal DMPlex_DotRealD_Internal(PetscInt dim, const PetscReal *x, const PetscReal *y) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += x[d]*y[d]; return sum;}

367: PETSC_STATIC_INLINE PetscReal DMPlex_NormD_Internal(PetscInt dim, const PetscReal *x) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += x[d]*x[d]; return PetscSqrtReal(sum);}

369: PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
370: {
372: #if defined(PETSC_USE_DEBUG)
373:   {
374:     PetscInt       dof;
376:     *start = *end = 0; /* Silence overzealous compiler warning */
377:     if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
378:     PetscSectionGetOffset(dm->defaultSection, point, start);
379:     PetscSectionGetDof(dm->defaultSection, point, &dof);
380:     *end = *start + dof;
381:   }
382: #else
383:   {
384:     const PetscSection s = dm->defaultSection;
385:     *start = s->atlasOff[point - s->pStart];
386:     *end   = *start + s->atlasDof[point - s->pStart];
387:   }
388: #endif
389:   return(0);
390: }

392: PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
393: {
395: #if defined(PETSC_USE_DEBUG)
396:   {
397:     PetscInt       dof;
399:     *start = *end = 0; /* Silence overzealous compiler warning */
400:     if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
401:     PetscSectionGetFieldOffset(dm->defaultSection, point, field, start);
402:     PetscSectionGetFieldDof(dm->defaultSection, point, field, &dof);
403:     *end = *start + dof;
404:   }
405: #else
406:   {
407:     const PetscSection s = dm->defaultSection->field[field];
408:     *start = s->atlasOff[point - s->pStart];
409:     *end   = *start + s->atlasDof[point - s->pStart];
410:   }
411: #endif
412:   return(0);
413: }

415: PETSC_STATIC_INLINE PetscErrorCode DMPlexGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
416: {
418: #if defined(PETSC_USE_DEBUG)
419:   {
421:     PetscInt       dof,cdof;
422:     *start = *end = 0; /* Silence overzealous compiler warning */
423:     if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
424:     if (!dm->defaultGlobalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default global Section. It will be created automatically by DMGetDefaultGlobalSection()");
425:     PetscSectionGetOffset(dm->defaultGlobalSection, point, start);
426:     PetscSectionGetDof(dm->defaultGlobalSection, point, &dof);
427:     PetscSectionGetConstraintDof(dm->defaultGlobalSection, point, &cdof);
428:     *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
429:   }
430: #else
431:   {
432:     const PetscSection s    = dm->defaultGlobalSection;
433:     const PetscInt     dof  = s->atlasDof[point - s->pStart];
434:     const PetscInt     cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
435:     *start = s->atlasOff[point - s->pStart];
436:     *end   = *start + dof - cdof + (dof < 0 ? 1 : 0);
437:   }
438: #endif
439:   return(0);
440: }

442: PETSC_STATIC_INLINE PetscErrorCode DMPlexGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
443: {
445: #if defined(PETSC_USE_DEBUG)
446:   {
447:     PetscInt       loff, lfoff, fdof, fcdof, ffcdof, f;
449:     *start = *end = 0; /* Silence overzealous compiler warning */
450:     if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
451:     if (!dm->defaultGlobalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default global Section. It will be crated automatically by DMGetDefaultGlobalSection()");
452:     PetscSectionGetOffset(dm->defaultGlobalSection, point, start);
453:     PetscSectionGetOffset(dm->defaultSection, point, &loff);
454:     PetscSectionGetFieldOffset(dm->defaultSection, point, field, &lfoff);
455:     PetscSectionGetFieldDof(dm->defaultSection, point, field, &fdof);
456:     PetscSectionGetFieldConstraintDof(dm->defaultSection, point, field, &fcdof);
457:     *start = *start < 0 ? *start - (lfoff-loff) : *start + lfoff-loff;
458:     for (f = 0; f < field; ++f) {
459:       PetscSectionGetFieldConstraintDof(dm->defaultSection, point, f, &ffcdof);
460:       *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
461:     }
462:     *end   = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
463:   }
464: #else
465:   {
466:     const PetscSection s     = dm->defaultSection;
467:     const PetscSection fs    = dm->defaultSection->field[field];
468:     const PetscSection gs    = dm->defaultGlobalSection;
469:     const PetscInt     loff  = s->atlasOff[point - s->pStart];
470:     const PetscInt     goff  = gs->atlasOff[point - s->pStart];
471:     const PetscInt     lfoff = fs->atlasOff[point - s->pStart];
472:     const PetscInt     fdof  = fs->atlasDof[point - s->pStart];
473:     const PetscInt     fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
474:     PetscInt           ffcdof = 0, f;

476:     for (f = 0; f < field; ++f) {
477:       const PetscSection ffs = dm->defaultSection->field[f];
478:       ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
479:     }
480:     *start = goff + (goff < 0 ? loff-lfoff + ffcdof : lfoff-loff - ffcdof);
481:     *end   = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
482:   }
483: #endif
484:   return(0);
485: }

487: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM,PetscInt,PetscInt,PetscDualSpace *);
488: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,const PetscInt[],PetscInt[]);
489: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection,PetscInt,PetscInt,PetscInt[],PetscBool,const PetscInt***,PetscInt,PetscInt[]);

491: #endif /* _PLEXIMPL_H */