xref: /petsc/include/petscdmmoab.h (revision 01314a76e9272c2e0f76df6c826eef91022873fd)
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 
8*01314a76SVijay 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)
17*01314a76SVijay 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)
18*01314a76SVijay Mahadevan #define MBERRVM(mbif,msg,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 :: %s",(PetscErrorCode)rval,msg,emsg.c_str());} } while(0)
191d72bce8STim Tautges 
201d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *moab);
216a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *moab);
221d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm);
231d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm);
241d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *iface);
256a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **iface);
266a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetRange(DM dm,moab::Range *range);
276a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetRange(DM dm,moab::Range **range);
281d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag);
291d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag);
306a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs);
316a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs);
321d72bce8STim Tautges 
33efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArrayRead(DM,Vec,void*);
34efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArrayRead(DM,Vec,void*);
35efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArray(DM,Vec,void*);
36efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArray(DM,Vec,void*);
37efd17f3eSVijay Mahadevan 
386a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateVector(DM dm, moab::Tag tag,PetscInt tag_size,moab::Range *range,PetscBool serial, PetscBool destroy_tag,Vec *X);
396a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateMatrix(DM dm, MatType mtype,Mat *J);
401d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecTag(Vec vec,moab::Tag *tag);
411d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecRange(Vec vec,moab::Range *range);
4251d15aeeSVijay Mahadevan 
4351d15aeeSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt nele, PetscInt nghost, DM *dm);
4451d15aeeSVijay Mahadevan 
451d72bce8STim Tautges #endif
46