Actual source code: dmpleximpl.h

petsc-master 2020-07-05
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: PETSC_EXTERN PetscLogEvent DMPLEX_Interpolate;
 11: PETSC_EXTERN PetscLogEvent DMPLEX_Partition;
 12: PETSC_EXTERN PetscLogEvent DMPLEX_PartSelf;
 13: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelInvert;
 14: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelCreateSF;
 15: PETSC_EXTERN PetscLogEvent DMPLEX_PartStratSF;
 16: PETSC_EXTERN PetscLogEvent DMPLEX_CreatePointSF;
 17: PETSC_EXTERN PetscLogEvent DMPLEX_Distribute;
 18: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeCones;
 19: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeLabels;
 20: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeSF;
 21: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeOverlap;
 22: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeField;
 23: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeData;
 24: PETSC_EXTERN PetscLogEvent DMPLEX_Migrate;
 25: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolateSF;
 26: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalBegin;
 27: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalEnd;
 28: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalBegin;
 29: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalEnd;
 30: PETSC_EXTERN PetscLogEvent DMPLEX_Stratify;
 31: PETSC_EXTERN PetscLogEvent DMPLEX_Symmetrize;
 32: PETSC_EXTERN PetscLogEvent DMPLEX_Preallocate;
 33: PETSC_EXTERN PetscLogEvent DMPLEX_ResidualFEM;
 34: PETSC_EXTERN PetscLogEvent DMPLEX_JacobianFEM;
 35: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolatorFEM;
 36: PETSC_EXTERN PetscLogEvent DMPLEX_InjectorFEM;
 37: PETSC_EXTERN PetscLogEvent DMPLEX_IntegralFEM;
 38: PETSC_EXTERN PetscLogEvent DMPLEX_CreateGmsh;
 39: PETSC_EXTERN PetscLogEvent DMPLEX_RebalanceSharedPoints;
 40: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromFile;
 41: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromCellList;
 42: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromCellList_Coordinates;

 44: PETSC_EXTERN PetscBool      PetscPartitionerRegisterAllCalled;
 45: PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterAll(void);

 47: typedef struct _PetscPartitionerOps *PetscPartitionerOps;
 48: struct _PetscPartitionerOps {
 49:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,PetscPartitioner);
 50:   PetscErrorCode (*setup)(PetscPartitioner);
 51:   PetscErrorCode (*view)(PetscPartitioner,PetscViewer);
 52:   PetscErrorCode (*destroy)(PetscPartitioner);
 53:   PetscErrorCode (*partition)(PetscPartitioner, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscSection, PetscSection, PetscSection, IS *);
 54: };

 56: struct _p_PetscPartitioner {
 57:   PETSCHEADER(struct _PetscPartitionerOps);
 58:   void        *data;            /* Implementation object */
 59:   PetscInt    height;           /* Height of points to partition into non-overlapping subsets */
 60:   PetscInt    edgeCut;          /* The number of edge cut by the partition */
 61:   PetscReal   balance;          /* The maximum partition size divided by the minimum size */
 62:   PetscViewer viewer;
 63:   PetscViewer viewerGraph;
 64:   PetscBool   viewGraph;
 65:   PetscBool   noGraph;          /* if true, the partitioner does not need the connectivity graph, only the number of local vertices */
 66:   PetscBool   usevwgt;          /* if true, the partitioner looks at the local section vertSection to weight the vertices of the graph */
 67: };

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

 73: typedef struct {
 74:   MPI_Comm  pcomm;
 75:   PetscInt  ptype;
 76:   PetscReal imbalanceRatio;
 77:   PetscInt  debugFlag;
 78:   PetscInt  randomSeed;
 79: } PetscPartitioner_ParMetis;

 81: typedef struct {
 82:   MPI_Comm  pcomm;
 83:   PetscInt  strategy;
 84:   PetscReal imbalance;
 85: } PetscPartitioner_PTScotch;

 87: static const char *const
 88: PTScotchStrategyList[] = {
 89:   "DEFAULT",
 90:   "QUALITY",
 91:   "SPEED",
 92:   "BALANCE",
 93:   "SAFETY",
 94:   "SCALABILITY",
 95:   "RECURSIVE",
 96:   "REMAP"
 97: };

 99: typedef struct {
100:   PetscSection section;   /* Sizes for each partition */
101:   IS           partition; /* Points in each partition */
102:   PetscBool    random;    /* Flag for a random partition */
103: } PetscPartitioner_Shell;

