Actual source code: petscdmmoab.h
1: #ifndef PETSCMOAB_H
2: #define PETSCMOAB_H
4: #include <petscvec.h>
5: #include <petscmat.h>
6: #include <petscdm.h>
7: #include <petscdt.h>
9: #include <string>
10: #include <moab/Core.hpp> /*I "moab/Core.hpp" I*/
11: #ifdef MOAB_HAVE_MPI
12: #include <moab/ParallelComm.hpp> /*I "moab/ParallelComm.hpp" I*/
13: #endif
15: /* The MBERR macro is used to save typing. It checks a MOAB error code
16: * (rval) and calls SETERRQ if not MB_SUCCESS. A message (msg) can
17: * also be passed in. */
20: #define MBERRV(mbif, rval) \
21: do { \
22: if (rval != moab::MB_SUCCESS) { \
23: std::string emsg; \
24: mbif->get_last_error(emsg); \
25: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "MOAB ERROR (%i): %s", (PetscErrorCode)rval, emsg.c_str()); \
26: } \
27: } while (0)
28: #define MBERRVM(mbif, msg, rval) \
29: do { \
30: if (rval != moab::MB_SUCCESS) { \
31: std::string emsg; \
32: mbif->get_last_error(emsg); \
33: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "MOAB ERROR (%i): %s :: %s", (PetscErrorCode)rval, msg, emsg.c_str()); \
34: } \
35: } while (0)
37: /* define enums for options to read and write MOAB files in parallel */
38: typedef enum {
39: READ_PART,
40: READ_DELETE,
41: BCAST_DELETE
42: } MoabReadMode;
43: static const char *const MoabReadModes[] = {"READ_PART", "READ_DELETE", "BCAST_DELETE", "MoabReadMode", "", 0};
44: typedef enum {
45: WRITE_PART,
46: FORMAT
47: } MoabWriteMode;
48: static const char *const MoabWriteModes[] = {"WRITE_PART", "FORMAT", "MoabWriteMode", "", 0};
50: PETSC_EXTERN PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *moab);
51: PETSC_EXTERN PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::Tag *ltog_tag, moab::Range *range, DM *moab);
52: PETSC_EXTERN PetscErrorCode DMMoabOutput(DM, const char *, const char *);
54: PETSC_EXTERN PetscErrorCode DMMoabSetInterface(DM, moab::Interface *);
55: PETSC_EXTERN PetscErrorCode DMMoabGetInterface(DM, moab::Interface **);
56: #ifdef MOAB_HAVE_MPI
57: PETSC_EXTERN PetscErrorCode DMMoabGetParallelComm(DM, moab::ParallelComm **);
58: #endif
60: PETSC_EXTERN PetscErrorCode DMMoabSetLocalVertices(DM, moab::Range *);
61: PETSC_EXTERN PetscErrorCode DMMoabGetAllVertices(DM, moab::Range *local);
62: PETSC_EXTERN PetscErrorCode DMMoabGetLocalVertices(DM, const moab::Range **, const moab::Range **);
63: PETSC_EXTERN PetscErrorCode DMMoabSetLocalElements(DM, moab::Range *);
64: PETSC_EXTERN PetscErrorCode DMMoabGetLocalElements(DM, const moab::Range **);
65: PETSC_EXTERN PetscErrorCode DMMoabSetLocalToGlobalTag(DM, moab::Tag);
66: PETSC_EXTERN PetscErrorCode DMMoabGetLocalToGlobalTag(DM, moab::Tag *);
67: PETSC_EXTERN PetscErrorCode DMMoabSetBlockSize(DM, PetscInt bs);
68: PETSC_EXTERN PetscErrorCode DMMoabGetBlockSize(DM, PetscInt *bs);
69: PETSC_EXTERN PetscErrorCode DMMoabSetBlockFills(DM, const PetscInt *, const PetscInt *);
70: PETSC_EXTERN PetscErrorCode DMMoabGetHierarchyLevel(DM, PetscInt *);
72: PETSC_EXTERN PetscErrorCode DMMoabGetDimension(DM dm, PetscInt *dim);
73: PETSC_EXTERN PetscErrorCode DMMoabGetBoundaryEntities(DM dm, moab::Range *, moab::Range *, moab::Range *);
74: PETSC_EXTERN PetscErrorCode DMMoabGetMaterialBlock(DM dm, const moab::EntityHandle, PetscInt *);
76: PETSC_EXTERN PetscErrorCode DMMoabGetSize(DM dm, PetscInt *, PetscInt *);
77: PETSC_EXTERN PetscErrorCode DMMoabGetLocalSize(DM dm, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
78: PETSC_EXTERN PetscErrorCode DMMoabGetOffset(DM dm, PetscInt *);
80: PETSC_EXTERN PetscErrorCode DMMoabVecGetArrayRead(DM, Vec, void *);
81: PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArrayRead(DM, Vec, void *);
82: PETSC_EXTERN PetscErrorCode DMMoabVecGetArray(DM, Vec, void *);
83: PETSC_EXTERN PetscErrorCode DMMoabVecRestoreArray(DM, Vec, void *);
85: PETSC_EXTERN PetscErrorCode DMMoabCreateVector(DM dm, moab::Tag tag, const moab::Range *range, PetscBool serial, PetscBool destroy_tag, Vec *X);
86: PETSC_EXTERN PetscErrorCode DMMoabGetVecTag(Vec vec, moab::Tag *tag);
87: PETSC_EXTERN PetscErrorCode DMMoabGetVecRange(Vec vec, moab::Range *range);
89: PETSC_EXTERN PetscErrorCode DMMoabSetFieldVector(DM, PetscInt, Vec);
90: PETSC_EXTERN PetscErrorCode DMMoabSetGlobalFieldVector(DM, Vec);
92: PETSC_EXTERN PetscErrorCode DMMoabCreateVertices(DM, const PetscReal *coords, PetscInt nverts, moab::Range *);
93: PETSC_EXTERN PetscErrorCode DMMoabCreateElement(DM, const moab::EntityType type, const moab::EntityHandle *conn, PetscInt nverts, moab::EntityHandle *elem);
94: PETSC_EXTERN PetscErrorCode DMMoabCreateSubmesh(DM dm, DM *newdm);
95: PETSC_EXTERN PetscErrorCode DMMoabRenumberMeshEntities(DM dm);
97: PETSC_EXTERN PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName);
98: PETSC_EXTERN PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName);
99: PETSC_EXTERN PetscErrorCode DMMoabSetFieldNames(DM dm, PetscInt nfields, const char *fields[]);
100: PETSC_EXTERN PetscErrorCode DMMoabGetFieldDof(DM dm, moab::EntityHandle point, PetscInt field, PetscInt *dof);
101: PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofs(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt field, PetscInt *dof);
102: PETSC_EXTERN PetscErrorCode DMMoabGetFieldDofsLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt field, PetscInt *dof);
103: PETSC_EXTERN PetscErrorCode DMMoabGetDofs(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof);
104: PETSC_EXTERN PetscErrorCode DMMoabGetDofsLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof);
105: PETSC_EXTERN PetscErrorCode DMMoabGetDofsBlocked(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof);
106: PETSC_EXTERN PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof);
108: PETSC_EXTERN PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm, PetscInt **dof);
109: PETSC_EXTERN PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm, PetscInt **dof);
111: /* discretization and assembly specific DMMoab interface functions */
112: PETSC_EXTERN PetscErrorCode DMMoabGetElementConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt *nconn, const moab::EntityHandle **conn);
113: PETSC_EXTERN PetscErrorCode DMMoabGetVertexConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt *nconn, moab::EntityHandle **conn);
114: PETSC_EXTERN PetscErrorCode DMMoabRestoreVertexConnectivity(DM dm, moab::EntityHandle ehandle, PetscInt *nconn, moab::EntityHandle **conn);
115: PETSC_EXTERN PetscErrorCode DMMoabGetVertexCoordinates(DM dm, PetscInt nconn, const moab::EntityHandle *conn, PetscReal *vpos);
116: PETSC_EXTERN PetscErrorCode DMMoabIsEntityOnBoundary(DM dm, const moab::EntityHandle ent, PetscBool *ent_on_boundary);
117: PETSC_EXTERN PetscErrorCode DMMoabCheckBoundaryVertices(DM dm, PetscInt nconn, const moab::EntityHandle *cnt, PetscBool *isbdvtx);
118: PETSC_EXTERN PetscErrorCode DMMoabGetBoundaryMarkers(DM dm, const moab::Range **bdvtx, const moab::Range **bdelems, const moab::Range **bdfaces);
120: /* TODO: Replace nverts/coords with just moab::EntityHandle -- can also eliminate dim */
121: /* TODO: Replace quad/npts with PetscDT */
122: PETSC_EXTERN PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature);
123: PETSC_EXTERN PetscErrorCode DMMoabFEMComputeBasis(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature, PetscReal *phypts, PetscReal *jxw, PetscReal *phi, PetscReal **dphi);
124: PETSC_EXTERN PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *xphy, PetscReal *natparam, PetscReal *phi);
126: /* DM utility creation interface */
127: PETSC_EXTERN PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscReal *, PetscInt, PetscInt, DM *);
128: PETSC_EXTERN PetscErrorCode DMMoabLoadFromFile(MPI_Comm, PetscInt, PetscInt, const char *, const char *, DM *);
130: /* Uniform refinement hierarchy interface */
131: PETSC_EXTERN PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees);
133: #endif