Actual source code: dmpleximpl.h

petsc-main 2021-04-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>

 10: #if defined(PETSC_HAVE_EXODUSII)
 11: #include <exodusII.h>
 12: #endif


 15: PETSC_EXTERN PetscLogEvent DMPLEX_Interpolate;
 16: PETSC_EXTERN PetscLogEvent DMPLEX_Partition;
 17: PETSC_EXTERN PetscLogEvent DMPLEX_PartSelf;
 18: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelInvert;
 19: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelCreateSF;
 20: PETSC_EXTERN PetscLogEvent DMPLEX_PartStratSF;
 21: PETSC_EXTERN PetscLogEvent DMPLEX_CreatePointSF;
 22: PETSC_EXTERN PetscLogEvent DMPLEX_Distribute;
 23: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeCones;
 24: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeLabels;
 25: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeSF;
 26: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeOverlap;
 27: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeField;
 28: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeData;
 29: PETSC_EXTERN PetscLogEvent DMPLEX_Migrate;
 30: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolateSF;
 31: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalBegin;
 32: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalEnd;
 33: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalBegin;
 34: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalEnd;
 35: PETSC_EXTERN PetscLogEvent DMPLEX_Stratify;
 36: PETSC_EXTERN PetscLogEvent DMPLEX_Symmetrize;
 37: PETSC_EXTERN PetscLogEvent DMPLEX_Preallocate;
 38: PETSC_EXTERN PetscLogEvent DMPLEX_ResidualFEM;
 39: PETSC_EXTERN PetscLogEvent DMPLEX_JacobianFEM;
 40: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolatorFEM;
 41: PETSC_EXTERN PetscLogEvent DMPLEX_InjectorFEM;
 42: PETSC_EXTERN PetscLogEvent DMPLEX_IntegralFEM;
 43: PETSC_EXTERN PetscLogEvent DMPLEX_CreateGmsh;
 44: PETSC_EXTERN PetscLogEvent DMPLEX_RebalanceSharedPoints;
 45: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromFile;
 46: PETSC_EXTERN PetscLogEvent DMPLEX_BuildFromCellList;
 47: PETSC_EXTERN PetscLogEvent DMPLEX_BuildCoordinatesFromCellList;
 48: PETSC_EXTERN PetscLogEvent DMPLEX_LocatePoints;

 50: typedef struct _DMPlexCellRefinerOps *DMPlexCellRefinerOps;
 51: struct _DMPlexCellRefinerOps {
 52:   PetscErrorCode (*refine)(DMPlexCellRefiner, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
 53:   PetscErrorCode (*mapsubcells)(DMPlexCellRefiner, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
 54:   PetscErrorCode (*getaffinetransforms)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
 55:   PetscErrorCode (*getaffinefacetransforms)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]);
 56:   PetscErrorCode (*getcellvertices)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[]);
 57:   PetscErrorCode (*getsubcellvertices)(DMPlexCellRefiner, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt *, PetscInt *[]);
 58:   PetscErrorCode (*mapcoords)(DMPlexCellRefiner, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
 59:   PetscErrorCode (*setup)(DMPlexCellRefiner);
 60:   PetscErrorCode (*destroy)(DMPlexCellRefiner);
 61: };

 63: typedef struct _n_PlexGeneratorFunctionList *PlexGeneratorFunctionList;
 64: struct _n_PlexGeneratorFunctionList {
 65:   PetscErrorCode    (*generate)(DM, PetscBool, DM*);
 66:   PetscErrorCode    (*refine)(DM, PetscReal*, DM*);
 67:   PetscErrorCode    (*adaptlabel)(DM, DMLabel, DM*);
 68:   char              *name;
 69:   PetscInt          dim;
 70:   PlexGeneratorFunctionList next;
 71: };

 73: struct _p_DMPlexCellRefiner {
 74:   PETSCHEADER(struct _DMPlexCellRefinerOps);
 75:   DM                    dm;          /* The original DM */
 76:   PetscBool             setupcalled;
 77:   DMPlexCellRefinerType type;
 78:   DMLabel               refineType;  /* The refinement type of each point, since there are multiple ways to refine a cell */
 79:   PetscInt              *ctOrder;    /* [i] = ct: An array with cell types in depth order */
 80:   PetscInt              *ctOrderInv; /* [ct] = i: An array with the ordinal numbers for each cell type */
 81:   PetscInt              *ctStart;    /* [ct]: The number for the first cell of each polytope type in the original mesh */
 82:   PetscInt              *ctStartNew; /* [ctNew]: The number for the first cell of each polytope type in the new mesh */
 83:   PetscInt              *offset;     /* [ct/rt][ctNew]: The offset in the new point numbering of a point of type ctNew produced from an old point of type ct or refine type rt */
 84:   PetscFE               *coordFE;    /* Finite element for each cell type, used for localized coordinate interpolation */
 85:   PetscFEGeom           **refGeom;   /* Geometry of the reference cell for each cell type */
 86:   DMLabel               adaptLabel;  /* Optional label indicating cells to be refined */
 87:   void                  *data;       /* refiner private data */
 88: };

 90: /* Utility struct to store the contents of a Fluent file in memory */
 91: typedef struct {
 92:   int          index;    /* Type of section */
 93:   unsigned int zoneID;
 94:   unsigned int first;
 95:   unsigned int last;
 96:   int          type;
 97:   int          nd;       /* Either ND or element-type */
 98:   void        *data;
 99: } FluentSection;

