1552f7358SJed Brown /* 2552f7358SJed Brown DMPlex, for parallel unstructured distributed mesh problems. 3552f7358SJed Brown */ 4552f7358SJed Brown #if !defined(__PETSCDMPLEX_H) 5552f7358SJed Brown #define __PETSCDMPLEX_H 6552f7358SJed Brown 7552f7358SJed Brown #include <petscsf.h> 8552f7358SJed Brown #include <petscdm.h> 9552f7358SJed Brown 10552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*); 11a89082eeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, const char[], DM*); 120fde5641SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*); 1306bbee04SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []); 14552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexClone(DM, DM*); 15552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetDimension(DM, PetscInt *); 16552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetDimension(DM, PetscInt); 17552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *); 18552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt); 19552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *); 20552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt); 21552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]); 22552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]); 239f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt); 2477c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt); 25552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]); 26552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]); 27552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *); 28552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt); 29552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]); 30552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]); 31552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt); 32552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *); 33552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]); 34552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]); 35552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *); 36552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM); 37552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexStratify(DM); 38552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCoordinateSection(DM, PetscSection *); 39552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCoordinateSection(DM, PetscSection); 40552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetPreallocationCenterDimension(DM,PetscInt); 41cb1e1211SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscSection, PetscSection, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool); 42552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*); 43552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*); 44552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,const void*); 45552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*); 46552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*); 47552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*); 48552f7358SJed Brown 49552f7358SJed Brown /*S 50552f7358SJed Brown DMLabel - Object which encapsulates a subset of the mesh from this DM 51552f7358SJed Brown 52552f7358SJed Brown Level: developer 53552f7358SJed Brown 54552f7358SJed Brown Concepts: grids, grid refinement 55552f7358SJed Brown 56552f7358SJed Brown .seealso: DM, DMPlexCreate(), DMPlexCreateLabel() 57552f7358SJed Brown S*/ 58552f7358SJed Brown typedef struct _n_DMLabel *DMLabel; 59552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelCreate(const char [], DMLabel *); 60552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelView(DMLabel, PetscViewer); 61552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelDestroy(DMLabel *); 62f807ecf7SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelDuplicate(DMLabel, DMLabel *); 6385de2ee3SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelGetName(DMLabel, const char **); 64552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetValue(DMLabel, PetscInt, PetscInt *); 65552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelSetValue(DMLabel, PetscInt, PetscInt); 66552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelClearValue(DMLabel, PetscInt, PetscInt); 67552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetNumValues(DMLabel, PetscInt *); 68552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetValueIS(DMLabel, IS *); 69552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetStratumSize(DMLabel, PetscInt, PetscInt *); 70552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetStratumIS(DMLabel, PetscInt, IS *); 71552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelClearStratum(DMLabel, PetscInt); 7291d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelCreateIndex(DMLabel, PetscInt, PetscInt); 7391d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelDestroyIndex(DMLabel); 7491d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelHasPoint(DMLabel, PetscInt, PetscBool *); 75552f7358SJed Brown 76552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateLabel(DM, const char []); 77552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelValue(DM, const char[], PetscInt, PetscInt *); 78552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetLabelValue(DM, const char[], PetscInt, PetscInt); 79552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexClearLabelValue(DM, const char[], PetscInt, PetscInt); 80552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelSize(DM, const char[], PetscInt *); 81552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelIdIS(DM, const char[], IS *); 82552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetStratumSize(DM, const char [], PetscInt, PetscInt *); 83552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetStratumIS(DM, const char [], PetscInt, IS *); 84552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexClearLabelStratum(DM, const char[], PetscInt); 85552f7358SJed Brown PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection, PetscSF, PetscBool, DMLabel, PetscInt, PetscSection *); 86552f7358SJed Brown 87552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetNumLabels(DM, PetscInt *); 88552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelName(DM, PetscInt, const char **); 89552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexHasLabel(DM, const char [], PetscBool *); 90552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabel(DM, const char *, DMLabel *); 91552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexAddLabel(DM, DMLabel); 92552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRemoveLabel(DM, const char [], DMLabel *); 93552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *); 94552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *); 95552f7358SJed Brown 96552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 97552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 98552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 99552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 100552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 101552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 102552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 103552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 104552f7358SJed Brown 10577c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **); 106552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *); 107552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *); 108552f7358SJed Brown 109552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *); 110552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *); 111552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal); 112552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *); 113552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool); 1149f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *); 11580330477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM); 116552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, const char[], PetscInt, DM*); 117552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexLoad(PetscViewer, DM); 118a89082eeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*); 119a89082eeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel); 120e6ccafaeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *); 121552f7358SJed Brown 122552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 123552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, DM *); 124552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm,PetscInt,const PetscInt[],DM *); 125552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *); 126552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *); 127552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 128552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, PetscInt, PetscInt,const PetscInt [],const PetscInt [], PetscInt,const PetscInt [],const IS [], PetscSection *); 129552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *); 1302be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel); 13109f723d9SJed Brown PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, DMLabel); 1322be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel); 133552f7358SJed Brown 134834e62ceSMatthew G. Knepley /* FVM Support */ 135*cc08537eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []); 136834e62ceSMatthew G. Knepley 137552f7358SJed Brown /* FEM Support */ 138552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometry(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 1397c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 1407c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 141552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode); 1427c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 143552f7358SJed Brown 144552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *); 1455bb3eff3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *); 146552f7358SJed Brown 147552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *); 148a89082eeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DM *); 149552f7358SJed Brown 150770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 151770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt); 152552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *); 153552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt); 154552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer); 155552f7358SJed Brown 156552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *); 157552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal); 158552f7358SJed Brown 159552f7358SJed Brown typedef struct { 160552f7358SJed Brown DM dm; 161552f7358SJed Brown Vec u; /* The base vector for the Jacbobian action J(u) x */ 162552f7358SJed Brown Mat J; /* Preconditioner for testing */ 163552f7358SJed Brown void *user; 164552f7358SJed Brown } JacActionCtx; 165552f7358SJed Brown 166552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexProjectFunction(DM, PetscInt, PetscScalar (**)(const PetscReal []), InsertMode, Vec); 167552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexProjectFunctionLocal(DM, PetscInt, PetscScalar (**)(const PetscReal []), InsertMode, Vec); 168552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeL2Diff(DM, PetscQuadrature[], PetscScalar (**)(const PetscReal []), Vec, PetscReal *); 169552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeResidualFEM(DM, Vec, Vec, void *); 170552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianActionFEM(DM, Mat, Vec, Vec, void *); 171552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianFEM(DM, Vec, Mat, Mat, MatStructure*,void *); 172552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetFEMIntegration(DM, 173552f7358SJed Brown PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 174552f7358SJed Brown const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 175552f7358SJed Brown void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 176552f7358SJed Brown void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]), 17702a80efeSMatthew G. Knepley PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 17802a80efeSMatthew G. Knepley const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 17902a80efeSMatthew G. Knepley void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 18002a80efeSMatthew G. Knepley void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), PetscScalar[]), 181552f7358SJed Brown PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], const PetscScalar[], 182552f7358SJed Brown const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 183552f7358SJed Brown void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 184552f7358SJed Brown void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 185552f7358SJed Brown void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 186552f7358SJed Brown void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]), 187552f7358SJed Brown PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 188552f7358SJed Brown const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 189552f7358SJed Brown void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 190552f7358SJed Brown void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 191552f7358SJed Brown void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 192552f7358SJed Brown void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[])); 193552f7358SJed Brown #endif 194