126bd1501SBarry Smith #if !defined(PETSCMOAB_H) 226bd1501SBarry Smith #define PETSCMOAB_H 31d72bce8STim Tautges 4b8ecf6d3SVijay Mahadevan #include <petscvec.h> /*I "petscvec.h" I*/ 5b8ecf6d3SVijay Mahadevan #include <petscmat.h> /*I "petscmat.h" I*/ 6b8ecf6d3SVijay Mahadevan #include <petscdm.h> /*I "petscdm.h" I*/ 763d025dbSVijay Mahadevan #include <petscdt.h> /*I "petscdt.h" I*/ 81d72bce8STim Tautges 901314a76SVijay Mahadevan #include <string> 10b8ecf6d3SVijay Mahadevan #include <moab/Core.hpp> /*I "moab/Core.hpp" I*/ 119daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 12b8ecf6d3SVijay Mahadevan #include <moab/ParallelComm.hpp> /*I "moab/ParallelComm.hpp" I*/ 139daf19fdSVijay Mahadevan #endif 141d72bce8STim Tautges 151d72bce8STim Tautges /* The MBERR macro is used to save typing. It checks a MOAB error code 161d72bce8STim Tautges * (rval) and calls SETERRQ if not MB_SUCCESS. A message (msg) can 171d72bce8STim Tautges * also be passed in. */ 182c71b3e2SJacob Faibussowitsch #define MBERR(msg, rval) PetscCheck(rval == moab::MB_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "MOAB ERROR (%i): %s", (PetscErrorCode)rval, msg) 192c71b3e2SJacob Faibussowitsch #define MBERRNM(rval) PetscCheck(rval == moab::MB_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "MOAB ERROR (%i)", rval) 20*9371c9d4SSatish Balay #define MBERRV(mbif, rval) \ 21*9371c9d4SSatish Balay do { \ 22*9371c9d4SSatish Balay if (rval != moab::MB_SUCCESS) { \ 23*9371c9d4SSatish Balay std::string emsg; \ 24*9371c9d4SSatish Balay mbif->get_last_error(emsg); \ 25*9371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "MOAB ERROR (%i): %s", (PetscErrorCode)rval, emsg.c_str()); \ 26*9371c9d4SSatish Balay } \ 27*9371c9d4SSatish Balay } while (0) 28*9371c9d4SSatish Balay #define MBERRVM(mbif, msg, rval) \ 29*9371c9d4SSatish Balay do { \ 30*9371c9d4SSatish Balay if (rval != moab::MB_SUCCESS) { \ 31*9371c9d4SSatish Balay std::string emsg; \ 32*9371c9d4SSatish Balay mbif->get_last_error(emsg); \ 33*9371c9d4SSatish Balay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "MOAB ERROR (%i): %s :: %s", (PetscErrorCode)rval, msg, emsg.c_str()); \ 34*9371c9d4SSatish Balay } \ 35*9371c9d4SSatish Balay } while (0) 361d72bce8STim Tautges 372e4e7c01SVijay Mahadevan /* define enums for options to read and write MOAB files in parallel */ 38*9371c9d4SSatish Balay typedef enum { 39*9371c9d4SSatish Balay READ_PART, 40*9371c9d4SSatish Balay READ_DELETE, 41*9371c9d4SSatish Balay BCAST_DELETE 42*9371c9d4SSatish Balay } MoabReadMode; 432e4e7c01SVijay Mahadevan static const char *const MoabReadModes[] = {"READ_PART", "READ_DELETE", "BCAST_DELETE", "MoabReadMode", "", 0}; 44*9371c9d4SSatish Balay typedef enum { 45*9371c9d4SSatish Balay WRITE_PART, 46*9371c9d4SSatish Balay FORMAT 47*9371c9d4SSatish Balay } MoabWriteMode; 482e4e7c01SVijay Mahadevan static const char *const MoabWriteModes[] = {"WRITE_PART", "FORMAT", "MoabWriteMode", "", 0}; 497023aa44SVijay Mahadevan 501d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *moab); 519daf19fdSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::Tag *ltog_tag, moab::Range *range, DM *moab); 529088682fSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabOutput(DM, const char *, const char *); 535eb88e9dSVijay Mahadevan 54b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetInterface(DM, moab::Interface *); 55b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetInterface(DM, moab::Interface **); 569daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 579daf19fdSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetParallelComm(DM, moab::ParallelComm **); 589daf19fdSVijay Mahadevan #endif 599daf19fdSVijay Mahadevan 60b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalVertices(DM, moab::Range *); 61b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetAllVertices(DM, moab::Range *local); 62b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalVertices(DM, const moab::Range **, const moab::Range **); 63b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalElements(DM, moab::Range *); 64b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalElements(DM, const moab::Range **); 65b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalToGlobalTag(DM, moab::Tag); 66b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalToGlobalTag(DM, moab::Tag *); 67b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetBlockSize(DM, PetscInt bs); 68b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBlockSize(DM, PetscInt *bs); 69b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetBlockFills(DM, const PetscInt *, const PetscInt *); 70755f3dfbSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetHierarchyLevel(DM, PetscInt *); 71b117cd09SVijay Mahadevan 724920ab11SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDimension(DM dm, PetscInt *dim); 73b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBoundaryEntities(DM dm, moab::Range *, moab::Range *, moab::Range *); 74b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetMaterialBlock(DM dm, const moab::EntityHandle, PetscInt *); 751d72bce8STim Tautges 7641dd5348SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetSize(DM dm, PetscInt *, PetscInt *); 7741dd5348SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalSize(DM dm, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 7800cc10feSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetOffset(DM dm, PetscInt *); 79212ad6d1SVijay Mahadevan 80efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArrayRead(DM, Vec, void *); 81efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArrayRead(DM, Vec, void *); 82efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArray(DM, Vec, void *); 83efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArray(DM, Vec, void *); 84efd17f3eSVijay Mahadevan 85351b8a77SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateVector(DM dm, moab::Tag tag, const moab::Range *range, PetscBool serial, PetscBool destroy_tag, Vec *X); 861d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecTag(Vec vec, moab::Tag *tag); 871d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecRange(Vec vec, moab::Range *range); 88212ad6d1SVijay Mahadevan 899088682fSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFieldVector(DM, PetscInt, Vec); 90212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetGlobalFieldVector(DM, Vec); 9151d15aeeSVijay Mahadevan 92304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateVertices(DM, const PetscReal *coords, PetscInt nverts, moab::Range *); 93304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateElement(DM, const moab::EntityType type, const moab::EntityHandle *conn, PetscInt nverts, moab::EntityHandle *elem); 94304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateSubmesh(DM dm, DM *newdm); 95304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabRenumberMeshEntities(DM dm); 96304006b3SVijay Mahadevan 9799fa7e03SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName); 985f7ae369SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName); 99f28b2503SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFieldNames(DM dm, PetscInt nfields, const char *fields[]); 1005eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDof(DM dm, moab::EntityHandle point, PetscInt field, PetscInt *dof); 1015eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofs(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt field, PetscInt *dof); 102212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofsLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt field, PetscInt *dof); 103212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofs(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof); 104212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof); 105212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsBlocked(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof); 106212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof); 1075eb88e9dSVijay Mahadevan 1088d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm, PetscInt **dof); 1098d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm, PetscInt **dof); 1108d8d51c8SVijay Mahadevan 1117023aa44SVijay Mahadevan /* discretization and assembly specific DMMoab interface functions */ 1127023aa44SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetElementConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt *nconn, const moab::EntityHandle **conn); 1138d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt *nconn, moab::EntityHandle **conn); 1148d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabRestoreVertexConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt *nconn, moab::EntityHandle **conn); 115cade3ad9SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexCoordinates(DM dm, PetscInt nconn, const moab::EntityHandle *conn, PetscReal *vpos); 11669263071SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabIsEntityOnBoundary(DM dm, const moab::EntityHandle ent, PetscBool *ent_on_boundary); 11769263071SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCheckBoundaryVertices(DM dm, PetscInt nconn, const moab::EntityHandle *cnt, PetscBool *isbdvtx); 1186d9eb265SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBoundaryMarkers(DM dm, const moab::Range **bdvtx, const moab::Range **bdelems, const moab::Range **bdfaces); 1195eb88e9dSVijay Mahadevan 12063d025dbSVijay Mahadevan /* TODO: Replace nverts/coords with just moab::EntityHandle -- can also eliminate dim */ 12163d025dbSVijay Mahadevan /* TODO: Replace quad/npts with PetscDT */ 122181f196bSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature); 123181f196bSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabFEMComputeBasis(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature, PetscReal *phypts, PetscReal *jxw, PetscReal *phi, PetscReal **dphi); 124181f196bSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *xphy, PetscReal *natparam, PetscReal *phi); 12563d025dbSVijay Mahadevan 1267023aa44SVijay Mahadevan /* DM utility creation interface */ 127c3dd00c7SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscReal *, PetscInt, PetscInt, DM *); 12849d66b22SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabLoadFromFile(MPI_Comm, PetscInt, PetscInt, const char *, const char *, DM *); 1295eb88e9dSVijay Mahadevan 130b117cd09SVijay Mahadevan /* Uniform refinement hierarchy interface */ 131b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees); 132b117cd09SVijay Mahadevan 1331d72bce8STim Tautges #endif 134