105: typedef struct {
106:   PetscInt dummy;
107: } PetscPartitioner_Simple;

109: typedef struct {
110:   PetscInt dummy;
111: } PetscPartitioner_Gather;

113: typedef struct _DMPlexCellRefinerOps *DMPlexCellRefinerOps;
114: struct _DMPlexCellRefinerOps {
115:   PetscErrorCode (*refine)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
116:   PetscErrorCode (*mapsubcells)(DMPlexCellRefiner, DMPolytopeType, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
117:   PetscErrorCode (*getaffinetransforms)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
118:   PetscErrorCode (*getaffinefacetransforms)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]);
119:   PetscErrorCode (*getcellvertices)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[]);
120:   PetscErrorCode (*getsubcellvertices)(DMPlexCellRefiner, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt *, PetscInt *[]);
121:   PetscErrorCode (*mapcoords)(DMPlexCellRefiner, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
122:   PetscErrorCode (*setup)(DMPlexCellRefiner);
123:   PetscErrorCode (*destroy)(DMPlexCellRefiner);
124: };

126: struct _p_DMPlexCellRefiner {
127:   PETSCHEADER(struct _DMPlexCellRefinerOps);
128:   DM                    dm;          /* The original DM */
129:   PetscBool             setupcalled;
130:   DMPlexCellRefinerType type;
131:   PetscInt              *ctOrder;    /* [i] = ct: An array with cell types in depth order */
132:   PetscInt              *ctOrderInv; /* [ct] = i: An array with the ordinal numbers for each cell type */
133:   PetscInt              *ctStart;    /* The number for the first cell of each polytope type in the original mesh, indexed by cell type */
134:   PetscInt              *ctStartNew; /* The number for the first cell of each polytope type in the new mesh, indexed by cell type */
135:   PetscInt              *offset;     /* [ct][ctNew]: The offset in the new point numbering of a point of type ctNew produced from an old point of type ct */
136:   PetscFE               *coordFE;    /* Finite element for each cell type, used for localized coordinate interpolation */
137:   PetscFEGeom           **refGeom;   /* Geometry of the reference cell for each cell type */
138:   void                  *data;       /* refiner private data */
139: };

141: /* Utility struct to store the contents of a Fluent file in memory */
142: typedef struct {
143:   int          index;    /* Type of section */
144:   unsigned int zoneID;
145:   unsigned int first;
146:   unsigned int last;
147:   int          type;
148:   int          nd;       /* Either ND or element-type */
149:   void        *data;
150: } FluentSection;

152: struct _PetscGridHash {
153:   PetscInt     dim;
154:   PetscReal    lower[3];    /* The lower-left corner */
155:   PetscReal    upper[3];    /* The upper-right corner */
156:   PetscReal    extent[3];   /* The box size */
157:   PetscReal    h[3];        /* The subbox size */
158:   PetscInt     n[3];        /* The number of subboxes */
159:   PetscSection cellSection; /* Offsets for cells in each subbox*/
160:   IS           cells;       /* List of cells in each subbox */
161:   DMLabel      cellsSparse; /* Sparse storage for cell map */
162: };