101: struct _PetscGridHash {
102:   PetscInt     dim;
103:   PetscReal    lower[3];    /* The lower-left corner */
104:   PetscReal    upper[3];    /* The upper-right corner */
105:   PetscReal    extent[3];   /* The box size */
106:   PetscReal    h[3];        /* The subbox size */
107:   PetscInt     n[3];        /* The number of subboxes */
108:   PetscSection cellSection; /* Offsets for cells in each subbox*/
109:   IS           cells;       /* List of cells in each subbox */
110:   DMLabel      cellsSparse; /* Sparse storage for cell map */
111: };

113: /* Point Numbering in Plex:

115:    Points are numbered contiguously by stratum. Strate are organized as follows:

117:    First Stratum:  Cells [height 0]
118:    Second Stratum: Vertices [depth 0]
119:    Third Stratum:  Faces [height 1]
120:    Fourth Stratum: Edges [depth 1]

122:    We do this so that the numbering of a cell-vertex mesh does not change after interpolation. Within a given stratum,
123:    we allow additional segregation of by cell type.
124: */
125: typedef struct {
126:   PetscInt             refct;

128:   PetscSection         coneSection;       /* Layout of cones (inedges for DAG) */
129:   PetscInt             maxConeSize;       /* Cached for fast lookup */
130:   PetscInt            *cones;             /* Cone for each point */
131:   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 */
132:   PetscSection         supportSection;    /* Layout of cones (inedges for DAG) */
133:   PetscInt             maxSupportSize;    /* Cached for fast lookup */
134:   PetscInt            *supports;          /* Cone for each point */
135:   PetscBool            refinementUniform; /* Flag for uniform cell refinement */
136:   PetscReal            refinementLimit;   /* Maximum volume for refined cell */
137:   PetscErrorCode     (*refinementFunc)(const PetscReal [], PetscReal *); /* Function giving the maximum volume for refined cell */
138:   PetscInt             overlap;           /* Overlap of the partitions as passed to DMPlexDistribute() or DMPlexDistributeOverlap() */
139:   DMPlexInterpolatedFlag interpolated;
140:   DMPlexInterpolatedFlag interpolatedCollective;

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

144:   /* Hierarchy */
145:   DMPlexCellRefinerType cellRefiner;       /* Strategy for refining cells */
146:   PetscBool             regularRefinement; /* This flag signals that we are a regular refinement of coarseMesh */

148:   /* Generation */
149:   char                *tetgenOpts;
150:   char                *triangleOpts;
151:   PetscPartitioner     partitioner;
152:   PetscBool            partitionBalance;  /* Evenly divide partition overlap when distributing */
153:   PetscBool            remeshBd;

155:   /* Submesh */
156:   DMLabel              subpointMap;       /* Label each original mesh point in the submesh with its depth, subpoint are the implicit numbering */
157:   IS                   subpointIS;        /* IS holding point number in the enclosing mesh of every point in the submesh chart */
158:   PetscObjectState     subpointState;     /* The state of subpointMap when the subpointIS was last created */

160:   /* Labels and numbering */
161:   PetscObjectState     depthState;        /* State of depth label, so that we can determine if a user changes it */
162:   PetscObjectState     celltypeState;     /* State of celltype label, so that we can determine if a user changes it */
163:   IS                   globalVertexNumbers;
164:   IS                   globalCellNumbers;

166:   /* Constraints */
167:   PetscSection         anchorSection;      /* maps constrained points to anchor points */
168:   IS                   anchorIS;           /* anchors indexed by the above section */
169:   PetscErrorCode     (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
170:   PetscErrorCode     (*computeanchormatrix)(DM,PetscSection,PetscSection,Mat);

172:   /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
173:   PetscSection         parentSection;     /* dof == 1 if point has parent */
174:   PetscInt            *parents;           /* point to parent */
175:   PetscInt            *childIDs;          /* point to child ID */
176:   PetscSection         childSection;      /* inverse of parent section */
177:   PetscInt            *children;          /* point to children */
178:   DM                   referenceTree;     /* reference tree to which child ID's refer */
179:   PetscErrorCode      (*getchildsymmetry)(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*);

181:   /* MATIS support */
182:   PetscSection         subdomainSection;

184:   /* Adjacency */
185:   PetscBool            useAnchors;        /* Replace constrained points with their anchors in adjacency lists */
186:   PetscErrorCode      (*useradjacency)(DM,PetscInt,PetscInt*,PetscInt[],void*); /* User callback for adjacency */
187:   void                *useradjacencyctx;  /* User context for callback */

189:   /* Projection */
190:   PetscInt             maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */
191:   PetscInt             activePoint;         /* current active point in iteration */

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

197:   /* Geometry */
198:   PetscReal            minradius;         /* Minimum distance from cell centroid to face */
199:   PetscBool            useHashLocation;   /* Use grid hashing for point location */
200:   PetscGridHash        lbox;              /* Local box for searching */
201:   void               (*coordFunc)(PetscInt, PetscInt, PetscInt, /* Function used to remap newly introduced vertices */
202:                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
203:                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
204:                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);

206:   /* Neighbors */
207:   PetscMPIInt*         neighbors;

209:   /* Debugging */
210:   PetscBool            printSetValues;
211:   PetscInt             printFEM;
212:   PetscInt             printL2;
213:   PetscReal            printTol;
214: } DM_Plex;

216: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM,PetscViewer);
217: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec,PetscViewer);
218: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec,PetscViewer);
219: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec,PetscViewer);
220: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec,PetscViewer);
221: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec,PetscViewer);
222: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec,PetscViewer);
223: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
224: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM,PetscViewer);
225: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject,PetscViewer);
226: #if defined(PETSC_HAVE_HDF5)
227: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
228: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
229: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
230: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
231: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
232: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
233: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
234: #endif

