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