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