236: PETSC_INTERN PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM, PetscSection, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
237: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM,PetscInt,const PetscInt[],IS*);
238: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *, DM);
239: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
240: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM []);
241: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
242: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM []);
243: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, DMLabel, DM *);
244: PETSC_INTERN PetscErrorCode DMAdaptMetric_Plex(DM, Vec, DMLabel, DM *);
245: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
246: PETSC_INTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
247: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
248: 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);
249: 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);
250: 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);
251: PETSC_INTERN PetscErrorCode DMProjectBdFieldLabelLocal_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[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
252: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
253: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *);
254: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
255: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);

257: PETSC_INTERN PetscErrorCode DMPlexLoadLabels_HDF5_Internal(DM, PetscViewer);
258: PETSC_INTERN PetscErrorCode DMPlexView_HDF5_Internal(DM, PetscViewer);
259: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Internal(DM, PetscViewer);
260: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM, PetscViewer);
261: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Internal(Vec, PetscViewer);
262: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec, PetscViewer);
263: PETSC_INTERN PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec, PetscViewer);
264: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec, PetscViewer);
265: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec, PetscViewer);
266: #if defined(PETSC_HAVE_EXODUSII)
267: PETSC_EXTERN PetscErrorCode DMView_PlexExodusII(DM, PetscViewer);
268: PETSC_INTERN PetscErrorCode VecView_PlexExodusII_Internal(Vec, PetscViewer);
269: PETSC_INTERN PetscErrorCode VecLoad_PlexExodusII_Internal(Vec, PetscViewer);
270: PETSC_INTERN PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec, int, int, int);
271: PETSC_INTERN PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec, int, int, int);
272: PETSC_INTERN PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec, int, int, int);
273: PETSC_INTERN PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec, int, int, int);
274: PETSC_INTERN PetscErrorCode EXOGetVarIndex_Internal(int, ex_entity_type, const char[], int*);
275: #endif
276: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM,PetscInt,PetscInt,PetscInt*);
277: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM,PetscInt,PetscBool,PetscBool,PetscBool,PetscInt*,PetscInt*[]);
278: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM,DMPolytopeType,const PetscInt[],PetscInt*,const DMPolytopeType*[],const PetscInt*[],const PetscInt*[]);
279: PETSC_INTERN PetscErrorCode DMPlexRestoreRawFaces_Internal(DM,DMPolytopeType,const PetscInt[],PetscInt*,const DMPolytopeType*[],const PetscInt*[],const PetscInt*[]);
280: PETSC_INTERN PetscErrorCode CellRefinerInCellTest_Internal(DMPolytopeType, const PetscReal[], PetscBool *);
281: PETSC_INTERN PetscErrorCode DMPlexCellRefinerAdaptLabel(DM, DMLabel, DM *);
282: PETSC_INTERN PetscErrorCode DMPlexComputeCellType_Internal(DM, PetscInt, PetscInt, DMPolytopeType *);
283: PETSC_INTERN PetscErrorCode DMPlexCreateCellTypeOrder_Internal(DMPolytopeType, PetscInt *[], PetscInt *[]);
284: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscScalar[], InsertMode);
285: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
286: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_SetTree(DM, PetscSection, PetscInt[], PetscInt[]);
287: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM,DM,const char *,DM*);
288: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM,DM,PetscSF,PetscInt *,Mat);
289: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM,DM,PetscSF,PetscInt *,Mat);
290: PETSC_EXTERN PetscErrorCode DMPlexAnchorsModifyMat(DM,PetscSection,PetscInt,PetscInt,const PetscInt[],const PetscInt ***,const PetscScalar[],PetscInt*,PetscInt*,PetscInt*[],PetscScalar*[],PetscInt[],PetscBool);
291: PETSC_EXTERN PetscErrorCode indicesPoint_private(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,PetscInt,PetscInt []);
292: PETSC_EXTERN PetscErrorCode indicesPointFields_private(PetscSection,PetscInt,PetscInt,PetscInt [],PetscBool,PetscInt,PetscInt []);
293: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM,PetscInt,const PetscScalar [],PetscInt,PetscInt *);
294: /* these two are PETSC_EXTERN just because of src/dm/impls/plex/tests/ex18.c */
295: PETSC_EXTERN PetscErrorCode DMPlexOrientCell_Internal(DM,PetscInt,PetscInt,PetscBool);
296: PETSC_EXTERN PetscErrorCode DMPlexOrientInterface_Internal(DM);

