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. */ 18*2c71b3e2SJacob Faibussowitsch #define MBERR(msg,rval) PetscCheck(rval == moab::MB_SUCCESS,PETSC_COMM_SELF,PETSC_ERR_LIB,"MOAB ERROR (%i): %s",(PetscErrorCode)rval,msg) 19*2c71b3e2SJacob Faibussowitsch #define MBERRNM(rval) PetscCheck(rval == moab::MB_SUCCESS,PETSC_COMM_SELF,PETSC_ERR_LIB,"MOAB ERROR (%i)",rval) 2098921bdaSJacob Faibussowitsch #define MBERRV(mbif,rval) do{if (rval != moab::MB_SUCCESS) { std::string emsg; mbif->get_last_error(emsg); SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MOAB ERROR (%i): %s",(PetscErrorCode)rval,emsg.c_str());} } while (0) 2198921bdaSJacob Faibussowitsch #define MBERRVM(mbif,msg,rval) do{if (rval != moab::MB_SUCCESS) { std::string emsg; mbif->get_last_error(emsg); SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"MOAB ERROR (%i): %s :: %s",(PetscErrorCode)rval,msg,emsg.c_str());} } while (0) 221d72bce8STim Tautges 232e4e7c01SVijay Mahadevan /* define enums for options to read and write MOAB files in parallel */ 242e4e7c01SVijay Mahadevan typedef enum {READ_PART, READ_DELETE, BCAST_DELETE} MoabReadMode; 252e4e7c01SVijay Mahadevan static const char *const MoabReadModes[] = {"READ_PART", "READ_DELETE", "BCAST_DELETE", "MoabReadMode", "", 0}; 262e4e7c01SVijay Mahadevan typedef enum {WRITE_PART, FORMAT} MoabWriteMode; 272e4e7c01SVijay Mahadevan static const char *const MoabWriteModes[] = {"WRITE_PART", "FORMAT", "MoabWriteMode", "", 0}; 287023aa44SVijay Mahadevan 291d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *moab); 309daf19fdSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::Tag *ltog_tag, moab::Range *range, DM *moab); 319088682fSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabOutput(DM, const char*, const char*); 325eb88e9dSVijay Mahadevan 33b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetInterface(DM, moab::Interface *); 34b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetInterface(DM, moab::Interface **); 359daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 369daf19fdSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetParallelComm(DM, moab::ParallelComm**); 379daf19fdSVijay Mahadevan #endif 389daf19fdSVijay Mahadevan 39b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalVertices(DM, moab::Range *); 40b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetAllVertices(DM, moab::Range *local); 41b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalVertices(DM, const moab::Range**, const moab::Range**); 42b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalElements(DM, moab::Range*); 43b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalElements(DM, const moab::Range**); 44b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalToGlobalTag(DM, moab::Tag); 45b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalToGlobalTag(DM, moab::Tag*); 46b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetBlockSize(DM, PetscInt bs); 47b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBlockSize(DM, PetscInt *bs); 48b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetBlockFills(DM, const PetscInt*, const PetscInt*); 49755f3dfbSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetHierarchyLevel(DM, PetscInt *); 50b117cd09SVijay Mahadevan 514920ab11SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDimension(DM dm, PetscInt *dim); 52b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBoundaryEntities(DM dm, moab::Range*, moab::Range*, moab::Range*); 53b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetMaterialBlock(DM dm, const moab::EntityHandle, PetscInt*); 541d72bce8STim Tautges 5541dd5348SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetSize(DM dm, PetscInt*, PetscInt*); 5641dd5348SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalSize(DM dm, PetscInt*, PetscInt*, PetscInt*, PetscInt*); 5700cc10feSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetOffset(DM dm, PetscInt*); 58212ad6d1SVijay Mahadevan 59efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArrayRead(DM, Vec, void*); 60efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArrayRead(DM, Vec, void*); 61efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArray(DM, Vec, void*); 62efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArray(DM, Vec, void*); 63efd17f3eSVijay Mahadevan 64351b8a77SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateVector(DM dm, moab::Tag tag, const moab::Range *range, PetscBool serial, PetscBool destroy_tag, Vec *X); 651d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecTag(Vec vec, moab::Tag *tag); 661d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecRange(Vec vec, moab::Range *range); 67212ad6d1SVijay Mahadevan 689088682fSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFieldVector(DM, PetscInt, Vec); 69212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetGlobalFieldVector(DM, Vec); 7051d15aeeSVijay Mahadevan 71304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateVertices(DM, const PetscReal* coords, PetscInt nverts, moab::Range*); 72304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateElement(DM, const moab::EntityType type, const moab::EntityHandle* conn, PetscInt nverts, moab::EntityHandle* elem); 73304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateSubmesh(DM dm, DM *newdm); 74304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabRenumberMeshEntities(DM dm); 75304006b3SVijay Mahadevan 7699fa7e03SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName); 775f7ae369SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName); 78f28b2503SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFieldNames(DM dm, PetscInt nfields, const char* fields[]); 795eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDof(DM dm, moab::EntityHandle point, PetscInt field, PetscInt* dof); 805eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofs(DM dm, PetscInt npoints, const moab::EntityHandle* points, PetscInt field, PetscInt* dof); 81212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofsLocal(DM dm, PetscInt npoints, const moab::EntityHandle* points, PetscInt field, PetscInt* dof); 82212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofs(DM dm, PetscInt npoints, const moab::EntityHandle* points, PetscInt* dof); 83212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsLocal(DM dm, PetscInt npoints, const moab::EntityHandle* points, PetscInt* dof); 84212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsBlocked(DM dm, PetscInt npoints, const moab::EntityHandle* points, PetscInt* dof); 85212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm, PetscInt npoints, const moab::EntityHandle* points, PetscInt* dof); 865eb88e9dSVijay Mahadevan 878d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm, PetscInt** dof); 888d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm, PetscInt** dof); 898d8d51c8SVijay Mahadevan 907023aa44SVijay Mahadevan /* discretization and assembly specific DMMoab interface functions */ 917023aa44SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetElementConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt* nconn, const moab::EntityHandle **conn); 928d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt* nconn, moab::EntityHandle **conn); 938d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabRestoreVertexConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt* nconn, moab::EntityHandle **conn); 94cade3ad9SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexCoordinates(DM dm, PetscInt nconn, const moab::EntityHandle *conn, PetscReal *vpos); 9569263071SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabIsEntityOnBoundary(DM dm, const moab::EntityHandle ent, PetscBool* ent_on_boundary); 9669263071SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCheckBoundaryVertices(DM dm, PetscInt nconn, const moab::EntityHandle *cnt, PetscBool* isbdvtx); 976d9eb265SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBoundaryMarkers(DM dm, const moab::Range **bdvtx, const moab::Range** bdelems, const moab::Range** bdfaces); 985eb88e9dSVijay Mahadevan 9963d025dbSVijay Mahadevan /* TODO: Replace nverts/coords with just moab::EntityHandle -- can also eliminate dim */ 10063d025dbSVijay Mahadevan /* TODO: Replace quad/npts with PetscDT */ 101181f196bSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature); 102181f196bSVijay 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); 103181f196bSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal* xphy, PetscReal* natparam, PetscReal* phi); 10463d025dbSVijay Mahadevan 1057023aa44SVijay Mahadevan /* DM utility creation interface */ 106c3dd00c7SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscReal*, PetscInt, PetscInt, DM*); 10749d66b22SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabLoadFromFile(MPI_Comm, PetscInt, PetscInt, const char*, const char*, DM*); 1085eb88e9dSVijay Mahadevan 109b117cd09SVijay Mahadevan /* Uniform refinement hierarchy interface */ 110b117cd09SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees); 111b117cd09SVijay Mahadevan 1121d72bce8STim Tautges #endif 113