164: /* Point Numbering in Plex:

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

168:    First Stratum:  Cells [height 0]
169:    Second Stratum: Vertices [depth 0]
170:    Third Stratum:  Faces [height 1]
171:    Fourth Stratum: Edges [depth 1]

173:    We do this so that the numbering of a cell-vertex mesh does not change after interpolation. Within a given stratum,
174:    we allow additional segregation of by cell type.
175: */
176: typedef struct {
177:   PetscInt             refct;

179:   PetscSection         coneSection;       /* Layout of cones (inedges for DAG) */
180:   PetscInt             maxConeSize;       /* Cached for fast lookup */
181:   PetscInt            *cones;             /* Cone for each point */
182:   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 */
183:   PetscSection         supportSection;    /* Layout of cones (inedges for DAG) */
184:   PetscInt             maxSupportSize;    /* Cached for fast lookup */
185:   PetscInt            *supports;          /* Cone for each point */
186:   PetscBool            refinementUniform; /* Flag for uniform cell refinement */
187:   PetscReal            refinementLimit;   /* Maximum volume for refined cell */
188:   PetscErrorCode     (*refinementFunc)(const PetscReal [], PetscReal *); /* Function giving the maximum volume for refined cell */
189:   PetscInt             overlap;           /* Overlap of the partitions as passed to DMPlexDistribute() or DMPlexDistributeOverlap() */
190:   DMPlexInterpolatedFlag interpolated;
191:   DMPlexInterpolatedFlag interpolatedCollective;

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

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

199:   /* Generation */
200:   char                *tetgenOpts;
201:   char                *triangleOpts;
202:   PetscPartitioner     partitioner;
203:   PetscBool            partitionBalance;  /* Evenly divide partition overlap when distributing */
204:   PetscBool            remeshBd;

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

211:   /* Labels and numbering */
212:   PetscObjectState     depthState;        /* State of depth label, so that we can determine if a user changes it */
213:   PetscObjectState     celltypeState;     /* State of celltype label, so that we can determine if a user changes it */
214:   IS                   globalVertexNumbers;
215:   IS                   globalCellNumbers;

217:   /* Constraints */
218:   PetscSection         anchorSection;      /* maps constrained points to anchor points */
219:   IS                   anchorIS;           /* anchors indexed by the above section */
220:   PetscErrorCode     (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
221:   PetscErrorCode     (*computeanchormatrix)(DM,PetscSection,PetscSection,Mat);

223:   /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
224:   PetscSection         parentSection;     /* dof == 1 if point has parent */
225:   PetscInt            *parents;           /* point to parent */
226:   PetscInt            *childIDs;          /* point to child ID */
227:   PetscSection         childSection;      /* inverse of parent section */
228:   PetscInt            *children;          /* point to children */
229:   DM                   referenceTree;     /* reference tree to which child ID's refer */
230:   PetscErrorCode      (*getchildsymmetry)(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*);

232:   /* MATIS support */
233:   PetscSection         subdomainSection;

235:   /* Adjacency */
236:   PetscBool            useAnchors;        /* Replace constrained points with their anchors in adjacency lists */
237:   PetscErrorCode      (*useradjacency)(DM,PetscInt,PetscInt*,PetscInt[],void*); /* User callback for adjacency */
238:   void                *useradjacencyctx;  /* User context for callback */

240:   /* Projection */
241:   PetscInt             maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */
242:   PetscInt             activePoint;         /* current active point in iteration */

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

248:   /* Geometry */
249:   PetscReal            minradius;         /* Minimum distance from cell centroid to face */
250:   PetscBool            useHashLocation;   /* Use grid hashing for point location */
251:   PetscGridHash        lbox;              /* Local box for searching */

253:   /* Neighbors */
254:   PetscMPIInt*         neighbors;

256:   /* Debugging */
257:   PetscBool            printSetValues;
258:   PetscInt             printFEM;
259:   PetscInt             printL2;
260:   PetscReal            printTol;
261: } DM_Plex;

263: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM,PetscViewer);
264: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec,PetscViewer);
265: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec,PetscViewer);
266: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec,PetscViewer);
267: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec,PetscViewer);
268: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec,PetscViewer);
269: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec,PetscViewer);
270: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
271: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM,PetscViewer);
272: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject,PetscViewer);
273: #if defined(PETSC_HAVE_HDF5)
274: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
275: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
276: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
277: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
278: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
279: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
280: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
281: #endif

283: PETSC_INTERN PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM, PetscSection, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
284: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM,PetscInt,const PetscInt[],IS*);
285: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *, DM);
286: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
287: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM []);
288: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
289: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM []);
290: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, DMLabel, DM *);
291: PETSC_INTERN PetscErrorCode DMAdaptMetric_Plex(DM, Vec, DMLabel, DM *);
292: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
293: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
294: 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);
295: 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);
296: 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);
297: 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);
298: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
299: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *);
300: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
301: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);