298: /* Applications may use this function */
299: PETSC_EXTERN PetscErrorCode DMPlexCreateNumbering_Plex(DM, PetscInt, PetscInt, PetscInt, PetscInt *, PetscSF, IS *);

301: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
302: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
303: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, DMLabel, DM *);
304: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, DMLabel, DM *);
305: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat*);

307: PETSC_INTERN PetscErrorCode DMPlexGetOverlap_Plex(DM, PetscInt *);

309: /* invert dihedral symmetry: return a^-1,
310:  * using the representation described in
311:  * DMPlexGetConeOrientation() */
312: PETSC_STATIC_INLINE PetscInt DihedralInvert(PetscInt N, PetscInt a)
313: {
314:   return (a <= 0) ? a : (N - a);
315: }

317: /* invert dihedral symmetry: return b * a,
318:  * using the representation described in
319:  * DMPlexGetConeOrientation() */
320: PETSC_STATIC_INLINE PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
321: {
322:   if (!N) return 0;
323:   return  (a >= 0) ?
324:          ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) :
325:          ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
326: }

328: /* swap dihedral symmetries: return b * a^-1,
329:  * using the representation described in
330:  * DMPlexGetConeOrientation() */
331: PETSC_STATIC_INLINE PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
332: {
333:   return DihedralCompose(N,DihedralInvert(N,a),b);
334: }

