1 /* 2 DMPlex, for parallel unstructured distributed mesh problems. 3 */ 4 #if !defined(__PETSCDMPLEX_H) 5 #define __PETSCDMPLEX_H 6 7 #include <petscdm.h> 8 #include <petscdt.h> 9 #include <petscfe.h> 10 #include <petscfv.h> 11 #include <petscsftypes.h> 12 13 PETSC_EXTERN PetscClassId PETSCPARTITIONER_CLASSID; 14 15 /*J 16 PetscPartitionerType - String with the name of a PETSc graph partitioner 17 18 Level: beginner 19 20 .seealso: PetscPartitionerSetType(), PetscPartitioner 21 J*/ 22 typedef const char *PetscPartitionerType; 23 #define PETSCPARTITIONERCHACO "chaco" 24 #define PETSCPARTITIONERPARMETIS "parmetis" 25 #define PETSCPARTITIONERSHELL "shell" 26 #define PETSCPARTITIONERSIMPLE "simple" 27 #define PETSCPARTITIONERGATHER "gather" 28 29 PETSC_EXTERN PetscFunctionList PetscPartitionerList; 30 PETSC_EXTERN PetscErrorCode PetscPartitionerCreate(MPI_Comm, PetscPartitioner *); 31 PETSC_EXTERN PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *); 32 PETSC_EXTERN PetscErrorCode PetscPartitionerSetType(PetscPartitioner, PetscPartitionerType); 33 PETSC_EXTERN PetscErrorCode PetscPartitionerGetType(PetscPartitioner, PetscPartitionerType *); 34 PETSC_EXTERN PetscErrorCode PetscPartitionerSetUp(PetscPartitioner); 35 PETSC_EXTERN PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner); 36 PETSC_STATIC_INLINE PetscErrorCode PetscPartitionerViewFromOptions(PetscPartitioner A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);} 37 PETSC_EXTERN PetscErrorCode PetscPartitionerView(PetscPartitioner, PetscViewer); 38 PETSC_EXTERN PetscErrorCode PetscPartitionerRegister(const char [], PetscErrorCode (*)(PetscPartitioner)); 39 PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterDestroy(void); 40 41 PETSC_EXTERN PetscErrorCode PetscPartitionerPartition(PetscPartitioner, DM, PetscSection, IS *); 42 43 PETSC_EXTERN PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner, PetscInt, const PetscInt[], const PetscInt[]); 44 PETSC_EXTERN PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner, PetscBool); 45 PETSC_EXTERN PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner, PetscBool *); 46 47 PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*); 48 PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char [], PetscInt, DM *); 49 PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*); 50 PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], PetscSF *, DM *); 51 PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []); 52 PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm, PetscInt, PetscBool, DM*); 53 PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *); 54 PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt); 55 PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *); 56 PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt); 57 PETSC_EXTERN PetscErrorCode DMPlexAddConeSize(DM, PetscInt, PetscInt); 58 PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]); 59 PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]); 60 PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt); 61 PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt); 62 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]); 63 PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]); 64 PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *); 65 PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt); 66 PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]); 67 PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]); 68 PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt); 69 PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *); 70 PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *); 71 PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]); 72 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]); 73 PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *); 74 PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM); 75 PETSC_EXTERN PetscErrorCode DMPlexStratify(DM); 76 PETSC_EXTERN PetscErrorCode DMPlexEqual(DM,DM,PetscBool*); 77 PETSC_EXTERN PetscErrorCode DMPlexReverseCell(DM, PetscInt); 78 PETSC_EXTERN PetscErrorCode DMPlexOrient(DM); 79 PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *); 80 PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *); 81 PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool); 82 PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*); 83 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,void*); 84 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*); 85 PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*); 86 PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*); 87 PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*,void*); 88 PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*); 89 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*); 90 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*); 91 PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*); 92 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*); 93 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*, void*); 94 95 PETSC_EXTERN PetscErrorCode DMPlexFilter(DM, DMLabel, PetscInt, DM *); 96 PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *); 97 PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *); 98 PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *); 99 100 PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *); 101 PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *); 102 PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *); 103 PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 104 105 /* Topological Operations */ 106 PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 107 PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 108 PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 109 PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 110 PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 111 PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 112 PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 113 PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 114 115 /* Mesh Generation */ 116 PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *); 117 PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM); 118 PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscBool, PetscReal, DM *); 119 PETSC_EXTERN PetscErrorCode DMPlexCreateSquareBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 120 PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 121 PETSC_EXTERN PetscErrorCode DMPlexCreateSquareMesh(DM, const PetscReal [], const PetscReal [], const PetscInt [], DMBoundaryType, DMBoundaryType); 122 PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscInt, PetscBool, DM *); 123 PETSC_EXTERN PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm, PetscInt, const PetscInt[], DMBoundaryType, DMBoundaryType, DMBoundaryType, DM *); 124 PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *); 125 PETSC_EXTERN PetscErrorCode DMPlexInvertCell(PetscInt, PetscInt, int []); 126 PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM); 127 PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscBool, PetscInt); 128 PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscBool, PetscInt); 129 130 PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *); 131 PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *); 132 133 /* Mesh Partitioning and Distribution */ 134 PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **); 135 PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *); 136 PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner); 137 PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, const char[], PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *); 138 PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **, IS *); 139 PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *); 140 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel); 141 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel); 142 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel); 143 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelPropagate(DM, DMLabel); 144 PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *); 145 PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF*, DM*); 146 PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *); 147 PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM,PetscSF,PetscSection,Vec,PetscSection,Vec); 148 PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *); 149 PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM,PetscSF,PetscSection,MPI_Datatype,void*,PetscSection,void**); 150 PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM); 151 PETSC_EXTERN PetscErrorCode DMPlexGetRedundantDM(DM, DM*); 152 PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseCone(DM,PetscBool); 153 PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseCone(DM,PetscBool*); 154 PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseClosure(DM,PetscBool); 155 PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseClosure(DM,PetscBool*); 156 PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM,PetscBool); 157 PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM,PetscBool*); 158 PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]); 159 160 PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, DMLabel, IS *); 161 PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *); 162 163 PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *); 164 PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *); 165 PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *); 166 PETSC_EXTERN PetscErrorCode DMPlexCreateOverlap(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *); 167 PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *); 168 PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *); 169 170 /* Submesh Support */ 171 PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, DM*); 172 PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel *, DM *); 173 PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*); 174 PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel); 175 PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *); 176 177 PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel); 178 PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscBool, DM); 179 PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel); 180 PETSC_EXTERN PetscErrorCode DMPlexLabelClearCells(DM, DMLabel); 181 182 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *); 183 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal); 184 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *); 185 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool); 186 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementFunction(DM, PetscErrorCode (**)(const PetscReal [], PetscReal *)); 187 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementFunction(DM, PetscErrorCode (*)(const PetscReal [], PetscReal *)); 188 PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *); 189 PETSC_EXTERN PetscErrorCode DMPlexGetRegularRefinement(DM, PetscBool *); 190 PETSC_EXTERN PetscErrorCode DMPlexSetRegularRefinement(DM, PetscBool); 191 192 /* Support for cell-vertex meshes */ 193 PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *); 194 PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *); 195 196 /* Geometry Support */ 197 PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *); 198 PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal); 199 PETSC_EXTERN PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar[], PetscReal[]); 200 PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar[], PetscReal[]); 201 PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt, PetscScalar[], PetscReal[]); 202 203 /* Point Location */ 204 typedef struct _PetscGridHash *PetscGridHash; 205 PETSC_EXTERN PetscErrorCode PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *); 206 PETSC_EXTERN PetscErrorCode PetscGridHashEnlarge(PetscGridHash, const PetscScalar []); 207 PETSC_EXTERN PetscErrorCode PetscGridHashSetGrid(PetscGridHash, const PetscInt [], const PetscReal []); 208 PETSC_EXTERN PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscScalar [], PetscInt [], PetscInt []); 209 PETSC_EXTERN PetscErrorCode PetscGridHashDestroy(PetscGridHash *); 210 211 /* FVM Support */ 212 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []); 213 PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *); 214 PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *); 215 PETSC_EXTERN PetscErrorCode DMPlexGetDataFVM(DM, PetscFV, Vec *, Vec *, DM *); 216 217 /* FEM Support */ 218 PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec); 219 PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, DMLabel); 220 221 PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, PetscInt, PetscInt, const PetscInt [], const PetscInt [], PetscInt, const PetscInt [], const IS [], const IS [], IS, PetscSection *); 222 PETSC_EXTERN PetscErrorCode DMPlexGetSubdomainSection(DM, PetscSection*); 223 224 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 225 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 226 PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFEM(DM, Vec *); 227 PETSC_EXTERN PetscErrorCode DMPlexCoordinatesToReference(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]); 228 PETSC_EXTERN PetscErrorCode DMPlexReferenceToCoordinates(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]); 229 230 PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 231 PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 232 PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode); 233 PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 234 PETSC_EXTERN PetscErrorCode DMPlexGetClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscInt *); 235 PETSC_EXTERN PetscErrorCode DMPlexRestoreClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscInt *, PetscInt **,PetscInt *); 236 PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 237 PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]); 238 PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection); 239 PETSC_EXTERN PetscErrorCode DMPlexCreateSpectralClosurePermutation(DM, PetscSection); 240 241 PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], PetscBool, DM *); 242 PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *); 243 PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char [], PetscBool, DM *); 244 PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *); 245 PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *); 246 PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *); 247 PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char [], PetscBool, DM *); 248 PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *); 249 PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char [], PetscBool, DM *); 250 PETSC_EXTERN PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm, const char [], PetscBool, DM *); 251 252 PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *); 253 PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DM *); 254 255 PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 256 PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt); 257 PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *); 258 PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt); 259 PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer); 260 261 PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *); 262 PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal); 263 264 typedef struct { 265 DM dm; 266 Vec u; /* The base vector for the Jacbobian action J(u) x */ 267 Mat J; /* Preconditioner for testing */ 268 void *user; 269 } JacActionCtx; 270 271 PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM, Vec); 272 PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt); 273 PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt*); 274 PETSC_EXTERN PetscErrorCode DMPlexProjectFieldLocal(DM, Vec, 275 void (**)(PetscInt, PetscInt, PetscInt, 276 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 277 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 278 PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec); 279 PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]); 280 PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffVec(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, Vec); 281 PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscReal *, void *); 282 PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorNested(DM, DM, Mat, void *); 283 PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorGeneral(DM, DM, Mat, void *); 284 PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *); 285 286 PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, MatNullSpace *); 287 288 PETSC_EXTERN PetscErrorCode DMPlexSetSNESLocalFEM(DM,void *,void *,void *); 289 PETSC_EXTERN PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM, Vec, void *); 290 PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *); 291 PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat, void *); 292 PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianActionFEM(DM, Vec, Vec, Vec, void *); 293 294 PETSC_EXTERN PetscErrorCode DMPlexTSComputeBoundary(DM, PetscReal, Vec, Vec, void *); 295 PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *); 296 PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *); 297 PETSC_EXTERN PetscErrorCode DMPlexTSComputeIJacobianFEM(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 298 299 PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *); 300 PETSC_EXTERN PetscErrorCode DMPlexReconstructGradientsFVM(DM,Vec,Vec); 301 302 /* anchors */ 303 PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection*, IS*); 304 PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS); 305 /* tree */ 306 PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM); 307 PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM*); 308 PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *); 309 PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM*); 310 PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt [], PetscInt []); 311 PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]); 312 PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *); 313 PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]); 314 PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *); 315 PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorReferenceTree(DM, Mat *); 316 PETSC_EXTERN PetscErrorCode DMPlexTransferVecTree(DM,Vec,DM,Vec,PetscSF,PetscSF,PetscInt *,PetscInt *,PetscBool,PetscReal); 317 318 /* natural order */ 319 PETSC_EXTERN PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM, PetscSection, PetscSF, PetscSF *); 320 PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalBegin(DM, Vec, Vec); 321 PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalEnd(DM, Vec, Vec); 322 PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalBegin(DM, Vec, Vec); 323 PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalEnd(DM, Vec, Vec); 324 325 /* mesh adaptation */ 326 PETSC_EXTERN PetscErrorCode DMPlexAdapt(DM, Vec, const char [], DM *); 327 #endif 328