303: PETSC_INTERN PetscErrorCode DMPlexBuildFromCellList_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt, const int[], PetscBool);
304: PETSC_INTERN PetscErrorCode DMPlexBuildFromCellList_Parallel_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt, const int[], PetscBool, PetscSF *);
305: PETSC_INTERN PetscErrorCode DMPlexBuildCoordinates_Internal(DM, PetscInt, PetscInt, PetscInt, const double[]);
306: PETSC_INTERN PetscErrorCode DMPlexBuildCoordinates_Parallel_Internal(DM, PetscInt, PetscInt, PetscInt, PetscSF, const PetscReal[]);
307: PETSC_INTERN PetscErrorCode DMPlexLoadLabels_HDF5_Internal(DM, PetscViewer);
308: PETSC_INTERN PetscErrorCode DMPlexView_HDF5_Internal(DM, PetscViewer);
309: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Internal(DM, PetscViewer);
310: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM, PetscViewer);
311: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Internal(Vec, PetscViewer);
312: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec, PetscViewer);
313: PETSC_INTERN PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec, PetscViewer);
314: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec, PetscViewer);
315: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec, PetscViewer);
316: /* TODO Make these INTERN */
317: PETSC_EXTERN PetscErrorCode DMPlexView_ExodusII_Internal(DM, int, PetscInt);
318: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec, int, int);
319: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec, int, int);
320: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec, int, int);
321: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec, int, int);
322: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM,PetscInt,PetscInt,PetscInt*);
323: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM,PetscInt,PetscBool,PetscBool,PetscBool,PetscInt*,PetscInt*[]);
324: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM,DMPolytopeType,const PetscInt[],PetscInt*,const DMPolytopeType*[],const PetscInt*[],const PetscInt*[]);
325: PETSC_INTERN PetscErrorCode DMPlexRestoreRawFaces_Internal(DM,DMPolytopeType,const PetscInt[],PetscInt*,const DMPolytopeType*[],const PetscInt*[],const PetscInt*[]);
326: PETSC_INTERN PetscErrorCode CellRefinerInCellTest_Internal(DMPolytopeType, const PetscReal[], PetscBool *);
327: PETSC_INTERN PetscErrorCode DMPlexComputeCellType_Internal(DM, PetscInt, PetscInt, DMPolytopeType *);
328: PETSC_INTERN PetscErrorCode DMPlexCreateCellTypeOrder_Internal(DMPolytopeType, PetscInt *[], PetscInt *[]);
329: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscScalar[], InsertMode);
330: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
331: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_SetTree(DM, PetscSection, PetscInt[], PetscInt[]);
332: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM,DM,const char *,DM*);
333: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM,DM,PetscSF,PetscInt *,Mat);
334: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM,DM,PetscSF,PetscInt *,Mat);
335: PETSC_EXTERN PetscErrorCode DMPlexAnchorsModifyMat(DM,PetscSection,PetscInt,PetscInt,const PetscInt[],const PetscInt ***,const PetscScalar[],PetscInt*,PetscInt*,PetscInt*[],PetscScalar*[],PetscInt[],PetscBool);
336: PETSC_EXTERN PetscErrorCode indicesPoint_private(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,PetscInt,PetscInt []);
337: PETSC_EXTERN PetscErrorCode indicesPointFields_private(PetscSection,PetscInt,PetscInt,PetscInt [],PetscBool,PetscInt,PetscInt []);
338: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM,PetscInt,const PetscScalar [],PetscInt,PetscInt *);
339: /* these two are PETSC_EXTERN just because of src/dm/impls/plex/tests/ex18.c */
340: PETSC_EXTERN PetscErrorCode DMPlexOrientCell_Internal(DM,PetscInt,PetscInt,PetscBool);
341: PETSC_EXTERN PetscErrorCode DMPlexOrientInterface_Internal(DM);

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

346: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
347: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
348: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, DMLabel, DM *);
349: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, DMLabel, DM *);
350: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat*);

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