336: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, PetscHashFormKey, IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
337: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM, PetscHashFormKey[], IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
338: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, PetscHashFormKey, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
339: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
340: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);

342: /* Matvec with A in row-major storage, x and y can be aliased */
343: PETSC_STATIC_INLINE void DMPlex_Mult2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
344: {
345:   PetscScalar z[2];
346:   z[0] = x[0]; z[1] = x[ldx];
347:   y[0]   = A[0]*z[0] + A[1]*z[1];
348:   y[ldx] = A[2]*z[0] + A[3]*z[1];
349:   (void)PetscLogFlops(6.0);
350: }
351: PETSC_STATIC_INLINE void DMPlex_Mult3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
352: {
353:   PetscScalar z[3];
354:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
355:   y[0]     = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
356:   y[ldx]   = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
357:   y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
358:   (void)PetscLogFlops(15.0);
359: }
360: PETSC_STATIC_INLINE void DMPlex_MultTranspose2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
361: {
362:   PetscScalar z[2];
363:   z[0] = x[0]; z[1] = x[ldx];
364:   y[0]   = A[0]*z[0] + A[2]*z[1];
365:   y[ldx] = A[1]*z[0] + A[3]*z[1];
366:   (void)PetscLogFlops(6.0);
367: }
368: PETSC_STATIC_INLINE void DMPlex_MultTranspose3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
369: {
370:   PetscScalar z[3];
371:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
372:   y[0]     = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
373:   y[ldx]   = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
374:   y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
375:   (void)PetscLogFlops(15.0);
376: }
377: PETSC_STATIC_INLINE void DMPlex_Mult2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
378: {
379:   PetscScalar z[2];
380:   z[0] = x[0]; z[1] = x[ldx];
381:   y[0]   = A[0]*z[0] + A[1]*z[1];
382:   y[ldx] = A[2]*z[0] + A[3]*z[1];
383:   (void)PetscLogFlops(6.0);
384: }
385: PETSC_STATIC_INLINE void DMPlex_Mult3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
386: {
387:   PetscScalar z[3];
388:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
389:   y[0]     = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
390:   y[ldx]   = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
391:   y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
392:   (void)PetscLogFlops(15.0);
393: }
394: PETSC_STATIC_INLINE void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
395: {
396:   PetscScalar z[2];
397:   z[0] = x[0]; z[1] = x[ldx];
398:   y[0]   += A[0]*z[0] + A[1]*z[1];
399:   y[ldx] += A[2]*z[0] + A[3]*z[1];
400:   (void)PetscLogFlops(6.0);
401: }
402: PETSC_STATIC_INLINE void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
403: {
404:   PetscScalar z[3];
405:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
406:   y[0]     += A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
407:   y[ldx]   += A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
408:   y[ldx*2] += A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
409:   (void)PetscLogFlops(15.0);
410: }
411: PETSC_STATIC_INLINE void DMPlex_MultTransposeReal_Internal(const PetscReal A[], PetscInt m, PetscInt n, PetscInt ldx, const PetscScalar x[], PetscScalar y[])
412: {
413:   PetscScalar z[3];
414:   PetscInt    i, j;
415:   for (i = 0; i < m; ++i) z[i] = x[i*ldx];
416:   for (j = 0; j < n; ++j) {
417:     const PetscInt l = j*ldx;
418:     y[l] = 0;
419:     for (i = 0; i < m; ++i) {
420:       y[l] += A[j*n+i]*z[i];
421:     }
422:   }
423:   (void)PetscLogFlops(2*m*n);
424: }
425: PETSC_STATIC_INLINE void DMPlex_MultTranspose2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
426: {
427:   PetscScalar z[2];
428:   z[0] = x[0]; z[1] = x[ldx];
429:   y[0]   = A[0]*z[0] + A[2]*z[1];
430:   y[ldx] = A[1]*z[0] + A[3]*z[1];
431:   (void)PetscLogFlops(6.0);
432: }
433: PETSC_STATIC_INLINE void DMPlex_MultTranspose3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
434: {
435:   PetscScalar z[3];
436:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
437:   y[0]     = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
438:   y[ldx]   = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
439:   y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
440:   (void)PetscLogFlops(15.0);
441: }

