xref: /petsc/src/dm/impls/moab/dmmoab.cxx (revision fd349b4114aa6c18db3d975e4ea5f0790da309c7)
11d72bce8STim Tautges #include <petsc-private/dmimpl.h> /*I  "petscdm.h"   I*/
21d72bce8STim Tautges #include "../src/vec/vec/impls/moab/vecmoabimpl.h" /*I  "petscvec.h"   I*/
31d72bce8STim Tautges 
41d72bce8STim Tautges #include <petscdmmoab.h>
51d72bce8STim Tautges #include "MBTagConventions.hpp"
61d72bce8STim Tautges 
71d72bce8STim Tautges typedef struct {
81d72bce8STim Tautges   PetscInt bs; /* Number of degrees of freedom on each entity, aka tag size in moab */
91d72bce8STim Tautges   PetscBool icreatedinstance; /* true if DM created moab instance internally, will destroy instance in DMDestroy */
101d72bce8STim Tautges   moab::ParallelComm *pcomm;
111d72bce8STim Tautges   moab::Interface *mbint;
121d72bce8STim Tautges   moab::Tag ltog_tag; /* moab supports "global id" tags, which are usually local to global numbering */
131d72bce8STim Tautges   moab::Range range;
141d72bce8STim Tautges } DM_Moab;
151d72bce8STim Tautges 
161d72bce8STim Tautges #undef __FUNCT__
17*fd349b41STim Tautges #define __FUNCT__ "DMCreateGlobalVector_Moab"
18*fd349b41STim Tautges PetscErrorCode DMCreateGlobalVector_Moab(DM dm,Vec *gvec)
19*fd349b41STim Tautges {
20*fd349b41STim Tautges   PetscErrorCode  ierr;
21*fd349b41STim Tautges   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
22*fd349b41STim Tautges 
23*fd349b41STim Tautges   PetscFunctionBegin;
24*fd349b41STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
25*fd349b41STim Tautges   PetscValidPointer(gvec,2);
26*fd349b41STim Tautges   PetscInt block_size = ((DM_Moab*)dm->data)->bs;
27*fd349b41STim Tautges   ierr = VecMoabCreate(dmmoab->mbint,dmmoab->pcomm,block_size,dmmoab->ltog_tag,dmmoab->range,PETSC_FALSE,PETSC_TRUE,gvec);CHKERRQ(ierr);
28*fd349b41STim Tautges   PetscFunctionReturn(0);
29*fd349b41STim Tautges }
30*fd349b41STim Tautges 
31*fd349b41STim Tautges 
32*fd349b41STim Tautges #undef __FUNCT__
33*fd349b41STim Tautges #define __FUNCT__ "DMCreateLocalVector_Moab"
34*fd349b41STim Tautges PetscErrorCode DMCreateLocalVector_Moab(DM dm,Vec *gvec)
35*fd349b41STim Tautges {
36*fd349b41STim Tautges   PetscErrorCode  ierr;
37*fd349b41STim Tautges   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
38*fd349b41STim Tautges 
39*fd349b41STim Tautges   PetscFunctionBegin;
40*fd349b41STim Tautges   PetscInt bs = 1;
41*fd349b41STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
42*fd349b41STim Tautges   PetscValidPointer(gvec,2);
43*fd349b41STim Tautges   ierr = VecMoabCreate(dmmoab->mbint,dmmoab->pcomm,bs,dmmoab->ltog_tag,dmmoab->range,PETSC_TRUE,PETSC_TRUE,gvec);CHKERRQ(ierr);
44*fd349b41STim Tautges   PetscFunctionReturn(0);
45*fd349b41STim Tautges }
46*fd349b41STim Tautges 
47*fd349b41STim Tautges #undef __FUNCT__
48*fd349b41STim Tautges #define __FUNCT__ "DMDestroy_Moab"
49*fd349b41STim Tautges PetscErrorCode DMDestroy_Moab(DM dm)
50*fd349b41STim Tautges {
51*fd349b41STim Tautges   PetscFunctionBegin;
52*fd349b41STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
53*fd349b41STim Tautges 
54*fd349b41STim Tautges   // Delete the DM_Moab:
55*fd349b41STim Tautges   if(dm->data) {
56*fd349b41STim Tautges     if (((DM_Moab*)dm->data)->icreatedinstance) {
57*fd349b41STim Tautges       delete ((DM_Moab*)dm->data)->mbint;
58*fd349b41STim Tautges       ((DM_Moab*)dm->data)->mbint = NULL;
59*fd349b41STim Tautges       ((DM_Moab*)dm->data)->pcomm = NULL;
60*fd349b41STim Tautges     }
61*fd349b41STim Tautges     delete (DM_Moab*)dm->data;
62*fd349b41STim Tautges     dm->data = NULL;
63*fd349b41STim Tautges   }
64*fd349b41STim Tautges   PetscFunctionReturn(0);
65*fd349b41STim Tautges }
66*fd349b41STim Tautges 
67*fd349b41STim Tautges 
68*fd349b41STim Tautges #undef __FUNCT__
69*fd349b41STim Tautges #define __FUNCT__ "DMInitialize_Moab"
70*fd349b41STim Tautges PetscErrorCode DMInitialize_Moab(DM dm)
71*fd349b41STim Tautges {
72*fd349b41STim Tautges   PetscFunctionBegin;
73*fd349b41STim Tautges 
74*fd349b41STim Tautges   // Create the DM_Moab and set dm->data
75*fd349b41STim Tautges   DM_Moab *dmmoab = new DM_Moab;
76*fd349b41STim Tautges   dmmoab->bs       = 0;
77*fd349b41STim Tautges   dmmoab->pcomm    = NULL;
78*fd349b41STim Tautges   dmmoab->mbint    = NULL;
79*fd349b41STim Tautges   dmmoab->ltog_tag = (moab::Tag)0;
80*fd349b41STim Tautges   dmmoab->icreatedinstance = PETSC_FALSE;
81*fd349b41STim Tautges   dm->data      = dmmoab;
82*fd349b41STim Tautges 
83*fd349b41STim Tautges     // initialize various functions
84*fd349b41STim Tautges   dm->ops->createglobalvector              = DMCreateGlobalVector_Moab;
85*fd349b41STim Tautges   dm->ops->createlocalvector               = DMCreateLocalVector_Moab;
86*fd349b41STim Tautges   dm->ops->destroy                         = DMDestroy_Moab;
87*fd349b41STim Tautges   PetscFunctionReturn(0);
88*fd349b41STim Tautges }
89*fd349b41STim Tautges 
90*fd349b41STim Tautges #undef __FUNCT__
911d72bce8STim Tautges #define __FUNCT__ "DMMoabCreate"
921d72bce8STim Tautges /*@
931d72bce8STim Tautges   DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance
941d72bce8STim Tautges 
951d72bce8STim Tautges   Collective on MPI_Comm
961d72bce8STim Tautges 
971d72bce8STim Tautges   Input Parameter:
981d72bce8STim Tautges . comm - The communicator for the DMMoab object
991d72bce8STim Tautges 
1001d72bce8STim Tautges   Output Parameter:
1011d72bce8STim Tautges . moab  - The DMMoab object
1021d72bce8STim Tautges 
1031d72bce8STim Tautges   Level: beginner
1041d72bce8STim Tautges 
1051d72bce8STim Tautges .keywords: DMMoab, create
1061d72bce8STim Tautges @*/
1071d72bce8STim Tautges PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *moab)
1081d72bce8STim Tautges {
1091d72bce8STim Tautges   PetscErrorCode ierr;
1101d72bce8STim Tautges 
1111d72bce8STim Tautges   PetscFunctionBegin;
1121d72bce8STim Tautges   PetscValidPointer(moab,2);
1131d72bce8STim Tautges   ierr = DMCreate(comm, moab);CHKERRQ(ierr);
1141d72bce8STim Tautges   ierr = DMSetType(*moab, DMMOAB);CHKERRQ(ierr);
1151d72bce8STim Tautges   PetscFunctionReturn(0);
1161d72bce8STim Tautges }
1171d72bce8STim Tautges 
1181d72bce8STim Tautges EXTERN_C_BEGIN
1191d72bce8STim Tautges #undef __FUNCT__
1201d72bce8STim Tautges #define __FUNCT__ "DMCreate_Moab"
1211d72bce8STim Tautges PetscErrorCode DMCreate_Moab(DM dm)
1221d72bce8STim Tautges {
1231d72bce8STim Tautges   DM_Moab        *moab;
1241d72bce8STim Tautges   PetscErrorCode ierr;
1251d72bce8STim Tautges 
1261d72bce8STim Tautges   PetscFunctionBegin;
1271d72bce8STim Tautges   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1281d72bce8STim Tautges   ierr     = PetscNewLog(dm, DM_Moab, &moab);CHKERRQ(ierr);
1291d72bce8STim Tautges   dm->data = moab;
1301d72bce8STim Tautges 
1311d72bce8STim Tautges   PetscFunctionReturn(0);
1321d72bce8STim Tautges }
1331d72bce8STim Tautges EXTERN_C_END
1341d72bce8STim Tautges 
1351d72bce8STim Tautges /*@
1361d72bce8STim Tautges   DMMoabCreateFromInstance - Creates a DMMoab object from an instance and (optionally) other data
1371d72bce8STim Tautges 
1381d72bce8STim Tautges   Collective on MPI_Comm
1391d72bce8STim Tautges 
1401d72bce8STim Tautges   Input Parameter:
1411d72bce8STim Tautges . comm - The communicator for the DMMoab object
1421d72bce8STim Tautges . moab - The MOAB Instance
1431d72bce8STim Tautges . pcomm - A ParallelComm; if NULL, creates one internally for the whole communicator
1441d72bce8STim Tautges . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag
1451d72bce8STim Tautges . range - If non-NULL, contains range of entities to which DOFs will be assigned
1461d72bce8STim Tautges 
1471d72bce8STim Tautges   Output Parameter:
1481d72bce8STim Tautges . moab  - The DMMoab object
1491d72bce8STim Tautges 
1501d72bce8STim Tautges   Level: beginner
1511d72bce8STim Tautges 
1521d72bce8STim Tautges .keywords: DMMoab, create
1531d72bce8STim Tautges @*/
1541d72bce8STim Tautges #undef __FUNCT__
1551d72bce8STim Tautges #define __FUNCT__ "DMMoabCreateFromInstance"
1561d72bce8STim Tautges PetscErrorCode DMMoabCreateFromInstance(MPI_Comm comm, moab::Interface *iface, moab::ParallelComm *pcomm, moab::Tag ltog_tag, moab::Range *range, DM *moab)
1571d72bce8STim Tautges {
1581d72bce8STim Tautges   PetscErrorCode ierr;
1591d72bce8STim Tautges 
1601d72bce8STim Tautges   PetscFunctionBegin;
1611d72bce8STim Tautges   PetscValidPointer(moab,2);
1621d72bce8STim Tautges   ierr = DMCreate(comm, moab);CHKERRQ(ierr);
1631d72bce8STim Tautges   ierr = DMSetType(*moab, DMMOAB);CHKERRQ(ierr);
1641d72bce8STim Tautges   ierr = DMInitialize_Moab(*moab);CHKERRQ(ierr);
1651d72bce8STim Tautges   ierr = DMMoabSetInterface(*moab, iface);CHKERRQ(ierr);
1661d72bce8STim Tautges   if (!pcomm) pcomm = new moab::ParallelComm(iface, comm);
1671d72bce8STim Tautges   ierr = DMMoabSetParallelComm(*moab, pcomm);CHKERRQ(ierr);
1681d72bce8STim Tautges   if (!ltog_tag) {
1691d72bce8STim Tautges     moab::ErrorCode merr = iface->tag_get_handle(GLOBAL_ID_TAG_NAME, ltog_tag);MBERRNM(merr);
1701d72bce8STim Tautges   }
1711d72bce8STim Tautges   ierr = DMMoabSetLocalToGlobalTag(*moab, ltog_tag);CHKERRQ(ierr);
1721d72bce8STim Tautges   if (range) {
1731d72bce8STim Tautges     ierr = DMMoabSetRange(*moab, *range);CHKERRQ(ierr);
1741d72bce8STim Tautges   }
1751d72bce8STim Tautges   PetscFunctionReturn(0);
1761d72bce8STim Tautges }
1771d72bce8STim Tautges 
1781d72bce8STim Tautges /*@
1791d72bce8STim Tautges   DMMoabCreateDMAndInstance - Creates a DMMoab object and a MOAB instance
1801d72bce8STim Tautges 
1811d72bce8STim Tautges   Collective on MPI_Comm
1821d72bce8STim Tautges 
1831d72bce8STim Tautges   Input Parameter:
1841d72bce8STim Tautges . comm - The communicator for the DMMoab object
1851d72bce8STim Tautges 
1861d72bce8STim Tautges   Output Parameter:
1871d72bce8STim Tautges . moab  - The DMMoab object
1881d72bce8STim Tautges 
1891d72bce8STim Tautges   Level: beginner
1901d72bce8STim Tautges 
1911d72bce8STim Tautges .keywords: DMMoab, create
1921d72bce8STim Tautges @*/
1931d72bce8STim Tautges #undef __FUNCT__
1941d72bce8STim Tautges #define __FUNCT__ "DMMoabCreateDMAndInstance"
1951d72bce8STim Tautges PetscErrorCode DMMoabCreateDMAndInstance(MPI_Comm comm, DM *dm)
1961d72bce8STim Tautges {
1971d72bce8STim Tautges   PetscErrorCode ierr;
1981d72bce8STim Tautges 
1991d72bce8STim Tautges   PetscFunctionBegin;
2001d72bce8STim Tautges   PetscValidPointer(dm,2);
2011d72bce8STim Tautges   PetscInt rank, nprocs;
2021d72bce8STim Tautges   MPI_Comm_rank(comm, &rank);
2031d72bce8STim Tautges   MPI_Comm_size(comm, &nprocs);
2041d72bce8STim Tautges   moab::Interface *iface = new moab::Core();
2051d72bce8STim Tautges   moab::ParallelComm *pcomm = new moab::ParallelComm(iface, comm);
2061d72bce8STim Tautges   ierr = DMMoabCreateFromInstance(comm, iface, pcomm, 0, NULL, dm);CHKERRQ(ierr);
2071d72bce8STim Tautges   ((DM_Moab*)(*dm)->data)->icreatedinstance = PETSC_TRUE;
2081d72bce8STim Tautges   PetscFunctionReturn(0);
2091d72bce8STim Tautges }
2101d72bce8STim Tautges 
2111d72bce8STim Tautges #undef __FUNCT__
2121d72bce8STim Tautges #define __FUNCT__ "DMMoabSetParallelComm"
2131d72bce8STim Tautges PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm)
2141d72bce8STim Tautges {
2151d72bce8STim Tautges   PetscFunctionBegin;
2161d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2171d72bce8STim Tautges   ((DM_Moab*)dm->data)->pcomm = pcomm;
2181d72bce8STim Tautges   ((DM_Moab*)dm->data)->mbint = pcomm->get_moab();
2191d72bce8STim Tautges   PetscFunctionReturn(0);
2201d72bce8STim Tautges }
2211d72bce8STim Tautges 
2221d72bce8STim Tautges 
2231d72bce8STim Tautges #undef __FUNCT__
2241d72bce8STim Tautges #define __FUNCT__ "DMMoabGetParallelComm"
2251d72bce8STim Tautges PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm)
2261d72bce8STim Tautges {
2271d72bce8STim Tautges   PetscFunctionBegin;
2281d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2291d72bce8STim Tautges   *pcomm = ((DM_Moab*)dm->data)->pcomm;
2301d72bce8STim Tautges   PetscFunctionReturn(0);
2311d72bce8STim Tautges }
2321d72bce8STim Tautges 
2331d72bce8STim Tautges 
2341d72bce8STim Tautges #undef __FUNCT__
2351d72bce8STim Tautges #define __FUNCT__ "DMMoabSetInterface"
2361d72bce8STim Tautges PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *iface)
2371d72bce8STim Tautges {
2381d72bce8STim Tautges   PetscFunctionBegin;
2391d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2401d72bce8STim Tautges   ((DM_Moab*)dm->data)->pcomm = NULL;
2411d72bce8STim Tautges   ((DM_Moab*)dm->data)->mbint = iface;
2421d72bce8STim Tautges   PetscFunctionReturn(0);
2431d72bce8STim Tautges }
2441d72bce8STim Tautges 
2451d72bce8STim Tautges 
2461d72bce8STim Tautges #undef __FUNCT__
2471d72bce8STim Tautges #define __FUNCT__ "DMMoabGetInterface"
2481d72bce8STim Tautges PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbint)
2491d72bce8STim Tautges {
2501d72bce8STim Tautges   PetscFunctionBegin;
2511d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2521d72bce8STim Tautges   *mbint = ((DM_Moab*)dm->data)->mbint;
2531d72bce8STim Tautges   PetscFunctionReturn(0);
2541d72bce8STim Tautges }
2551d72bce8STim Tautges 
2561d72bce8STim Tautges 
2571d72bce8STim Tautges #undef __FUNCT__
2581d72bce8STim Tautges #define __FUNCT__ "DMMoabSetRange"
2591d72bce8STim Tautges PetscErrorCode DMMoabSetRange(DM dm,moab::Range range)
2601d72bce8STim Tautges {
2611d72bce8STim Tautges   PetscFunctionBegin;
2621d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2631d72bce8STim Tautges   ((DM_Moab*)dm->data)->range = range;
2641d72bce8STim Tautges   PetscFunctionReturn(0);
2651d72bce8STim Tautges }
2661d72bce8STim Tautges 
2671d72bce8STim Tautges 
2681d72bce8STim Tautges #undef __FUNCT__
2691d72bce8STim Tautges #define __FUNCT__ "DMMoabGetRange"
2701d72bce8STim Tautges PetscErrorCode DMMoabGetRange(DM dm,moab::Range *range)
2711d72bce8STim Tautges {
2721d72bce8STim Tautges   PetscFunctionBegin;
2731d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2741d72bce8STim Tautges   *range = ((DM_Moab*)dm->data)->range;
2751d72bce8STim Tautges   PetscFunctionReturn(0);
2761d72bce8STim Tautges }
2771d72bce8STim Tautges 
2781d72bce8STim Tautges #undef __FUNCT__
2791d72bce8STim Tautges #define __FUNCT__ "DMMoabSetLocalToGlobalTag"
2801d72bce8STim Tautges PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag)
2811d72bce8STim Tautges {
2821d72bce8STim Tautges   PetscFunctionBegin;
2831d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2841d72bce8STim Tautges   ((DM_Moab*)dm->data)->ltog_tag = ltogtag;
2851d72bce8STim Tautges   PetscFunctionReturn(0);
2861d72bce8STim Tautges }
2871d72bce8STim Tautges 
2881d72bce8STim Tautges 
2891d72bce8STim Tautges #undef __FUNCT__
2901d72bce8STim Tautges #define __FUNCT__ "DMMoabGetLocalToGlobalTag"
2911d72bce8STim Tautges PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag)
2921d72bce8STim Tautges {
2931d72bce8STim Tautges   PetscFunctionBegin;
2941d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2951d72bce8STim Tautges   *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag;
2961d72bce8STim Tautges   PetscFunctionReturn(0);
2971d72bce8STim Tautges }
2981d72bce8STim Tautges 
2991d72bce8STim Tautges 
3001d72bce8STim Tautges #undef __FUNCT__
3011d72bce8STim Tautges #define __FUNCT__ "DMMoabSetBlockSize"
3021d72bce8STim Tautges PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs)
3031d72bce8STim Tautges {
3041d72bce8STim Tautges   PetscFunctionBegin;
3051d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3061d72bce8STim Tautges   ((DM_Moab*)dm->data)->bs = bs;
3071d72bce8STim Tautges   PetscFunctionReturn(0);
3081d72bce8STim Tautges }
3091d72bce8STim Tautges 
3101d72bce8STim Tautges 
3111d72bce8STim Tautges #undef __FUNCT__
3121d72bce8STim Tautges #define __FUNCT__ "DMMoabGetBlockSize"
3131d72bce8STim Tautges PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs)
3141d72bce8STim Tautges {
3151d72bce8STim Tautges   PetscFunctionBegin;
3161d72bce8STim Tautges   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3171d72bce8STim Tautges   *bs = ((DM_Moab*)dm->data)->bs;
3181d72bce8STim Tautges   PetscFunctionReturn(0);
3191d72bce8STim Tautges }
3201d72bce8STim Tautges 
3211d72bce8STim Tautges 
3221d72bce8STim Tautges #undef __FUNCT__
3231d72bce8STim Tautges #define __FUNCT__ "DMMoabCreateVectorFromTag"
3241d72bce8STim Tautges PetscErrorCode DMMoabCreateVectorFromTag(DM dm,moab::Tag tag,moab::Range range,PetscBool serial, PetscBool destroy_tag,Vec *X)
3251d72bce8STim Tautges {
3261d72bce8STim Tautges   PetscErrorCode ierr;
3271d72bce8STim Tautges   PetscFunctionBegin;
3281d72bce8STim Tautges   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
3291d72bce8STim Tautges   ierr = VecMoabCreateFromTag(dmmoab->mbint, dmmoab->pcomm, tag, dmmoab->ltog_tag,range, serial, destroy_tag, X);CHKERRQ(ierr);
3301d72bce8STim Tautges   PetscFunctionReturn(0);
3311d72bce8STim Tautges }
3321d72bce8STim Tautges 
3331d72bce8STim Tautges #undef __FUNCT__
3341d72bce8STim Tautges #define __FUNCT__ "DMMoabCreateVector"
3351d72bce8STim Tautges PetscErrorCode DMMoabCreateVector(DM dm,PetscInt tag_size,moab::Range range,PetscBool serial,PetscBool destroy_tag,Vec *vec)
3361d72bce8STim Tautges {
3371d72bce8STim Tautges   PetscErrorCode ierr;
3381d72bce8STim Tautges   PetscFunctionBegin;
3391d72bce8STim Tautges   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
3401d72bce8STim Tautges   ierr = VecMoabCreate(dmmoab->mbint, dmmoab->pcomm, tag_size, dmmoab->ltog_tag, range, serial, destroy_tag, vec);CHKERRQ(ierr);
3411d72bce8STim Tautges   PetscFunctionReturn(0);
3421d72bce8STim Tautges }
3431d72bce8STim Tautges 
3441d72bce8STim Tautges 
3451d72bce8STim Tautges #undef __FUNCT__
3461d72bce8STim Tautges #define __FUNCT__ "DMMoabGetVecTag"
3471d72bce8STim Tautges PetscErrorCode DMMoabGetVecTag(Vec vec,moab::Tag *tag)
3481d72bce8STim Tautges {
3491d72bce8STim Tautges   PetscErrorCode ierr;
3501d72bce8STim Tautges   PetscFunctionBegin;
3511d72bce8STim Tautges   ierr = VecMoabGetTag(vec,tag);CHKERRQ(ierr);
3521d72bce8STim Tautges   PetscFunctionReturn(0);
3531d72bce8STim Tautges }
3541d72bce8STim Tautges 
3551d72bce8STim Tautges 
3561d72bce8STim Tautges #undef __FUNCT__
3571d72bce8STim Tautges #define __FUNCT__ "DMMoabGetVecRange"
3581d72bce8STim Tautges PetscErrorCode DMMoabGetVecRange(Vec vec,moab::Range *range)
3591d72bce8STim Tautges {
3601d72bce8STim Tautges   PetscErrorCode ierr;
3611d72bce8STim Tautges   PetscFunctionBegin;
362*fd349b41STim Tautges   ierr = VecMoabGetRange(vec,range);CHKERRQ(ierr);
3631d72bce8STim Tautges   PetscFunctionReturn(0);
3641d72bce8STim Tautges }
3651d72bce8STim Tautges 
3661d72bce8STim Tautges 
367