354: /* invert dihedral symmetry: return a^-1,
355:  * using the representation described in
356:  * DMPlexGetConeOrientation() */
357: PETSC_STATIC_INLINE PetscInt DihedralInvert(PetscInt N, PetscInt a)
358: {
359:   return (a <= 0) ? a : (N - a);
360: }

362: /* invert dihedral symmetry: return b * a,
363:  * using the representation described in
364:  * DMPlexGetConeOrientation() */
365: PETSC_STATIC_INLINE PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
366: {
367:   if (!N) return 0;
368:   return  (a >= 0) ?
369:          ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) :
370:          ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
371: }

373: /* swap dihedral symmetries: return b * a^-1,
374:  * using the representation described in
375:  * DMPlexGetConeOrientation() */
376: PETSC_STATIC_INLINE PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
377: {
378:   return DihedralCompose(N,DihedralInvert(N,a),b);
379: }

381: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, IS , PetscReal, Vec, Vec, PetscReal, Vec, void *);
382: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM, IS , PetscReal, Vec, Vec, PetscReal, Vec, void *);
383: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
384: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
385: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);

387: /* Matvec with A in row-major storage, x and y can be aliased */
388: PETSC_STATIC_INLINE void DMPlex_Mult2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
389: {
390:   PetscScalar z[2];
391:   z[0] = x[0]; z[1] = x[ldx];
392:   y[0]   = A[0]*z[0] + A[1]*z[1];
393:   y[ldx] = A[2]*z[0] + A[3]*z[1];
394:   (void)PetscLogFlops(6.0);
395: }
396: PETSC_STATIC_INLINE void DMPlex_Mult3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
397: {
398:   PetscScalar z[3];
399:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
400:   y[0]     = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
401:   y[ldx]   = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
402:   y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
403:   (void)PetscLogFlops(15.0);
404: }
405: PETSC_STATIC_INLINE void DMPlex_MultTranspose2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
406: {
407:   PetscScalar z[2];
408:   z[0] = x[0]; z[1] = x[ldx];
409:   y[0]   = A[0]*z[0] + A[2]*z[1];
410:   y[ldx] = A[1]*z[0] + A[3]*z[1];
411:   (void)PetscLogFlops(6.0);
412: }
413: PETSC_STATIC_INLINE void DMPlex_MultTranspose3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
414: {
415:   PetscScalar z[3];
416:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
417:   y[0]     = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
418:   y[ldx]   = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
419:   y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
420:   (void)PetscLogFlops(15.0);
421: }
422: PETSC_STATIC_INLINE void DMPlex_Mult2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
423: {
424:   PetscScalar z[2];
425:   z[0] = x[0]; z[1] = x[ldx];
426:   y[0]   = A[0]*z[0] + A[1]*z[1];
427:   y[ldx] = A[2]*z[0] + A[3]*z[1];
428:   (void)PetscLogFlops(6.0);
429: }
430: PETSC_STATIC_INLINE void DMPlex_Mult3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
431: {
432:   PetscScalar z[3];
433:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
434:   y[0]     = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
435:   y[ldx]   = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
436:   y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
437:   (void)PetscLogFlops(15.0);
438: }
439: PETSC_STATIC_INLINE void DMPlex_MultTransposeReal_Internal(const PetscReal A[], PetscInt m, PetscInt n, PetscInt ldx, const PetscScalar x[], PetscScalar y[])
440: {
441:   PetscScalar z[3];
442:   PetscInt    i, j;
443:   for (i = 0; i < m; ++i) z[i] = x[i*ldx];
444:   for (j = 0; j < n; ++j) {
445:     const PetscInt l = j*ldx;
446:     y[l] = 0;
447:     for (i = 0; i < m; ++i) {
448:       y[l] += A[j*n+i]*z[i];
449:     }
450:   }
451:   (void)PetscLogFlops(2*m*n);
452: }
453: PETSC_STATIC_INLINE void DMPlex_MultTranspose2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
454: {
455:   PetscScalar z[2];
456:   z[0] = x[0]; z[1] = x[ldx];
457:   y[0]   = A[0]*z[0] + A[2]*z[1];
458:   y[ldx] = A[1]*z[0] + A[3]*z[1];
459:   (void)PetscLogFlops(6.0);
460: }
461: PETSC_STATIC_INLINE void DMPlex_MultTranspose3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
462: {
463:   PetscScalar z[3];
464:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
465:   y[0]     = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
466:   y[ldx]   = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
467:   y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
468:   (void)PetscLogFlops(15.0);
469: }