443: PETSC_STATIC_INLINE void DMPlex_MatMult2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
444: {
445:   PetscInt j;
446:   for (j = 0; j < n; ++j) {
447:     PetscScalar z[2];
448:     z[0] = B[0+j]; z[1] = B[1*ldb+j];
449:     DMPlex_Mult2D_Internal(A, 1, z, z);
450:     C[0+j] = z[0]; C[1*ldb+j] = z[1];
451:   }
452:   (void)PetscLogFlops(8.0*n);
453: }
454: PETSC_STATIC_INLINE void DMPlex_MatMult3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
455: {
456:   PetscInt j;
457:   for (j = 0; j < n; ++j) {
458:     PetscScalar z[3];
459:     z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
460:     DMPlex_Mult3D_Internal(A, 1, z, z);
461:     C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
462:   }
463:   (void)PetscLogFlops(8.0*n);
464: }
465: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
466: {
467:   PetscInt j;
468:   for (j = 0; j < n; ++j) {
469:     PetscScalar z[2];
470:     z[0] = B[0+j]; z[1] = B[1*ldb+j];
471:     DMPlex_MultTranspose2D_Internal(A, 1, z, z);
472:     C[0+j] = z[0]; C[1*ldb+j] = z[1];
473:   }
474:   (void)PetscLogFlops(8.0*n);
475: }
476: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
477: {
478:   PetscInt j;
479:   for (j = 0; j < n; ++j) {
480:     PetscScalar z[3];
481:     z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
482:     DMPlex_MultTranspose3D_Internal(A, 1, z, z);
483:     C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
484:   }
485:   (void)PetscLogFlops(8.0*n);
486: }

488: PETSC_STATIC_INLINE void DMPlex_MatMultLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
489: {
490:   PetscInt j;
491:   for (j = 0; j < m; ++j) {
492:     DMPlex_MultTranspose2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
493:   }
494:   (void)PetscLogFlops(8.0*m);
495: }
496: PETSC_STATIC_INLINE void DMPlex_MatMultLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
497: {
498:   PetscInt j;
499:   for (j = 0; j < m; ++j) {
500:     DMPlex_MultTranspose3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
501:   }
502:   (void)PetscLogFlops(8.0*m);
503: }
504: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
505: {
506:   PetscInt j;
507:   for (j = 0; j < m; ++j) {
508:     DMPlex_Mult2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
509:   }
510:   (void)PetscLogFlops(8.0*m);
511: }
512: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
513: {
514:   PetscInt j;
515:   for (j = 0; j < m; ++j) {
516:     DMPlex_Mult3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
517:   }
518:   (void)PetscLogFlops(8.0*m);
519: }

