xref: /petsc/include/petscdmmoab.h (revision fc4180132621dce802f5cfc3099d116b2d2b65bb)
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 
24*fc418013SVijay Mahadevan #define MOAB_PARWOPTS_WRITE_PART "WRITE_PART"
25*fc418013SVijay Mahadevan #define MOAB_PARWOPTS_WRITE_FORMAT "FORMAT"
26*fc418013SVijay Mahadevan 
277023aa44SVijay Mahadevan typedef const char* MoabReadMode;
28*fc418013SVijay 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);
395eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost);
405eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range);
415eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range);
421d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag);
431d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag);
446a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs);
456a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs);
469088682fSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetBoundaryEntities(DM dm,moab::Range *bdvtx,moab::Range* bdfaces);
471d72bce8STim Tautges 
48efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArrayRead(DM,Vec,void*);
49efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArrayRead(DM,Vec,void*);
50efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecGetArray(DM,Vec,void*);
51efd17f3eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArray(DM,Vec,void*);
52efd17f3eSVijay Mahadevan 
536a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateVector(DM dm, moab::Tag tag,PetscInt tag_size,moab::Range *range,PetscBool serial, PetscBool destroy_tag,Vec *X);
546a36134cSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateMatrix(DM dm, MatType mtype,Mat *J);
551d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecTag(Vec vec,moab::Tag *tag);
561d72bce8STim Tautges PETSC_EXTERN PetscErrorCode DMMoabGetVecRange(Vec vec,moab::Range *range);
579088682fSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFieldVector(DM, PetscInt, Vec);
5851d15aeeSVijay Mahadevan 
595eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields);
605eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point, PetscInt field,PetscInt* dof);
615eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points, PetscInt field,PetscInt* dof);
625eb88e9dSVijay Mahadevan 
637023aa44SVijay Mahadevan /* discretization and assembly specific DMMoab interface functions */
647023aa44SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn);
657023aa44SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos);
667023aa44SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx,PetscBool* elem_on_boundary);
675eb88e9dSVijay Mahadevan 
687023aa44SVijay Mahadevan /* DM utility creation interface */
695eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm,PetscInt,PetscInt,PetscInt,DM*);
705eb88e9dSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoabLoadFromFile(MPI_Comm,PetscInt,const char*,const char*,DM*);
715eb88e9dSVijay Mahadevan 
721d72bce8STim Tautges #endif
73