471: PETSC_STATIC_INLINE void DMPlex_MatMult2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
472: {
473:   PetscInt j;
474:   for (j = 0; j < n; ++j) {
475:     PetscScalar z[2];
476:     z[0] = B[0+j]; z[1] = B[1*ldb+j];
477:     DMPlex_Mult2D_Internal(A, 1, z, z);
478:     C[0+j] = z[0]; C[1*ldb+j] = z[1];
479:   }
480:   (void)PetscLogFlops(8.0*n);
481: }
482: PETSC_STATIC_INLINE void DMPlex_MatMult3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
483: {
484:   PetscInt j;
485:   for (j = 0; j < n; ++j) {
486:     PetscScalar z[3];
487:     z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
488:     DMPlex_Mult3D_Internal(A, 1, z, z);
489:     C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
490:   }
491:   (void)PetscLogFlops(8.0*n);
492: }
493: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
494: {
495:   PetscInt j;
496:   for (j = 0; j < n; ++j) {
497:     PetscScalar z[2];
498:     z[0] = B[0+j]; z[1] = B[1*ldb+j];
499:     DMPlex_MultTranspose2D_Internal(A, 1, z, z);
500:     C[0+j] = z[0]; C[1*ldb+j] = z[1];
501:   }
502:   (void)PetscLogFlops(8.0*n);
503: }
504: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
505: {
506:   PetscInt j;
507:   for (j = 0; j < n; ++j) {
508:     PetscScalar z[3];
509:     z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
510:     DMPlex_MultTranspose3D_Internal(A, 1, z, z);
511:     C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
512:   }
513:   (void)PetscLogFlops(8.0*n);
514: }

516: PETSC_STATIC_INLINE void DMPlex_MatMultLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
517: {
518:   PetscInt j;
519:   for (j = 0; j < m; ++j) {
520:     DMPlex_MultTranspose2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
521:   }
522:   (void)PetscLogFlops(8.0*m);
523: }
524: PETSC_STATIC_INLINE void DMPlex_MatMultLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
525: {
526:   PetscInt j;
527:   for (j = 0; j < m; ++j) {
528:     DMPlex_MultTranspose3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
529:   }
530:   (void)PetscLogFlops(8.0*m);
531: }
532: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
533: {
534:   PetscInt j;
535:   for (j = 0; j < m; ++j) {
536:     DMPlex_Mult2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
537:   }
538:   (void)PetscLogFlops(8.0*m);
539: }
540: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
541: {
542:   PetscInt j;
543:   for (j = 0; j < m; ++j) {
544:     DMPlex_Mult3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
545:   }
546:   (void)PetscLogFlops(8.0*m);
547: }

549: PETSC_STATIC_INLINE void DMPlex_Transpose2D_Internal(PetscScalar A[])
550: {
551:   PetscScalar tmp;
552:   tmp = A[1]; A[1] = A[2]; A[2] = tmp;
553: }
554: PETSC_STATIC_INLINE void DMPlex_Transpose3D_Internal(PetscScalar A[])
555: {
556:   PetscScalar tmp;
557:   tmp = A[1]; A[1] = A[3]; A[3] = tmp;
558:   tmp = A[2]; A[2] = A[6]; A[6] = tmp;
559:   tmp = A[5]; A[5] = A[7]; A[7] = tmp;
560: }

562: PETSC_STATIC_INLINE void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
563: {
564:   const PetscReal invDet = 1.0/detJ;

566:   invJ[0] =  invDet*J[3];
567:   invJ[1] = -invDet*J[1];
568:   invJ[2] = -invDet*J[2];
569:   invJ[3] =  invDet*J[0];
570:   (void)PetscLogFlops(5.0);
571: }