521: PETSC_STATIC_INLINE void DMPlex_PTAP2DReal_Internal(const PetscReal P[], const PetscScalar A[], PetscScalar C[])
522: {
523:   PetscScalar out[4];
524:   PetscInt i, j, k, l;
525:   for (i = 0; i < 2; ++i) {
526:     for (j = 0; j < 2; ++j) {
527:       out[i*2+j] = 0.;
528:       for (k = 0; k < 2; ++k) {
529:         for (l = 0; l < 2; ++l) {
530:           out[i*2+j] += P[k*2+i]*A[k*2+l]*P[l*2+j];
531:         }
532:       }
533:     }
534:   }
535:   for (i = 0; i < 2*2; ++i) C[i] = out[i];
536:   (void)PetscLogFlops(48.0);
537: }
538: PETSC_STATIC_INLINE void DMPlex_PTAP3DReal_Internal(const PetscReal P[], const PetscScalar A[], PetscScalar C[])
539: {
540:   PetscScalar out[9];
541:   PetscInt i, j, k, l;
542:   for (i = 0; i < 3; ++i) {
543:     for (j = 0; j < 3; ++j) {
544:       out[i*3+j] = 0.;
545:       for (k = 0; k < 3; ++k) {
546:         for (l = 0; l < 3; ++l) {
547:           out[i*3+j] += P[k*3+i]*A[k*3+l]*P[l*3+j];
548:         }
549:       }
550:     }
551:   }
552:   for (i = 0; i < 2*2; ++i) C[i] = out[i];
553:   (void)PetscLogFlops(243.0);
554: }
555: /* TODO Fix for aliasing of A and C */
556: PETSC_STATIC_INLINE void DMPlex_PTAPReal_Internal(const PetscReal P[], PetscInt m, PetscInt n, const PetscScalar A[], PetscScalar C[])
557: {
558:   PetscInt i, j, k, l;
559:   for (i = 0; i < n; ++i) {
560:     for (j = 0; j < n; ++j) {
561:       C[i*n+j] = 0.;
562:       for (k = 0; k < m; ++k) {
563:         for (l = 0; l < m; ++l) {
564:           C[i*n+j] += P[k*n+i]*A[k*m+l]*P[l*n+j];
565:         }
566:       }
567:     }
568:   }
569:   (void)PetscLogFlops(243.0);
570: }

572: PETSC_STATIC_INLINE void DMPlex_Transpose2D_Internal(PetscScalar A[])
573: {
574:   PetscScalar tmp;
575:   tmp = A[1]; A[1] = A[2]; A[2] = tmp;
576: }
577: PETSC_STATIC_INLINE void DMPlex_Transpose3D_Internal(PetscScalar A[])
578: {
579:   PetscScalar tmp;
580:   tmp = A[1]; A[1] = A[3]; A[3] = tmp;
581:   tmp = A[2]; A[2] = A[6]; A[6] = tmp;
582:   tmp = A[5]; A[5] = A[7]; A[7] = tmp;
583: }

585: PETSC_STATIC_INLINE void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
586: {
587:   const PetscReal invDet = 1.0/detJ;

589:   invJ[0] =  invDet*J[3];
590:   invJ[1] = -invDet*J[1];
591:   invJ[2] = -invDet*J[2];
592:   invJ[3] =  invDet*J[0];
593:   (void)PetscLogFlops(5.0);
594: }

596: PETSC_STATIC_INLINE void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
597: {
598:   const PetscReal invDet = 1.0/detJ;

600:   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
601:   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
602:   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
603:   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
604:   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
605:   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
606:   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
607:   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
608:   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
609:   (void)PetscLogFlops(37.0);
610: }

612: PETSC_STATIC_INLINE void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
613: {
614:   *detJ = J[0]*J[3] - J[1]*J[2];
615:   (void)PetscLogFlops(3.0);
616: }

