1552f7358SJed Brown /* 2552f7358SJed Brown DMPlex, for parallel unstructured distributed mesh problems. 3552f7358SJed Brown */ 426bd1501SBarry Smith #if !defined(PETSCDMPLEX_H) 526bd1501SBarry Smith #define PETSCDMPLEX_H 6552f7358SJed Brown 7224ef0b1SMatthew Knepley #include <petscsection.h> 8552f7358SJed Brown #include <petscdm.h> 986fe8405SMatthew G. Knepley #include <petscdmplextypes.h> 10a0845e3aSMatthew G. Knepley #include <petscdt.h> 11a0845e3aSMatthew G. Knepley #include <petscfe.h> 120bacfa0aSMatthew G. Knepley #include <petscfv.h> 13799f573fSMatthew G. Knepley #include <petscsftypes.h> 14c330f8ffSToby Isaac #include <petscdmfield.h> 15552f7358SJed Brown 1677623264SMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCPARTITIONER_CLASSID; 1777623264SMatthew G. Knepley 1877623264SMatthew G. Knepley /*J 1977623264SMatthew G. Knepley PetscPartitionerType - String with the name of a PETSc graph partitioner 2077623264SMatthew G. Knepley 2177623264SMatthew G. Knepley Level: beginner 2277623264SMatthew G. Knepley 2377623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitioner 2477623264SMatthew G. Knepley J*/ 2577623264SMatthew G. Knepley typedef const char *PetscPartitionerType; 2677623264SMatthew G. Knepley #define PETSCPARTITIONERCHACO "chaco" 2777623264SMatthew G. Knepley #define PETSCPARTITIONERPARMETIS "parmetis" 28137cd93aSLisandro Dalcin #define PETSCPARTITIONERPTSCOTCH "ptscotch" 2977623264SMatthew G. Knepley #define PETSCPARTITIONERSHELL "shell" 30555a9cf8SMatthew G. Knepley #define PETSCPARTITIONERSIMPLE "simple" 316462276dSToby Isaac #define PETSCPARTITIONERGATHER "gather" 32de68236aSVaclav Hapla #define PETSCPARTITIONERMATPARTITIONING "matpartitioning" 3377623264SMatthew G. Knepley 3477623264SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscPartitionerList; 3577623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate(MPI_Comm, PetscPartitioner *); 3677623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *); 3777623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetType(PetscPartitioner, PetscPartitionerType); 3877623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerGetType(PetscPartitioner, PetscPartitionerType *); 3977623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetUp(PetscPartitioner); 4077623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner); 418aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscPartitionerViewFromOptions(PetscPartitioner A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);} 4277623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerView(PetscPartitioner, PetscViewer); 4377623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerRegister(const char [], PetscErrorCode (*)(PetscPartitioner)); 4477623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterDestroy(void); 4577623264SMatthew G. Knepley 46f8987ae8SMichael Lange PETSC_EXTERN PetscErrorCode PetscPartitionerPartition(PetscPartitioner, DM, PetscSection, IS *); 4777623264SMatthew G. Knepley 485680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner, PetscInt, const PetscInt[], const PetscInt[]); 49df8115dbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner, PetscBool); 50df8115dbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner, PetscBool *); 5177623264SMatthew G. Knepley 52de68236aSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscPartitionerMatPartitioningGetMatPartitioning(PetscPartitioner part, MatPartitioning *mp); 53de68236aSVaclav Hapla 54552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*); 5527c04023SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char [], PetscInt, DM *); 560fde5641SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*); 57fba955ccSBarry Smith PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListParallel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const PetscReal[], PetscSF *, DM *); 5806bbee04SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []); 598415267dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm, PetscInt, PetscBool, DM*); 60a9074c1eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetOptionsPrefix(DM, const char []); 61552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *); 62552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt); 63552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *); 64552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt); 65feb68dd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexAddConeSize(DM, PetscInt, PetscInt); 66552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]); 670ce7577fSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeTuple(DM, IS, PetscSection *, IS *); 68af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]); 69af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexRestoreConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]); 70af9eab45SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursiveVertices(DM, IS, IS *); 71552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]); 729f074e33SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt); 7377c88f5bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt); 74552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]); 75552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]); 76552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *); 77552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt); 78552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]); 79552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]); 80552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt); 81552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *); 828cb4d582SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *); 83552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]); 84552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]); 85552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *); 86552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM); 87552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexStratify(DM); 884e3744c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexEqual(DM,DM,PetscBool*); 89cb11e54dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexReverseCell(DM, PetscInt); 90fa01033eSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexOrientCell(DM,PetscInt,PetscInt,const PetscInt[]); 9127d0e99aSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCompareOrientations(DM, PetscInt, PetscInt, const PetscInt [], PetscInt *, PetscBool *); 92d27a0f52SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexOrient(DM); 93e101f074SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool); 94552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*); 95081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,void*); 96a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*); 97a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*); 98a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*); 99081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*,void*); 100552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*); 101552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*); 102a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*); 103a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM,PetscInt,PetscInt,PetscInt*,PetscInt*); 104a89cf0ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM,PetscInt,PetscInt,PetscScalar*,void*); 105081a2d76SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM,PetscInt,PetscInt,const PetscScalar*, void*); 106552f7358SJed Brown 107*a7748839SVaclav Hapla /* Topological interpolation */ 108*a7748839SVaclav Hapla PETSC_EXTERN const char * const DMPlexInterpolatedFlags[]; 109*a7748839SVaclav Hapla typedef enum { 110*a7748839SVaclav Hapla DMPLEX_INTERPOLATED_UNINITIALIZED = -1, 111*a7748839SVaclav Hapla DMPLEX_INTERPOLATED_NONE = 0, 112*a7748839SVaclav Hapla DMPLEX_INTERPOLATED_PARTIAL, 113*a7748839SVaclav Hapla DMPLEX_INTERPOLATED_MIXED, 114*a7748839SVaclav Hapla DMPLEX_INTERPOLATED_FULL 115*a7748839SVaclav Hapla } DMPlexInterpolatedFlag; 116*a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *); 117*a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *); 118*a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexInterpolatePointSF(DM, PetscSF); 119*a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsInterpolated(DM, DMPlexInterpolatedFlag *); 120*a7748839SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexIsInterpolatedCollective(DM, DMPlexInterpolatedFlag *); 121*a7748839SVaclav Hapla 1228b51c812SSatish Balay PETSC_EXTERN PetscErrorCode DMPlexFilter(DM, DMLabel, PetscInt, DM *); 123552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *); 124552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *); 125ef48cebcSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *); 12608a22f4bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRankField(DM, Vec *); 12718e14f0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateLabelField(DM, DMLabel, Vec *); 128552f7358SJed Brown 12975d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *); 1303ded2ed9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *); 13175d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *); 13275d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *); 13375d3a19aSMatthew G. Knepley 13475d3a19aSMatthew G. Knepley /* Topological Operations */ 135552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 136552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 137552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 138552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 139552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 140552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **); 141552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 142552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]); 143552f7358SJed Brown 14475d3a19aSMatthew G. Knepley /* Mesh Generation */ 145552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *); 1463a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegister(const char[],PetscErrorCode (*)(DM,PetscBool,DM*),PetscErrorCode (*)(DM,double*,DM*),PetscInt); 1473a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegisterAll(void); 1483a074057SBarry Smith PETSC_EXTERN PetscErrorCode DMPlexGenerateRegisterDestroy(void); 149930319cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM); 1501df5d5c5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscBool, PetscReal, DM *); 15126492d91SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSquareBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 152552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []); 153768d5fceSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *); 15465a81367SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm, PetscInt, PetscBool, DM *); 155dbc1dc17SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm, PetscInt, DMBoundaryType, DM *); 15624119c2aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm, PetscInt, PetscBool, DM *); 15700dabe28SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, PetscBool, DM *); 15800dabe28SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexExtrude(DM, PetscInt, PetscReal, PetscBool, PetscBool, DM *); 159552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *); 160459e96c1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInvertCell(PetscInt, PetscInt, int []); 161068a5610SStefano Zampini 162ca8062c8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM); 16325c50c26SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscInt); 16425c50c26SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscInt); 165bb6a34a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckGeometry(DM); 16603da9461SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckPointSF(DM); 167f84a5eb8SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCheckConesConformOnInterfaces(DM); 16843fa8764SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCheckCellShape(DM, PetscBool, PetscReal); 16975d3a19aSMatthew G. Knepley 170d9deefdfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *); 171776fd405SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *); 172d9deefdfSMatthew G. Knepley 173d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], PetscBool, DM *); 174d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *); 175d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char [], PetscBool, DM *); 176d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *); 177d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *); 178d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *); 179d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char [], PetscBool, DM *); 180d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *); 181d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char [], PetscBool, DM *); 182d6218766SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm, const char [], PetscBool, DM *); 183f2801cd6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm, const char [], PetscBool, DM *); 184d6218766SMatthew G. Knepley 18575d3a19aSMatthew G. Knepley /* Mesh Partitioning and Distribution */ 186552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **); 1875680f57bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *); 18871bb2955SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner); 189ac17c9edSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, const char[], PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *); 1903fa7752dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **, IS *); 191552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *); 1921fd9873aSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel); 1935abbe4feSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel); 19424d039d7SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel); 195be200f8dSMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelPropagate(DM, DMLabel); 196aa3148a8SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *); 19798ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetPartitionBalance(DM, PetscBool); 19898ba2d7fSLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetPartitionBalance(DM, PetscBool *); 19980cf41d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF*, DM*); 200b9f40539SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *); 201cb54e036SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexGetOverlap(DM, PetscInt *); 202552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM,PetscSF,PetscSection,Vec,PetscSection,Vec); 203b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *); 204552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM,PetscSF,PetscSection,MPI_Datatype,void*,PetscSection,void**); 2058b879b41SFlorian Wechsung PETSC_EXTERN PetscErrorCode DMPlexRebalanceSharedPoints(DM, PetscInt, PetscBool, PetscBool, PetscBool*); 20645800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM); 207a13df41bSStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetGatherDM(DM, PetscSF*, DM*); 208a13df41bSStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetRedundantDM(DM, PetscSF*, DM*); 2093c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUser(DM, PetscErrorCode (*)(DM,PetscInt,PetscInt*,PetscInt[],void*),void*); 2103c1f0c11SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUser(DM, PetscErrorCode (**)(DM,PetscInt,PetscInt*,PetscInt[],void*),void**); 21170034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseCone(DM,PetscBool); 21270034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseCone(DM,PetscBool*); 21370034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseClosure(DM,PetscBool); 21470034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseClosure(DM,PetscBool*); 215a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM,PetscBool); 216a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM,PetscBool*); 21770034214SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]); 218f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetMigrationSF(DM, PetscSF); 219f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetMigrationSF(DM, PetscSF *); 22075d3a19aSMatthew G. Knepley 221c99f2573SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, DMLabel, IS *); 222f5cedc29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *); 223f5cedc29SMatthew G. Knepley 224b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *); 225b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *); 226b0a623aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *); 227dccdeccaSVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabel(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *); 22846f9b1c3SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *); 22945800214SMichael Lange PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *); 230b0a623aaSMatthew G. Knepley 23175d3a19aSMatthew G. Knepley /* Submesh Support */ 232158acfadSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, PetscBool, DM*); 2337db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel, DMLabel *, DMLabel *, DM *, DM *); 234552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*); 235552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel); 236552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *); 237559a1558SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetSubpoint(DM, PetscInt, PetscInt *); 23844171101SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetAuxiliaryPoint(DM, DM, PetscInt, PetscInt *); 239552f7358SJed Brown 2402be2b188SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel); 2410f66a230SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscBool, DM); 2426cf0e42fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel); 243f402d5e4SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexLabelClearCells(DM, DMLabel); 244552f7358SJed Brown 24575d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *); 24675d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal); 24775d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *); 24875d3a19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool); 2491e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRefinementFunction(DM, PetscErrorCode (**)(const PetscReal [], PetscReal *)); 2501e4cba6aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRefinementFunction(DM, PetscErrorCode (*)(const PetscReal [], PetscReal *)); 2512389894bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *); 2520aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetRegularRefinement(DM, PetscBool *); 2530aef6b92SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetRegularRefinement(DM, PetscBool); 254e5337592SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexRefineSimplexToTensor(DM, DM*); 255552f7358SJed Brown 25618ad9376SMatthew G. Knepley /* Support for cell-vertex meshes */ 25718ad9376SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *); 258f5d69827SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt [], PetscInt, PetscInt [], PetscInt [], PetscInt [], PetscBool *); 25918ad9376SMatthew G. Knepley 2600f0bf202SMatthew G. Knepley /* Geometry Support */ 2610f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *); 2620f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal); 263741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar[], PetscReal[]); 264741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar[], PetscReal[]); 265741bfc07SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt, PetscScalar[], PetscReal[]); 2660f0bf202SMatthew G. Knepley 267c4eade1cSMatthew G. Knepley /* Point Location */ 268c4eade1cSMatthew G. Knepley typedef struct _PetscGridHash *PetscGridHash; 269c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *); 270c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashEnlarge(PetscGridHash, const PetscScalar []); 271c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashSetGrid(PetscGridHash, const PetscInt [], const PetscReal []); 2721c6dfc3eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscScalar [], PetscInt [], PetscInt []); 273c4eade1cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscGridHashDestroy(PetscGridHash *); 2743985bb02SVaclav Hapla PETSC_EXTERN PetscErrorCode DMPlexFindVertices(DM, PetscInt, const PetscReal [], PetscReal, PetscInt []); 275c4eade1cSMatthew G. Knepley 276834e62ceSMatthew G. Knepley /* FVM Support */ 277cc08537eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal [], PetscReal []); 2780f0bf202SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *); 279856ac710SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *); 280b27d5b9eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetDataFVM(DM, PetscFV, Vec *, Vec *, DM *); 281834e62ceSMatthew G. Knepley 282552f7358SJed Brown /* FEM Support */ 283cf0b7c11SKarl Rupp PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFEM(DM, Vec *); 284bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec); 2851c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM, PetscReal, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], 286b278463cSMatthew G. Knepley PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, Vec ); 2871c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [], 288b278463cSMatthew G. Knepley void (*)(PetscInt, PetscInt, PetscInt, 289b278463cSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 290b278463cSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 291191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], 292b278463cSMatthew G. Knepley PetscScalar[]), 293b278463cSMatthew G. Knepley void *, Vec); 2941c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM, PetscReal, Vec, Vec, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt [], 295b278463cSMatthew G. Knepley PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *, Vec); 296e752be1aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, PetscInt, DMLabel); 297d7ddef95SMatthew G. Knepley 298baef929fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, DMLabel [], const PetscInt [], const PetscInt [], PetscInt, const PetscInt [], const IS [], const IS [], IS, PetscSection *); 299be36d101SStefano Zampini PETSC_EXTERN PetscErrorCode DMPlexGetSubdomainSection(DM, PetscSection*); 30075d3a19aSMatthew G. Knepley 3018e0841e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 302dfccc68fSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 303d6143a4eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCoordinatesToReference(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]); 3049d150b73SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceToCoordinates(DM, PetscInt, PetscInt, const PetscReal [], PetscReal[]); 305d6143a4eSToby Isaac 3067c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 3077c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]); 308552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode); 3097c1f9639SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 3106ecaa68aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscInt *); 31146bdb399SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexRestoreClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscInt *, PetscInt **,PetscInt *); 3120405ed22SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode); 3137c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]); 314ee72d991SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection); 315bc1eb3faSJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetClosurePermutationTensor(DM, PetscInt, PetscSection); 316552f7358SJed Brown 317552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *); 3187db7e0a7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DMLabel, DM *); 319552f7358SJed Brown 320770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 321770b213bSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt); 322552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *); 323552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt); 324552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer); 325552f7358SJed Brown 3262f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **); 3272f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **); 3282f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **); 3292f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **); 3302f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **); 3312f856554SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **); 3322f856554SMatthew G. Knepley 333552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *); 334552f7358SJed Brown PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal); 335552f7358SJed Brown 336552f7358SJed Brown typedef struct { 337552f7358SJed Brown DM dm; 33828d58a37SPierre Jolivet Vec u; /* The base vector for the Jacobian action J(u) x */ 339552f7358SJed Brown Mat J; /* Preconditioner for testing */ 340552f7358SJed Brown void *user; 341552f7358SJed Brown } JacActionCtx; 342552f7358SJed Brown 3433351dd3dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM, Vec); 344b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt); 345b29cfa1cSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt*); 346e29c0834SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexProjectFieldLocal(DM, Vec, 347e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 348e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 349e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 350e29c0834SMatthew G. Knepley PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec); 351574a98acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *); 3520163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]); 3530163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffVec(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, Vec); 354338f77d5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM, Vec, Vec, void *); 355b8feb594SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscScalar *, void *); 3569b6f715bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdIntegral(DM, Vec, DMLabel, PetscInt, const PetscInt[], 3579b6f715bSMatthew G. Knepley void (*)(PetscInt, PetscInt, PetscInt, 3589b6f715bSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 3599b6f715bSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 3609b6f715bSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 3619b6f715bSMatthew G. Knepley PetscScalar *, void *); 36268132eb9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorNested(DM, DM, Mat, void *); 36368132eb9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorGeneral(DM, DM, Mat, void *); 3641555c271SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeGradientClementInterpolant(DM, Vec, Vec); 3657c927364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *); 3662cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixNested(DM, DM, Mat, void *); 3672cb13a97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixGeneral(DM, DM, Mat, void *); 3681c41a8caSMatthew G. Knepley 369026175e5SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, MatNullSpace *); 370fd86e548SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBodies(DM, PetscInt, DMLabel, const PetscInt[], const PetscInt[], MatNullSpace *); 371026175e5SToby Isaac 3729f520fc2SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetSNESLocalFEM(DM,void *,void *,void *); 373bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM, Vec, void *); 3740f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *); 3750f2d7e86SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat, void *); 3767a73cf09SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianAction(DM, IS, PetscReal, PetscReal, Vec, Vec, Vec, Vec, void *); 377e3a53471SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdResidualSingle(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, Vec); 3781624f205SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeBdJacobianSingle(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, PetscReal, Mat, Mat); 3790f2d7e86SMatthew G. Knepley 380ef68eab9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeBoundary(DM, PetscReal, Vec, Vec, void *); 381924a1b8fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *); 382adbe6fbbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *); 383756a1f44SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTSComputeIJacobianFEM(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 384adbe6fbbSMatthew G. Knepley 3851c41a8caSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *); 3865d16530eSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReconstructGradientsFVM(DM,Vec,Vec); 387a68b90caSToby Isaac 388a17985deSToby Isaac /* anchors */ 389a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection*, IS*); 390a17985deSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS); 391d6a7ad0dSToby Isaac /* tree */ 392d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM); 393d6a7ad0dSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM*); 394dcbd3bf7SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *); 395da43764aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM*); 396b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt [], PetscInt []); 397b2f41788SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]); 398d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *); 399d961a43aSToby Isaac PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]); 4006f5f1567SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *); 4013b490a17SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorReferenceTree(DM, Mat *); 402bbcf1e96SToby Isaac PETSC_EXTERN PetscErrorCode DMPlexTransferVecTree(DM,Vec,DM,Vec,PetscSF,PetscSF,PetscInt *,PetscInt *,PetscBool,PetscReal); 403fa534816SMatthew G. Knepley 404fa534816SMatthew G. Knepley /* natural order */ 405fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM, PetscSection, PetscSF, PetscSF *); 406f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexSetGlobalToNaturalSF(DM, PetscSF); 407f94b4a02SBlaise Bourdin PETSC_EXTERN PetscErrorCode DMPlexGetGlobalToNaturalSF(DM, PetscSF *); 408fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalBegin(DM, Vec, Vec); 409fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalEnd(DM, Vec, Vec); 410fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalBegin(DM, Vec, Vec); 411fa534816SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalEnd(DM, Vec, Vec); 4120bbe5d31Sbarral 4130bbe5d31Sbarral /* mesh adaptation */ 4140f7e110dSbarral PETSC_EXTERN PetscErrorCode DMPlexAdapt(DM, Vec, const char [], DM *); 415ca3d3a14SMatthew G. Knepley 416ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGlobalToLocalBasis(DM, Vec); 417ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexLocalToGlobalBasis(DM, Vec); 418ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexCreateBasisRotation(DM, PetscReal, PetscReal, PetscReal); 419552f7358SJed Brown #endif 420