573: PETSC_STATIC_INLINE void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
574: {
575:   const PetscReal invDet = 1.0/detJ;

577:   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
578:   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
579:   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
580:   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
581:   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
582:   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
583:   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
584:   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
585:   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
586:   (void)PetscLogFlops(37.0);
587: }

589: PETSC_STATIC_INLINE void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
590: {
591:   *detJ = J[0]*J[3] - J[1]*J[2];
592:   (void)PetscLogFlops(3.0);
593: }

595: PETSC_STATIC_INLINE void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
596: {
597:   *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
598:            J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
599:            J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
600:   (void)PetscLogFlops(12.0);
601: }

603: PETSC_STATIC_INLINE void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
604: {
605:   *detJ = PetscRealPart(J[0])*PetscRealPart(J[3]) - PetscRealPart(J[1])*PetscRealPart(J[2]);
606:   (void)PetscLogFlops(3.0);
607: }

609: PETSC_STATIC_INLINE void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
610: {
611:   *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])) +
612:            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])) +
613:            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])));
614:   (void)PetscLogFlops(12.0);
615: }

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

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

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

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

625: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Translate_Private(PetscInt ornt, PetscInt *start, PetscBool *reverse)
626: {
628:   *reverse = (ornt < 0) ? PETSC_TRUE : PETSC_FALSE;
629:   *start = *reverse ? -(ornt+1) : ornt;
630:   return(0);
631: }

633: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Combine_Private(PetscInt coneSize, PetscInt origStart, PetscBool origReverse, PetscInt rotateStart, PetscBool rotateReverse, PetscInt *newStart, PetscBool *newReverse)
634: {
636:   *newReverse = (origReverse == rotateReverse) ? PETSC_FALSE : PETSC_TRUE;
637:   *newStart = rotateReverse ? (coneSize + rotateStart - origStart) : (coneSize + origStart - rotateStart);
638:   *newStart %= coneSize;
639:   return(0);
640: }

642: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_TranslateBack_Private(PetscInt coneSize, PetscInt start, PetscBool reverse, PetscInt *ornt)
643: {
645:   if (coneSize < 3) {
646:     /* edges just get flipped if start == 1 regardless direction */
647:     *ornt = start ? -2 : 0;
648:   } else {
649:     *ornt = reverse ? -(start+1) : start;
650:   }
651:   return(0);
652: }

654: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Permute_Private(PetscInt n, const PetscInt arr[], PetscInt start, PetscBool reverse, PetscInt newarr[])
655: {
656:   PetscInt i;

659:   if (reverse) {for (i=0; i<n; i++) newarr[i] = arr[(n+start-i)%n];}
660:   else         {for (i=0; i<n; i++) newarr[i] = arr[(start+i)%n];}
661:   return(0);
662: }

664: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM,PetscInt,PetscInt,PetscDualSpace *);
665: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection,PetscBool,PetscInt,PetscInt,PetscInt *,PetscBool,const PetscInt[],const PetscInt[],PetscInt[]);
666: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection,PetscBool,PetscInt,PetscInt,PetscInt[],PetscBool,const PetscInt***,PetscInt,const PetscInt[],PetscInt[]);
667: PETSC_INTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
668: PETSC_INTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);

670: PETSC_EXTERN PetscErrorCode DMSNESGetFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
671: PETSC_EXTERN PetscErrorCode DMSNESRestoreFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
672: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM, PetscSection, IS, PetscReal, Vec, Vec, Vec, void *);
673: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM, PetscSection, PetscSection, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
674: PETSC_INTERN PetscErrorCode DMCreateSubDomainDM_Plex(DM,DMLabel,PetscInt,IS*,DM*);
675: PETSC_INTERN PetscErrorCode DMPlexBasisTransformPoint_Internal(DM, DM, Vec, PetscInt, PetscBool[], PetscBool, PetscScalar *);
676: PETSC_EXTERN PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM, DM, Vec, PetscInt, PetscBool, PetscInt, PetscScalar *);
677: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscReal *, PetscReal *, void *);
678: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApply_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscScalar *, PetscScalar *, void *);
679: PETSC_INTERN PetscErrorCode DMCreateNeumannOverlap_Plex(DM, IS*, Mat*, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void **);

681: #endif /* _PLEXIMPL_H */