618: PETSC_STATIC_INLINE void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
619: {
620:   *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
621:            J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
622:            J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
623:   (void)PetscLogFlops(12.0);
624: }

626: PETSC_STATIC_INLINE void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
627: {
628:   *detJ = PetscRealPart(J[0])*PetscRealPart(J[3]) - PetscRealPart(J[1])*PetscRealPart(J[2]);
629:   (void)PetscLogFlops(3.0);
630: }

632: PETSC_STATIC_INLINE void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
633: {
634:   *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])) +
635:            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])) +
636:            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])));
637:   (void)PetscLogFlops(12.0);
638: }

640: 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];}

642: 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;}

644: 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;}

646: 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);}

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

650: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Translate_Private(PetscInt ornt, PetscInt *start, PetscBool *reverse)
651: {
653:   *reverse = (ornt < 0) ? PETSC_TRUE : PETSC_FALSE;
654:   *start = *reverse ? -(ornt+1) : ornt;
655:   return(0);
656: }

658: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Combine_Private(PetscInt coneSize, PetscInt origStart, PetscBool origReverse, PetscInt rotateStart, PetscBool rotateReverse, PetscInt *newStart, PetscBool *newReverse)
659: {
661:   *newReverse = (origReverse == rotateReverse) ? PETSC_FALSE : PETSC_TRUE;
662:   *newStart = rotateReverse ? (coneSize + rotateStart - origStart) : (coneSize + origStart - rotateStart);
663:   *newStart %= coneSize;
664:   return(0);
665: }

667: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_TranslateBack_Private(PetscInt coneSize, PetscInt start, PetscBool reverse, PetscInt *ornt)
668: {
670:   if (coneSize < 3) {
671:     /* edges just get flipped if start == 1 regardless direction */
672:     *ornt = start ? -2 : 0;
673:   } else {
674:     *ornt = reverse ? -(start+1) : start;
675:   }
676:   return(0);
677: }

679: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Permute_Private(PetscInt n, const PetscInt arr[], PetscInt start, PetscBool reverse, PetscInt newarr[])
680: {
681:   PetscInt i;

684:   if (reverse) {for (i=0; i<n; i++) newarr[i] = arr[(n+start-i)%n];}
685:   else         {for (i=0; i<n; i++) newarr[i] = arr[(start+i)%n];}
686:   return(0);
687: }

689: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM,PetscInt,PetscInt,PetscDualSpace *);
690: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection,PetscBool,PetscInt,PetscInt,PetscInt *,PetscBool,const PetscInt[],const PetscInt[],PetscInt[]);
691: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection,PetscBool,PetscInt,PetscInt,PetscInt[],PetscBool,const PetscInt***,PetscInt,const PetscInt[],PetscInt[]);
692: PETSC_INTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
693: PETSC_INTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);

695: PETSC_EXTERN PetscErrorCode DMPlexGetAllCells_Internal(DM, IS *);
696: PETSC_EXTERN PetscErrorCode DMSNESGetFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
697: PETSC_EXTERN PetscErrorCode DMSNESRestoreFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
698: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM, PetscSection, IS, PetscReal, Vec, Vec, Vec, void *);
699: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM, PetscSection, PetscSection, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
700: PETSC_INTERN PetscErrorCode DMCreateSubDomainDM_Plex(DM,DMLabel,PetscInt,IS*,DM*);
701: PETSC_INTERN PetscErrorCode DMPlexBasisTransformPoint_Internal(DM, DM, Vec, PetscInt, PetscBool[], PetscBool, PetscScalar *);
702: PETSC_EXTERN PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM, DM, Vec, PetscInt, PetscBool, PetscInt, PetscScalar *);
703: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscReal *, PetscReal *, void *);
704: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApply_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscScalar *, PetscScalar *, void *);
705: PETSC_INTERN PetscErrorCode DMCreateNeumannOverlap_Plex(DM, IS*, Mat*, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void **);

707: #endif /* _PLEXIMPL_H */