1 /* 2 DMPlex, for parallel unstructured distributed mesh problems. 3 */ 4 #if !defined(__PETSCDMPLEX_H) 5 #define __PETSCDMPLEX_H 6 7 #include <petscsftypes.h> 8 #include <petscdm.h> 9 10 PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*); 11 PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, const char[], PetscInt, DM*); 12 PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char [], PetscInt, DM *); 13 PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*); 14 PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []); 15 PETSC_EXTERN PetscErrorCode DMPlexGetDimension(DM, PetscInt *); 16 PETSC_EXTERN PetscErrorCode DMPlexSetDimension(DM, PetscInt); 17 PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *); 18 PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt); 19 PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *); 20 PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt); 21 PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]); 22 PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]); 23 PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt); 24 PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt); 25 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]); 26 PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]); 27 PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *); 28 PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt); 29 PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]); 30 PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]); 31 PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt); 32 PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *); 33 PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *); 34 PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]); 35 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]); 36 PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *); 37 PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM); 38 PETSC_EXTERN PetscErrorCode DMPlexStratify(DM); 39 PETSC_EXTERN PetscErrorCode DMPlexEqual(DM,DM,PetscBool*); 40 PETSC_EXTERN PetscErrorCode DMPlexOrient(DM); 41 PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *); 42 PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *); 43 PETSC_EXTERN PetscErrorCode DMPlexLoad(PetscViewer, DM); 44 PETSC_EXTERN PetscErrorCode DMPlexGetCoordinateSection(DM, PetscSection *); 45 PETSC_EXTERN PetscErrorCode DMPlexSetCoordinateSection(DM, PetscSection); 46 PETSC_EXTERN PetscErrorCode DMPlexSetPreallocationCenterDimension(DM,PetscInt); 47 PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscSection, PetscSection, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool); 48 PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*); 49 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*); 50 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,const void*); 51 PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*); 52 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*); 53 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*); 54 55 /*S 56 DMLabel - Object which encapsulates a subset of the mesh from this DM 57 58 Level: developer 59 60 Concepts: grids, grid refinement 61 62 .seealso: DM, DMPlexCreate(), DMPlexCreateLabel() 63 S*/ 64 typedef struct _n_DMLabel *DMLabel; 65 PETSC_EXTERN PetscErrorCode DMLabelCreate(const char [], DMLabel *); 66 PETSC_EXTERN PetscErrorCode DMLabelView(DMLabel, PetscViewer); 67 PETSC_EXTERN PetscErrorCode DMLabelDestroy(DMLabel *); 68 PETSC_EXTERN PetscErrorCode DMLabelDuplicate(DMLabel, DMLabel *); 69 PETSC_EXTERN PetscErrorCode DMLabelGetName(DMLabel, const char **); 70 PETSC_EXTERN PetscErrorCode DMLabelGetValue(DMLabel, PetscInt, PetscInt *); 71 PETSC_EXTERN PetscErrorCode DMLabelSetValue(DMLabel, PetscInt, PetscInt); 72 PETSC_EXTERN PetscErrorCode DMLabelClearValue(DMLabel, PetscInt, PetscInt); 73 PETSC_EXTERN PetscErrorCode DMLabelGetNumValues(DMLabel, PetscInt *); 74 PETSC_EXTERN PetscErrorCode DMLabelGetValueIS(DMLabel, IS *); 75 PETSC_EXTERN PetscErrorCode DMLabelGetStratumSize(DMLabel, PetscInt, PetscInt *); 76 PETSC_EXTERN PetscErrorCode DMLabelGetStratumIS(DMLabel, PetscInt, IS *); 77 PETSC_EXTERN PetscErrorCode DMLabelClearStratum(DMLabel, PetscInt); 78 PETSC_EXTERN PetscErrorCode DMLabelCreateIndex(DMLabel, PetscInt, PetscInt); 79 PETSC_EXTERN PetscErrorCode DMLabelDestroyIndex(DMLabel); 80 PETSC_EXTERN PetscErrorCode DMLabelHasPoint(DMLabel, PetscInt, PetscBool *); 81 82 PETSC_EXTERN PetscErrorCode DMPlexCreateLabel(DM, const char []); 83 PETSC_EXTERN PetscErrorCode DMPlexGetLabelValue(DM, const char[], PetscInt, PetscInt *); 84 PETSC_EXTERN PetscErrorCode DMPlexSetLabelValue(DM, const char[], PetscInt, PetscInt); 85 PETSC_EXTERN PetscErrorCode DMPlexClearLabelValue(DM, const char[], PetscInt, PetscInt); 86 PETSC_EXTERN PetscErrorCode DMPlexGetLabelSize(DM, const char[], PetscInt *); 87 PETSC_EXTERN PetscErrorCode DMPlexGetLabelIdIS(DM, const char[], IS *); 88 PETSC_EXTERN PetscErrorCode DMPlexGetStratumSize(DM, const char [], PetscInt, PetscInt *); 89 PETSC_EXTERN PetscErrorCode DMPlexGetStratumIS(DM, const char [], PetscInt, IS *); 90 PETSC_EXTERN PetscErrorCode DMPlexClearLabelStratum(DM, const char[], PetscInt); 91 PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection, PetscSF, PetscBool, DMLabel, PetscInt, PetscSection *); 92 93 PETSC_EXTERN PetscErrorCode DMPlexGetNumLabels(DM, PetscInt *); 94 PETSC_EXTERN PetscErrorCode DMPlexGetLabelName(DM, PetscInt, const char **); 95 PETSC_EXTERN PetscErrorCode DMPlexHasLabel(DM, const char [], PetscBool *); 96 PETSC_EXTERN PetscErrorCode DMPlexGetLabel(DM, const char *, DMLabel *); 97 PETSC_EXTERN PetscErrorCode DMPlexAddLabel(DM, DMLabel); 98 PETSC_EXTERN PetscErrorCode DMPlexRemoveLabel(DM, const char [], DMLabel *); 99 PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *); 100 PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *); 101 102 PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *); 103 PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *); 104 PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *); 105 PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 106 107 /* Topological Operations */ 108 PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 109 PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 110 PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 111 PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 112 PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 113 PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 114 PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 115 PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 116 117 /* Mesh Generation */ 118 PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *); 119 PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM); 120 PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 121 PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, DM *); 122 PETSC_EXTERN PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm,PetscInt,const PetscInt[],DM *); 123 PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *); 124 PETSC_EXTERN PetscErrorCode DMPlexInvertCell(PetscInt, PetscInt, int []); 125 126 /* Mesh Partitioning and Distribution */ 127 PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **); 128 PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, const char[], PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *); 129 PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *); 130 PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, const char[], PetscInt, DM*); 131 132 /* Submesh Support */ 133 PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*); 134 PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel); 135 PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *); 136 137 PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, DMLabel); 138 PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel); 139 PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, PetscBool, DM); 140 141 /* Mesh Refinement */ 142 typedef PetscInt CellRefiner; 143 144 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *); 145 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal); 146 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *); 147 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool); 148 149 /* Support for cell-vertex meshes */ 150 PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *); 151 PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *); 152 153 /* FVM Support */ 154 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []); 155 156 /* FEM Support */ 157 PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, PetscInt, PetscInt,const PetscInt [],const PetscInt [], PetscInt,const PetscInt [],const IS [], PetscSection *); 158 159 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometry(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 160 PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 161 PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 162 PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode); 163 PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 164 165 PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *); 166 PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *); 167 168 PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *); 169 PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DM *); 170 171 PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 172 PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt); 173 PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *); 174 PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt); 175 PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer); 176 177 PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *); 178 PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal); 179 180 typedef struct { 181 DM dm; 182 Vec u; /* The base vector for the Jacbobian action J(u) x */ 183 Mat J; /* Preconditioner for testing */ 184 void *user; 185 } JacActionCtx; 186 187 PETSC_EXTERN PetscErrorCode DMPlexProjectFunction(DM, PetscInt, void (**)(const PetscReal [], PetscScalar *), InsertMode, Vec); 188 PETSC_EXTERN PetscErrorCode DMPlexProjectFunctionLocal(DM, PetscInt, void (**)(const PetscReal [], PetscScalar *), InsertMode, Vec); 189 PETSC_EXTERN PetscErrorCode DMPlexComputeL2Diff(DM, PetscQuadrature[], void (**)(const PetscReal [], PetscScalar *), Vec, PetscReal *); 190 PETSC_EXTERN PetscErrorCode DMPlexComputeResidualFEM(DM, Vec, Vec, void *); 191 PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianActionFEM(DM, Mat, Vec, Vec, void *); 192 PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianFEM(DM, Vec, Mat, Mat, MatStructure*,void *); 193 PETSC_EXTERN PetscErrorCode DMPlexSetFEMIntegration(DM, 194 PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 195 const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 196 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 197 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]), 198 PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 199 const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 200 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 201 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), PetscScalar[]), 202 PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], const PetscScalar[], 203 const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 204 void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 205 void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 206 void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 207 void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]), 208 PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 209 const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 210 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 211 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 212 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 213 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[])); 214 #endif 215