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 <petscdm.h> 8a0845e3aSMatthew G. Knepley #include <petscdt.h> 9a0845e3aSMatthew G. Knepley #include <petscfe.h> 100bacfa0aSMatthew G. Knepley #include <petscfv.h> 11799f573fSMatthew G. Knepley #include <petscsftypes.h> 12552f7358SJed Brown 1377623264SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCPARTITIONER_CLASSID; 1477623264SMatthew G. Knepley 1577623264SMatthew G. Knepley /*J 1677623264SMatthew G. Knepley PetscPartitionerType - String with the name of a PETSc graph partitioner 1777623264SMatthew G. Knepley 1877623264SMatthew G. Knepley Level: beginner 1977623264SMatthew G. Knepley 2077623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitioner 2177623264SMatthew G. Knepley J*/ 2277623264SMatthew G. Knepley typedef const char *PetscPartitionerType; 2377623264SMatthew G. Knepley #define PETSCPARTITIONERCHACO "chaco" 2477623264SMatthew G. Knepley #define PETSCPARTITIONERPARMETIS "parmetis" 2577623264SMatthew G. Knepley #define PETSCPARTITIONERSHELL "shell" 26555a9cf8SMatthew G. Knepley #define PETSCPARTITIONERSIMPLE "simple" 2777623264SMatthew G. Knepley 2877623264SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscPartitionerList; 2977623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate(MPI_Comm, PetscPartitioner *); 3077623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *); 3177623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetType(PetscPartitioner, PetscPartitionerType); 3277623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerGetType(PetscPartitioner, PetscPartitionerType *); 3377623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetUp(PetscPartitioner); 3477623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner); 358aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscPartitionerViewFromOptions(PetscPartitioner A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);} 3677623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerView(PetscPartitioner, PetscViewer); 3777623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerRegister(const char [], PetscErrorCode (*)(PetscPartitioner)); 3877623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterDestroy(void); 3977623264SMatthew G. Knepley 40f8987ae8SMichael Lange PETSC_EXTERN PetscErrorCode PetscPartitionerPartition(PetscPartitioner, DM, PetscSection, IS *); 4177623264SMatthew G. Knepley 425680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner, PetscInt, const PetscInt[], const PetscInt[]); 4377623264SMatthew G. Knepley 44552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*); 4527c04023SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char [], PetscInt, DM *); 460fde5641SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*); 4706bbee04SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []); 488415267dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm, PetscInt, PetscBool, DM*); 49552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *); 50552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt); 51552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *); 52552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt); 53feb68dd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexAddConeSize(DM, PetscInt, PetscInt); 54552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]); 55552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]); 569f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt); 5777c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt); 58552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]); 59552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]); 60552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *); 61552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt); 62552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]); 63552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]); 64552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt); 65552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *); 668cb4d582SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *); 67552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]); 68552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]); 69552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *); 70552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM); 71552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexStratify(DM); 724e3744c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexEqual(DM,DM,PetscBool*); 73cb11e54dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexReverseCell(DM, PetscInt); 74d27a0f52SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrient(DM); 7575d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *); 764e3744c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *); 7775d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLoad(PetscViewer, DM); 78e101f074SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool); 79552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*); 80552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,const void*); 81a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*); 82a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*); 83a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*); 841ce3176fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*,const void*); 85552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*); 86552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*); 87a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*); 88a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*); 89a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*); 9033879625SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*,const void*); 91552f7358SJed Brown 92552f7358SJed Brown /*S 93552f7358SJed Brown DMLabel - Object which encapsulates a subset of the mesh from this DM 94552f7358SJed Brown 95552f7358SJed Brown Level: developer 96552f7358SJed Brown 97552f7358SJed Brown Concepts: grids, grid refinement 98552f7358SJed Brown 99552f7358SJed Brown .seealso: DM, DMPlexCreate(), DMPlexCreateLabel() 100552f7358SJed Brown S*/ 101552f7358SJed Brown typedef struct _n_DMLabel *DMLabel; 102552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelCreate(const char [], DMLabel *); 103552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelView(DMLabel, PetscViewer); 104552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelDestroy(DMLabel *); 105f807ecf7SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelDuplicate(DMLabel, DMLabel *); 10685de2ee3SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelGetName(DMLabel, const char **); 107552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetValue(DMLabel, PetscInt, PetscInt *); 108552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelSetValue(DMLabel, PetscInt, PetscInt); 109552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelClearValue(DMLabel, PetscInt, PetscInt); 110c52d6403SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLabelInsertIS(DMLabel, IS, PetscInt); 111552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetNumValues(DMLabel, PetscInt *); 112a9e27379SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLabelGetStratumBounds(DMLabel, PetscInt, PetscInt *, PetscInt *); 113552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetValueIS(DMLabel, IS *); 1140225b034SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLabelStratumHasPoint(DMLabel, PetscInt, PetscInt, PetscBool *); 115552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetStratumSize(DMLabel, PetscInt, PetscInt *); 116552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelGetStratumIS(DMLabel, PetscInt, IS *); 117552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMLabelClearStratum(DMLabel, PetscInt); 11891d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelCreateIndex(DMLabel, PetscInt, PetscInt); 11991d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelDestroyIndex(DMLabel); 12046b85170SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLabelHasValue(DMLabel, PetscInt, PetscBool *); 12191d07b70SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLabelHasPoint(DMLabel, PetscInt, PetscBool *); 122ae932cb5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLabelFilter(DMLabel, PetscInt, PetscInt); 123f186335cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLabelPermute(DMLabel, IS, DMLabel *); 1242eb1fa14SMichael Lange PETSC_EXTERN PetscErrorCode DMLabelDistribute(DMLabel, PetscSF, DMLabel *); 12533aee3d3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLabelConvertToSection(DMLabel, PetscSection *, IS *); 126552f7358SJed Brown 127552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateLabel(DM, const char []); 128552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelValue(DM, const char[], PetscInt, PetscInt *); 129552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetLabelValue(DM, const char[], PetscInt, PetscInt); 130552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexClearLabelValue(DM, const char[], PetscInt, PetscInt); 131552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelSize(DM, const char[], PetscInt *); 132552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelIdIS(DM, const char[], IS *); 133552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetStratumSize(DM, const char [], PetscInt, PetscInt *); 134552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetStratumIS(DM, const char [], PetscInt, IS *); 135552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexClearLabelStratum(DM, const char[], PetscInt); 136a85475f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetLabelOutput(DM, const char[], PetscBool *); 137a85475f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetLabelOutput(DM, const char[], PetscBool); 138552f7358SJed Brown PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection, PetscSF, PetscBool, DMLabel, PetscInt, PetscSection *); 139552f7358SJed Brown 140552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetNumLabels(DM, PetscInt *); 141552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabelName(DM, PetscInt, const char **); 142552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexHasLabel(DM, const char [], PetscBool *); 143552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetLabel(DM, const char *, DMLabel *); 144a85475f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetLabelByNum(DM, PetscInt, DMLabel *); 145552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexAddLabel(DM, DMLabel); 146552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRemoveLabel(DM, const char [], DMLabel *); 147552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *); 148552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *); 149ef48cebcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *); 150552f7358SJed Brown 15175d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *); 1523ded2ed9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *); 15375d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *); 15475d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 15575d3a19aSMatthew G. Knepley 15675d3a19aSMatthew G. Knepley /* Topological Operations */ 157552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 158552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 159552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 160552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 161552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 162552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 163552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 164552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 165552f7358SJed Brown 16675d3a19aSMatthew G. Knepley /* Mesh Generation */ 167552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *); 168930319cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM); 1695c386225SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyLabels(DM, DM); 1701df5d5c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscBool, PetscReal, DM *); 17126492d91SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSquareBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 172552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 1731a77d578SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSquareMesh(DM, const PetscReal [], const PetscReal [], const PetscInt [], DMBoundaryType, DMBoundaryType); 174552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, DM *); 175bff4a2f0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm, PetscInt, const PetscInt[], DMBoundaryType, DMBoundaryType, DMBoundaryType, DM *); 176552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *); 177459e96c1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInvertCell(PetscInt, PetscInt, int []); 178f05c90a3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLocalizeCoordinate(DM, const PetscScalar[], PetscScalar[]); 179d49fcaecSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLocalizeCoordinates(DM); 180ca8062c8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM); 18158723a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscBool, PetscInt); 1829bf0dad6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscBool, PetscInt); 18375d3a19aSMatthew G. Knepley 184d9deefdfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *); 185776fd405SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *); 186d9deefdfSMatthew G. Knepley 18775d3a19aSMatthew G. Knepley /* Mesh Partitioning and Distribution */ 188552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **); 1895680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *); 19071bb2955SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner); 191ac17c9edSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, const char[], PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *); 19245800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **); 193552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *); 1941fd9873aSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel); 1955abbe4feSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel); 19624d039d7SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel); 197aa3148a8SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *); 19880cf41d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF*, DM*); 199b9f40539SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *); 200552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM,PetscSF,PetscSection,Vec,PetscSection,Vec); 201b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *); 202552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM,PetscSF,PetscSection,MPI_Datatype,void*,PetscSection,void**); 20345800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM); 20470034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseCone(DM,PetscBool); 20570034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseCone(DM,PetscBool*); 20670034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseClosure(DM,PetscBool); 20770034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseClosure(DM,PetscBool*); 208a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM,PetscBool); 209a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM,PetscBool*); 21070034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]); 21175d3a19aSMatthew G. Knepley 212f5cedc29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, IS *); 213f5cedc29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *); 214f5cedc29SMatthew G. Knepley 215b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *); 216b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *); 217b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *); 218a9f1d5b2SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateOverlap(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *); 21946f9b1c3SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *); 22045800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *); 221b0a623aaSMatthew G. Knepley 22275d3a19aSMatthew G. Knepley /* Submesh Support */ 2233cf72582SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, DM*); 2243cf72582SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel *, DM *); 225552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*); 226552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel); 227552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *); 228552f7358SJed Brown 22909f723d9SJed Brown PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, DMLabel); 2302be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel); 2310f66a230SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscBool, DM); 2326cf0e42fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel); 233552f7358SJed Brown 23475d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *); 23575d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal); 23675d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *); 23775d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool); 2383913d7c8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCoarseDM(DM, DM *); 2393913d7c8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetCoarseDM(DM, DM); 2402389894bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *); 241552f7358SJed Brown 24218ad9376SMatthew G. Knepley /* Support for cell-vertex meshes */ 24318ad9376SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *); 244f5d69827SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *); 24518ad9376SMatthew G. Knepley 2460f0bf202SMatthew G. Knepley /* Geometry Support */ 2470f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *); 2480f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal); 2490f0bf202SMatthew G. Knepley 250*c4eade1cSMatthew G. Knepley /* Point Location */ 251*c4eade1cSMatthew G. Knepley typedef struct _PetscGridHash *PetscGridHash; 252*c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *); 253*c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashEnlarge(PetscGridHash, const PetscScalar []); 254*c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashSetGrid(PetscGridHash, const PetscInt [], const PetscReal []); 255*c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscReal [], PetscInt [], PetscInt []); 256*c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashDestroy(PetscGridHash *); 257*c4eade1cSMatthew G. Knepley 258834e62ceSMatthew G. Knepley /* FVM Support */ 259cc08537eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []); 2600f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *); 261856ac710SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *); 262834e62ceSMatthew G. Knepley 263552f7358SJed Brown /* FEM Support */ 264d7ddef95SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, Vec, PetscReal, Vec, Vec, Vec); 265d7ddef95SMatthew G. Knepley 266ab1d0545SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, PetscInt, PetscInt, const PetscInt [], const PetscInt [], PetscInt, const PetscInt [], const IS [], const IS [], IS, PetscSection *); 26775d3a19aSMatthew G. Knepley 2688e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 2698e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscFE, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 270c0d900a5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFEM(DM, Vec *); 2717c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 2727c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 273552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode); 2747c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 2750405ed22SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 2767c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]); 277ee72d991SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection); 278552f7358SJed Brown 279ca522641SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], PetscBool, DM *); 280552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *); 28133751fbdSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char [], PetscBool, DM *); 2825bb3eff3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *); 28344cd5272SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *); 284331830f3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *); 2857d282ae0SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char [], PetscBool, DM *); 2862f0bd6dcSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *); 2872f0bd6dcSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char [], PetscBool, DM *); 288552f7358SJed Brown 289552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *); 290a89082eeSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DM *); 291552f7358SJed Brown 292770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 293770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt); 294552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *); 295552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt); 296552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer); 297552f7358SJed Brown 298552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *); 299552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal); 300552f7358SJed Brown 3016c73c22cSMatthew G. Knepley typedef struct _n_Boundary *DMBoundary; 3020e0f25f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexAddBoundary(DM, PetscBool, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(), PetscInt, const PetscInt *, void *); 3036c73c22cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetNumBoundary(DM, PetscInt *); 3040e0f25f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetBoundary(DM, PetscInt, PetscBool *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(), PetscInt *, const PetscInt **, void **); 3050225b034SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexIsBoundaryPoint(DM, PetscInt, PetscBool *); 3068e136ac0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyBoundary(DM, DM); 3076c73c22cSMatthew G. Knepley 308552f7358SJed Brown typedef struct { 309552f7358SJed Brown DM dm; 310552f7358SJed Brown Vec u; /* The base vector for the Jacbobian action J(u) x */ 311552f7358SJed Brown Mat J; /* Preconditioner for testing */ 312552f7358SJed Brown void *user; 313552f7358SJed Brown } JacActionCtx; 314552f7358SJed Brown 3153351dd3dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM, Vec); 316b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt); 317b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt*); 318ad917190SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFunction(DM, PetscErrorCode (**)(PetscInt, const PetscReal [], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec); 319ad917190SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFunctionLocal(DM, PetscErrorCode (**)(PetscInt, const PetscReal [], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec); 3201ee2f23fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFieldLocal(DM, Vec, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode, Vec); 321ad917190SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2Diff(DM, PetscErrorCode (**)(PetscInt, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *); 322ad917190SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2GradientDiff(DM, PetscErrorCode (**)(PetscInt, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *); 323ad917190SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscErrorCode (**)(PetscInt, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]); 3240f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscReal *, void *); 325bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorFEM(DM, DM, Mat, void *); 3267c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *); 3271c41a8caSMatthew G. Knepley 328026175e5SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, MatNullSpace *); 329026175e5SToby Isaac 3300f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *); 3310f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat,void *); 3320f2d7e86SMatthew G. Knepley 333924a1b8fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *); 334adbe6fbbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *); 335adbe6fbbSMatthew G. Knepley 3361c41a8caSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *); 337a68b90caSToby Isaac 338a17985deSToby Isaac /* anchors */ 339a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection*, IS*); 340a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS); 341d6a7ad0dSToby Isaac /* tree */ 342d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM); 343d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM*); 344dcbd3bf7SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *); 345da43764aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM*); 346b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt [], PetscInt []); 347b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]); 348d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *); 349d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]); 3506f5f1567SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *); 351552f7358SJed Brown #endif 352