11d72bce8STim Tautges #if !defined(__PETSCMOAB_H) 21d72bce8STim Tautges #define __PETSCMOAB_H 31d72bce8STim Tautges 41d72bce8STim Tautges #include <petscvec.h> 51d72bce8STim Tautges #include <petscmat.h> 61d72bce8STim Tautges #include <petscdm.h> 71d72bce8STim Tautges 801314a76SVijay Mahadevan #include <string> 91d72bce8STim Tautges #include <moab/Core.hpp> 101d72bce8STim Tautges #include <moab/ParallelComm.hpp> 111d72bce8STim Tautges 121d72bce8STim Tautges /* The MBERR macro is used to save typing. It checks a MOAB error code 131d72bce8STim Tautges * (rval) and calls SETERRQ if not MB_SUCCESS. A message (msg) can 141d72bce8STim Tautges * also be passed in. */ 156a36134cSVijay Mahadevan #define MBERR(msg,rval) do{if(rval != moab::MB_SUCCESS) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"MOAB ERROR (%i): %s",(PetscErrorCode)rval,msg);} while(0) 161d72bce8STim Tautges #define MBERRNM(rval) do{if(rval != moab::MB_SUCCESS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"MOAB ERROR (%i)",rval);} while(0) 1701314a76SVijay Mahadevan #define MBERRV(mbif,rval) do{if(rval != moab::MB_SUCCESS) { std::string emsg; mbif->get_last_error(emsg); SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"MOAB ERROR (%i): %s",(PetscErrorCode)rval,emsg.c_str());} } while(0) 185eb88e9dSVijay Mahadevan #define MBERRVM(mbif,msg,rval) do{if(rval != moab::MB_SUCCESS) { std::string emsg; mbif->get_last_error(emsg); SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"MOAB ERROR (%i): %s :: %s",(PetscErrorCode)rval,msg,emsg.c_str());} } while(0) 191d72bce8STim Tautges 207023aa44SVijay Mahadevan #define MOAB_PARROPTS_READ_PART "READ_PART" 217023aa44SVijay Mahadevan #define MOAB_PARROPTS_READ_DELETE "READ_DELETE" 227023aa44SVijay Mahadevan #define MOAB_PARROPTS_BCAST_DELETE "BCAST_DELETE" 237023aa44SVijay Mahadevan 24fc418013SVijay Mahadevan #define MOAB_PARWOPTS_WRITE_PART "WRITE_PART" 25fc418013SVijay Mahadevan #define MOAB_PARWOPTS_WRITE_FORMAT "FORMAT" 26fc418013SVijay Mahadevan 277023aa44SVijay Mahadevan typedef const char* MoabReadMode; 28fc418013SVijay Mahadevan typedef const char* MoabWriteMode; 297023aa44SVijay Mahadevan 301d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *moab); 316a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *moab); 329088682fSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabOutput(DM,const char*,const char*); 335eb88e9dSVijay Mahadevan 341d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm); 351d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm); 361d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *iface); 376a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **iface); 385eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalVertices(DM dm,moab::Range *range); 398d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetAllVertices(DM dm,moab::Range *local); 405eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost); 415eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range); 425eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range); 431d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag); 441d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag); 456a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs); 466a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs); 474920ab11SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDimension(DM dm,PetscInt *dim); 4869263071SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBoundaryEntities(DM dm,moab::Range *bdvtx,moab::Range* bdfaces,moab::Range* bdelems); 491d72bce8STim Tautges 50212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetSize(DM dm,PetscInt *ng); 51212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalSize(DM dm,PetscInt *nl,PetscInt *ng); 52212ad6d1SVijay Mahadevan 53efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArrayRead(DM,Vec,void*); 54efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArrayRead(DM,Vec,void*); 55efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArray(DM,Vec,void*); 56efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArray(DM,Vec,void*); 57efd17f3eSVijay Mahadevan 58*bb8f3634SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateVector(DM dm, moab::Tag tag,moab::Range *range,PetscBool serial, PetscBool destroy_tag,Vec *X); 596a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateMatrix(DM dm, MatType mtype,Mat *J); 601d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecTag(Vec vec,moab::Tag *tag); 611d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecRange(Vec vec,moab::Range *range); 62212ad6d1SVijay Mahadevan 639088682fSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFieldVector(DM, PetscInt, Vec); 64212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetGlobalFieldVector(DM, Vec); 6551d15aeeSVijay Mahadevan 665eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields); 675eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point, PetscInt field,PetscInt* dof); 685eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points, PetscInt field,PetscInt* dof); 69212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points, PetscInt field,PetscInt* dof); 70212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof); 71212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof); 72212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof); 73212ad6d1SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof); 745eb88e9dSVijay Mahadevan 758d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm,PetscInt** dof); 768d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm,PetscInt** dof); 778d8d51c8SVijay Mahadevan 787023aa44SVijay Mahadevan /* discretization and assembly specific DMMoab interface functions */ 797023aa44SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn); 808d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn, moab::EntityHandle **conn); 818d8d51c8SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabRestoreVertexConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn, moab::EntityHandle **conn); 827023aa44SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos); 8369263071SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabIsEntityOnBoundary(DM dm,const moab::EntityHandle ent,PetscBool* ent_on_boundary); 8469263071SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx); 855eb88e9dSVijay Mahadevan 867023aa44SVijay Mahadevan /* DM utility creation interface */ 875eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm,PetscInt,PetscInt,PetscInt,DM*); 885eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabLoadFromFile(MPI_Comm,PetscInt,const char*,const char*,DM*); 895eb88e9dSVijay Mahadevan 901d72bce8STim Tautges #endif 91