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 DMPlexGetCones(DM, PetscInt *[]); 34 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]); 35 PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *); 36 PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM); 37 PETSC_EXTERN PetscErrorCode DMPlexStratify(DM); 38 PETSC_EXTERN PetscErrorCode DMPlexOrient(DM); 39 PETSC_EXTERN PetscErrorCode DMPlexGetCoordinateSection(DM, PetscSection *); 40 PETSC_EXTERN PetscErrorCode DMPlexSetCoordinateSection(DM, PetscSection); 41 PETSC_EXTERN PetscErrorCode DMPlexSetPreallocationCenterDimension(DM,PetscInt); 42 PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscSection, PetscSection, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool); 43 PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*); 44 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*); 45 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,const void*); 46 PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*); 47 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*); 48 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*); 49 50 /*S 51 DMLabel - Object which encapsulates a subset of the mesh from this DM 52 53 Level: developer 54 55 Concepts: grids, grid refinement 56 57 .seealso: DM, DMPlexCreate(), DMPlexCreateLabel() 58 S*/ 59 typedef struct _n_DMLabel *DMLabel; 60 PETSC_EXTERN PetscErrorCode DMLabelCreate(const char [], DMLabel *); 61 PETSC_EXTERN PetscErrorCode DMLabelView(DMLabel, PetscViewer); 62 PETSC_EXTERN PetscErrorCode DMLabelDestroy(DMLabel *); 63 PETSC_EXTERN PetscErrorCode DMLabelDuplicate(DMLabel, DMLabel *); 64 PETSC_EXTERN PetscErrorCode DMLabelGetName(DMLabel, const char **); 65 PETSC_EXTERN PetscErrorCode DMLabelGetValue(DMLabel, PetscInt, PetscInt *); 66 PETSC_EXTERN PetscErrorCode DMLabelSetValue(DMLabel, PetscInt, PetscInt); 67 PETSC_EXTERN PetscErrorCode DMLabelClearValue(DMLabel, PetscInt, PetscInt); 68 PETSC_EXTERN PetscErrorCode DMLabelGetNumValues(DMLabel, PetscInt *); 69 PETSC_EXTERN PetscErrorCode DMLabelGetValueIS(DMLabel, IS *); 70 PETSC_EXTERN PetscErrorCode DMLabelGetStratumSize(DMLabel, PetscInt, PetscInt *); 71 PETSC_EXTERN PetscErrorCode DMLabelGetStratumIS(DMLabel, PetscInt, IS *); 72 PETSC_EXTERN PetscErrorCode DMLabelClearStratum(DMLabel, PetscInt); 73 PETSC_EXTERN PetscErrorCode DMLabelCreateIndex(DMLabel, PetscInt, PetscInt); 74 PETSC_EXTERN PetscErrorCode DMLabelDestroyIndex(DMLabel); 75 PETSC_EXTERN PetscErrorCode DMLabelHasPoint(DMLabel, PetscInt, PetscBool *); 76 77 PETSC_EXTERN PetscErrorCode DMPlexCreateLabel(DM, const char []); 78 PETSC_EXTERN PetscErrorCode DMPlexGetLabelValue(DM, const char[], PetscInt, PetscInt *); 79 PETSC_EXTERN PetscErrorCode DMPlexSetLabelValue(DM, const char[], PetscInt, PetscInt); 80 PETSC_EXTERN PetscErrorCode DMPlexClearLabelValue(DM, const char[], PetscInt, PetscInt); 81 PETSC_EXTERN PetscErrorCode DMPlexGetLabelSize(DM, const char[], PetscInt *); 82 PETSC_EXTERN PetscErrorCode DMPlexGetLabelIdIS(DM, const char[], IS *); 83 PETSC_EXTERN PetscErrorCode DMPlexGetStratumSize(DM, const char [], PetscInt, PetscInt *); 84 PETSC_EXTERN PetscErrorCode DMPlexGetStratumIS(DM, const char [], PetscInt, IS *); 85 PETSC_EXTERN PetscErrorCode DMPlexClearLabelStratum(DM, const char[], PetscInt); 86 PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection, PetscSF, PetscBool, DMLabel, PetscInt, PetscSection *); 87 88 PETSC_EXTERN PetscErrorCode DMPlexGetNumLabels(DM, PetscInt *); 89 PETSC_EXTERN PetscErrorCode DMPlexGetLabelName(DM, PetscInt, const char **); 90 PETSC_EXTERN PetscErrorCode DMPlexHasLabel(DM, const char [], PetscBool *); 91 PETSC_EXTERN PetscErrorCode DMPlexGetLabel(DM, const char *, DMLabel *); 92 PETSC_EXTERN PetscErrorCode DMPlexAddLabel(DM, DMLabel); 93 PETSC_EXTERN PetscErrorCode DMPlexRemoveLabel(DM, const char [], DMLabel *); 94 PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *); 95 PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *); 96 97 PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 98 PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 99 PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 100 PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 101 PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 102 PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 103 PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 104 PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 105 106 PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **); 107 PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *); 108 PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *); 109 110 PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *); 111 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *); 112 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal); 113 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *); 114 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool); 115 PETSC_EXTERN PetscErrorCode DMPlexInvertCell(PetscInt, PetscInt, int []); 116 PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *); 117 PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM); 118 PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, const char[], PetscInt, DM*); 119 PETSC_EXTERN PetscErrorCode DMPlexLoad(PetscViewer, DM); 120 PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*); 121 PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel); 122 PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *); 123 124 PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 125 PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, DM *); 126 PETSC_EXTERN PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm,PetscInt,const PetscInt[],DM *); 127 PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *); 128 PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *); 129 PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *); 130 PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 131 PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, PetscInt, PetscInt,const PetscInt [],const PetscInt [], PetscInt,const PetscInt [],const IS [], PetscSection *); 132 PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *); 133 PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, PetscBool, DM); 134 PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, DMLabel); 135 PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel); 136 137 /* Support for cell-vertex meshes */ 138 PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *); 139 PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *); 140 141 /* FVM Support */ 142 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []); 143 144 /* FEM Support */ 145 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometry(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 146 PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 147 PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 148 PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode); 149 PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 150 151 PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *); 152 PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *); 153 154 PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *); 155 PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DM *); 156 157 PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 158 PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt); 159 PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *); 160 PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt); 161 PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer); 162 163 PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *); 164 PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal); 165 166 typedef struct { 167 DM dm; 168 Vec u; /* The base vector for the Jacbobian action J(u) x */ 169 Mat J; /* Preconditioner for testing */ 170 void *user; 171 } JacActionCtx; 172 173 PETSC_EXTERN PetscErrorCode DMPlexProjectFunction(DM, PetscInt, void (**)(const PetscReal [], PetscScalar *), InsertMode, Vec); 174 PETSC_EXTERN PetscErrorCode DMPlexProjectFunctionLocal(DM, PetscInt, void (**)(const PetscReal [], PetscScalar *), InsertMode, Vec); 175 PETSC_EXTERN PetscErrorCode DMPlexComputeL2Diff(DM, PetscQuadrature[], void (**)(const PetscReal [], PetscScalar *), Vec, PetscReal *); 176 PETSC_EXTERN PetscErrorCode DMPlexComputeResidualFEM(DM, Vec, Vec, void *); 177 PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianActionFEM(DM, Mat, Vec, Vec, void *); 178 PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianFEM(DM, Vec, Mat, Mat, MatStructure*,void *); 179 PETSC_EXTERN PetscErrorCode DMPlexSetFEMIntegration(DM, 180 PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 181 const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 182 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 183 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]), 184 PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 185 const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 186 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 187 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), PetscScalar[]), 188 PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], const PetscScalar[], 189 const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 190 void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 191 void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 192 void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 193 void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]), 194 PetscErrorCode (*)(PetscInt, 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[]), 198 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 199 void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[])); 200 #endif 201