1*1d72bce8STim Tautges #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/ 2*1d72bce8STim Tautges #include "../src/vec/vec/impls/moab/vecmoabimpl.h" /*I "petscvec.h" I*/ 3*1d72bce8STim Tautges 4*1d72bce8STim Tautges #include <petscdmmoab.h> 5*1d72bce8STim Tautges #include "MBTagConventions.hpp" 6*1d72bce8STim Tautges 7*1d72bce8STim Tautges typedef struct { 8*1d72bce8STim Tautges PetscInt bs; /* Number of degrees of freedom on each entity, aka tag size in moab */ 9*1d72bce8STim Tautges PetscBool icreatedinstance; /* true if DM created moab instance internally, will destroy instance in DMDestroy */ 10*1d72bce8STim Tautges moab::ParallelComm *pcomm; 11*1d72bce8STim Tautges moab::Interface *mbint; 12*1d72bce8STim Tautges moab::Tag ltog_tag; /* moab supports "global id" tags, which are usually local to global numbering */ 13*1d72bce8STim Tautges moab::Range range; 14*1d72bce8STim Tautges } DM_Moab; 15*1d72bce8STim Tautges 16*1d72bce8STim Tautges #undef __FUNCT__ 17*1d72bce8STim Tautges #define __FUNCT__ "DMMoabCreate" 18*1d72bce8STim Tautges /*@ 19*1d72bce8STim Tautges DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance 20*1d72bce8STim Tautges 21*1d72bce8STim Tautges Collective on MPI_Comm 22*1d72bce8STim Tautges 23*1d72bce8STim Tautges Input Parameter: 24*1d72bce8STim Tautges . comm - The communicator for the DMMoab object 25*1d72bce8STim Tautges 26*1d72bce8STim Tautges Output Parameter: 27*1d72bce8STim Tautges . moab - The DMMoab object 28*1d72bce8STim Tautges 29*1d72bce8STim Tautges Level: beginner 30*1d72bce8STim Tautges 31*1d72bce8STim Tautges .keywords: DMMoab, create 32*1d72bce8STim Tautges @*/ 33*1d72bce8STim Tautges PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *moab) 34*1d72bce8STim Tautges { 35*1d72bce8STim Tautges PetscErrorCode ierr; 36*1d72bce8STim Tautges 37*1d72bce8STim Tautges PetscFunctionBegin; 38*1d72bce8STim Tautges PetscValidPointer(moab,2); 39*1d72bce8STim Tautges ierr = DMCreate(comm, moab);CHKERRQ(ierr); 40*1d72bce8STim Tautges ierr = DMSetType(*moab, DMMOAB);CHKERRQ(ierr); 41*1d72bce8STim Tautges PetscFunctionReturn(0); 42*1d72bce8STim Tautges } 43*1d72bce8STim Tautges 44*1d72bce8STim Tautges EXTERN_C_BEGIN 45*1d72bce8STim Tautges #undef __FUNCT__ 46*1d72bce8STim Tautges #define __FUNCT__ "DMCreate_Moab" 47*1d72bce8STim Tautges PetscErrorCode DMCreate_Moab(DM dm) 48*1d72bce8STim Tautges { 49*1d72bce8STim Tautges DM_Moab *moab; 50*1d72bce8STim Tautges PetscErrorCode ierr; 51*1d72bce8STim Tautges 52*1d72bce8STim Tautges PetscFunctionBegin; 53*1d72bce8STim Tautges PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 54*1d72bce8STim Tautges ierr = PetscNewLog(dm, DM_Moab, &moab);CHKERRQ(ierr); 55*1d72bce8STim Tautges dm->data = moab; 56*1d72bce8STim Tautges 57*1d72bce8STim Tautges PetscFunctionReturn(0); 58*1d72bce8STim Tautges } 59*1d72bce8STim Tautges EXTERN_C_END 60*1d72bce8STim Tautges 61*1d72bce8STim Tautges /*@ 62*1d72bce8STim Tautges DMMoabCreateFromInstance - Creates a DMMoab object from an instance and (optionally) other data 63*1d72bce8STim Tautges 64*1d72bce8STim Tautges Collective on MPI_Comm 65*1d72bce8STim Tautges 66*1d72bce8STim Tautges Input Parameter: 67*1d72bce8STim Tautges . comm - The communicator for the DMMoab object 68*1d72bce8STim Tautges . moab - The MOAB Instance 69*1d72bce8STim Tautges . pcomm - A ParallelComm; if NULL, creates one internally for the whole communicator 70*1d72bce8STim Tautges . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag 71*1d72bce8STim Tautges . range - If non-NULL, contains range of entities to which DOFs will be assigned 72*1d72bce8STim Tautges 73*1d72bce8STim Tautges Output Parameter: 74*1d72bce8STim Tautges . moab - The DMMoab object 75*1d72bce8STim Tautges 76*1d72bce8STim Tautges Level: beginner 77*1d72bce8STim Tautges 78*1d72bce8STim Tautges .keywords: DMMoab, create 79*1d72bce8STim Tautges @*/ 80*1d72bce8STim Tautges #undef __FUNCT__ 81*1d72bce8STim Tautges #define __FUNCT__ "DMMoabCreateFromInstance" 82*1d72bce8STim Tautges PetscErrorCode DMMoabCreateFromInstance(MPI_Comm comm, moab::Interface *iface, moab::ParallelComm *pcomm, moab::Tag ltog_tag, moab::Range *range, DM *moab) 83*1d72bce8STim Tautges { 84*1d72bce8STim Tautges PetscErrorCode ierr; 85*1d72bce8STim Tautges 86*1d72bce8STim Tautges PetscFunctionBegin; 87*1d72bce8STim Tautges PetscValidPointer(moab,2); 88*1d72bce8STim Tautges ierr = DMCreate(comm, moab);CHKERRQ(ierr); 89*1d72bce8STim Tautges ierr = DMSetType(*moab, DMMOAB);CHKERRQ(ierr); 90*1d72bce8STim Tautges ierr = DMInitialize_Moab(*moab);CHKERRQ(ierr); 91*1d72bce8STim Tautges ierr = DMMoabSetInterface(*moab, iface);CHKERRQ(ierr); 92*1d72bce8STim Tautges if (!pcomm) pcomm = new moab::ParallelComm(iface, comm); 93*1d72bce8STim Tautges ierr = DMMoabSetParallelComm(*moab, pcomm);CHKERRQ(ierr); 94*1d72bce8STim Tautges if (!ltog_tag) { 95*1d72bce8STim Tautges moab::ErrorCode merr = iface->tag_get_handle(GLOBAL_ID_TAG_NAME, ltog_tag);MBERRNM(merr); 96*1d72bce8STim Tautges } 97*1d72bce8STim Tautges ierr = DMMoabSetLocalToGlobalTag(*moab, ltog_tag);CHKERRQ(ierr); 98*1d72bce8STim Tautges if (range) { 99*1d72bce8STim Tautges ierr = DMMoabSetRange(*moab, *range);CHKERRQ(ierr); 100*1d72bce8STim Tautges } 101*1d72bce8STim Tautges PetscFunctionReturn(0); 102*1d72bce8STim Tautges } 103*1d72bce8STim Tautges 104*1d72bce8STim Tautges /*@ 105*1d72bce8STim Tautges DMMoabCreateDMAndInstance - Creates a DMMoab object and a MOAB instance 106*1d72bce8STim Tautges 107*1d72bce8STim Tautges Collective on MPI_Comm 108*1d72bce8STim Tautges 109*1d72bce8STim Tautges Input Parameter: 110*1d72bce8STim Tautges . comm - The communicator for the DMMoab object 111*1d72bce8STim Tautges 112*1d72bce8STim Tautges Output Parameter: 113*1d72bce8STim Tautges . moab - The DMMoab object 114*1d72bce8STim Tautges 115*1d72bce8STim Tautges Level: beginner 116*1d72bce8STim Tautges 117*1d72bce8STim Tautges .keywords: DMMoab, create 118*1d72bce8STim Tautges @*/ 119*1d72bce8STim Tautges #undef __FUNCT__ 120*1d72bce8STim Tautges #define __FUNCT__ "DMMoabCreateDMAndInstance" 121*1d72bce8STim Tautges PetscErrorCode DMMoabCreateDMAndInstance(MPI_Comm comm, DM *dm) 122*1d72bce8STim Tautges { 123*1d72bce8STim Tautges PetscErrorCode ierr; 124*1d72bce8STim Tautges 125*1d72bce8STim Tautges PetscFunctionBegin; 126*1d72bce8STim Tautges PetscValidPointer(dm,2); 127*1d72bce8STim Tautges PetscInt rank, nprocs; 128*1d72bce8STim Tautges MPI_Comm_rank(comm, &rank); 129*1d72bce8STim Tautges MPI_Comm_size(comm, &nprocs); 130*1d72bce8STim Tautges moab::Interface *iface = new moab::Core(); 131*1d72bce8STim Tautges moab::ParallelComm *pcomm = new moab::ParallelComm(iface, comm); 132*1d72bce8STim Tautges ierr = DMMoabCreateFromInstance(comm, iface, pcomm, 0, NULL, dm);CHKERRQ(ierr); 133*1d72bce8STim Tautges ((DM_Moab*)(*dm)->data)->icreatedinstance = PETSC_TRUE; 134*1d72bce8STim Tautges PetscFunctionReturn(0); 135*1d72bce8STim Tautges } 136*1d72bce8STim Tautges 137*1d72bce8STim Tautges #undef __FUNCT__ 138*1d72bce8STim Tautges #define __FUNCT__ "DMInitialize_Moab" 139*1d72bce8STim Tautges PetscErrorCode DMInitialize_Moab(DM dm) 140*1d72bce8STim Tautges { 141*1d72bce8STim Tautges PetscFunctionBegin; 142*1d72bce8STim Tautges 143*1d72bce8STim Tautges // Create the DM_Moab and set dm->data 144*1d72bce8STim Tautges DM_Moab *dmmoab = new DM_Moab; 145*1d72bce8STim Tautges dmmoab->bs = 0; 146*1d72bce8STim Tautges dmmoab->pcomm = NULL; 147*1d72bce8STim Tautges dmmoab->mbint = NULL; 148*1d72bce8STim Tautges dmmoab->ltog_tag = (moab::Tag)0; 149*1d72bce8STim Tautges dmmoab->icreatedinstance = PETSC_FALSE; 150*1d72bce8STim Tautges dm->data = dmmoab; 151*1d72bce8STim Tautges 152*1d72bce8STim Tautges // initialize various functions 153*1d72bce8STim Tautges dm->ops->createglobalvector = DMCreateGlobalVector_Moab; 154*1d72bce8STim Tautges dm->ops->createlocalvector = DMCreateLocalVector_Moab; 155*1d72bce8STim Tautges dm->ops->destroy = DMDestroy_Moab; 156*1d72bce8STim Tautges PetscFunctionReturn(0); 157*1d72bce8STim Tautges } 158*1d72bce8STim Tautges 159*1d72bce8STim Tautges #undef __FUNCT__ 160*1d72bce8STim Tautges #define __FUNCT__ "DMMoabSetParallelComm" 161*1d72bce8STim Tautges PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm) 162*1d72bce8STim Tautges { 163*1d72bce8STim Tautges PetscFunctionBegin; 164*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 165*1d72bce8STim Tautges ((DM_Moab*)dm->data)->pcomm = pcomm; 166*1d72bce8STim Tautges ((DM_Moab*)dm->data)->mbint = pcomm->get_moab(); 167*1d72bce8STim Tautges PetscFunctionReturn(0); 168*1d72bce8STim Tautges } 169*1d72bce8STim Tautges 170*1d72bce8STim Tautges 171*1d72bce8STim Tautges #undef __FUNCT__ 172*1d72bce8STim Tautges #define __FUNCT__ "DMMoabGetParallelComm" 173*1d72bce8STim Tautges PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm) 174*1d72bce8STim Tautges { 175*1d72bce8STim Tautges PetscFunctionBegin; 176*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 177*1d72bce8STim Tautges *pcomm = ((DM_Moab*)dm->data)->pcomm; 178*1d72bce8STim Tautges PetscFunctionReturn(0); 179*1d72bce8STim Tautges } 180*1d72bce8STim Tautges 181*1d72bce8STim Tautges 182*1d72bce8STim Tautges #undef __FUNCT__ 183*1d72bce8STim Tautges #define __FUNCT__ "DMMoabSetInterface" 184*1d72bce8STim Tautges PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *iface) 185*1d72bce8STim Tautges { 186*1d72bce8STim Tautges PetscFunctionBegin; 187*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 188*1d72bce8STim Tautges ((DM_Moab*)dm->data)->pcomm = NULL; 189*1d72bce8STim Tautges ((DM_Moab*)dm->data)->mbint = iface; 190*1d72bce8STim Tautges PetscFunctionReturn(0); 191*1d72bce8STim Tautges } 192*1d72bce8STim Tautges 193*1d72bce8STim Tautges 194*1d72bce8STim Tautges #undef __FUNCT__ 195*1d72bce8STim Tautges #define __FUNCT__ "DMMoabGetInterface" 196*1d72bce8STim Tautges PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbint) 197*1d72bce8STim Tautges { 198*1d72bce8STim Tautges PetscFunctionBegin; 199*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 200*1d72bce8STim Tautges *mbint = ((DM_Moab*)dm->data)->mbint; 201*1d72bce8STim Tautges PetscFunctionReturn(0); 202*1d72bce8STim Tautges } 203*1d72bce8STim Tautges 204*1d72bce8STim Tautges 205*1d72bce8STim Tautges #undef __FUNCT__ 206*1d72bce8STim Tautges #define __FUNCT__ "DMMoabSetRange" 207*1d72bce8STim Tautges PetscErrorCode DMMoabSetRange(DM dm,moab::Range range) 208*1d72bce8STim Tautges { 209*1d72bce8STim Tautges PetscFunctionBegin; 210*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 211*1d72bce8STim Tautges ((DM_Moab*)dm->data)->range = range; 212*1d72bce8STim Tautges PetscFunctionReturn(0); 213*1d72bce8STim Tautges } 214*1d72bce8STim Tautges 215*1d72bce8STim Tautges 216*1d72bce8STim Tautges #undef __FUNCT__ 217*1d72bce8STim Tautges #define __FUNCT__ "DMMoabGetRange" 218*1d72bce8STim Tautges PetscErrorCode DMMoabGetRange(DM dm,moab::Range *range) 219*1d72bce8STim Tautges { 220*1d72bce8STim Tautges PetscFunctionBegin; 221*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 222*1d72bce8STim Tautges *range = ((DM_Moab*)dm->data)->range; 223*1d72bce8STim Tautges PetscFunctionReturn(0); 224*1d72bce8STim Tautges } 225*1d72bce8STim Tautges 226*1d72bce8STim Tautges #undef __FUNCT__ 227*1d72bce8STim Tautges #define __FUNCT__ "DMMoabSetLocalToGlobalTag" 228*1d72bce8STim Tautges PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag) 229*1d72bce8STim Tautges { 230*1d72bce8STim Tautges PetscFunctionBegin; 231*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 232*1d72bce8STim Tautges ((DM_Moab*)dm->data)->ltog_tag = ltogtag; 233*1d72bce8STim Tautges PetscFunctionReturn(0); 234*1d72bce8STim Tautges } 235*1d72bce8STim Tautges 236*1d72bce8STim Tautges 237*1d72bce8STim Tautges #undef __FUNCT__ 238*1d72bce8STim Tautges #define __FUNCT__ "DMMoabGetLocalToGlobalTag" 239*1d72bce8STim Tautges PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag) 240*1d72bce8STim Tautges { 241*1d72bce8STim Tautges PetscFunctionBegin; 242*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 243*1d72bce8STim Tautges *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag; 244*1d72bce8STim Tautges PetscFunctionReturn(0); 245*1d72bce8STim Tautges } 246*1d72bce8STim Tautges 247*1d72bce8STim Tautges 248*1d72bce8STim Tautges #undef __FUNCT__ 249*1d72bce8STim Tautges #define __FUNCT__ "DMMoabSetBlockSize" 250*1d72bce8STim Tautges PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs) 251*1d72bce8STim Tautges { 252*1d72bce8STim Tautges PetscFunctionBegin; 253*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 254*1d72bce8STim Tautges ((DM_Moab*)dm->data)->bs = bs; 255*1d72bce8STim Tautges PetscFunctionReturn(0); 256*1d72bce8STim Tautges } 257*1d72bce8STim Tautges 258*1d72bce8STim Tautges 259*1d72bce8STim Tautges #undef __FUNCT__ 260*1d72bce8STim Tautges #define __FUNCT__ "DMMoabGetBlockSize" 261*1d72bce8STim Tautges PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs) 262*1d72bce8STim Tautges { 263*1d72bce8STim Tautges PetscFunctionBegin; 264*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 265*1d72bce8STim Tautges *bs = ((DM_Moab*)dm->data)->bs; 266*1d72bce8STim Tautges PetscFunctionReturn(0); 267*1d72bce8STim Tautges } 268*1d72bce8STim Tautges 269*1d72bce8STim Tautges 270*1d72bce8STim Tautges #undef __FUNCT__ 271*1d72bce8STim Tautges #define __FUNCT__ "DMCreateGlobalVector_Moab" 272*1d72bce8STim Tautges PetscErrorCode DMCreateGlobalVector_Moab(DM dm,Vec *gvec) 273*1d72bce8STim Tautges { 274*1d72bce8STim Tautges PetscErrorCode ierr; 275*1d72bce8STim Tautges DM_Moab *dmmoab = (DM_Moab*)dm->data; 276*1d72bce8STim Tautges 277*1d72bce8STim Tautges PetscFunctionBegin; 278*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 279*1d72bce8STim Tautges PetscValidPointer(gvec,2); 280*1d72bce8STim Tautges PetscInt block_size = ((DM_Moab*)dm->data)->bs; 281*1d72bce8STim Tautges ierr = VecMoabCreate(dmmoab->mbint,dmmoab->pcomm,block_size,dmmoab->ltog_tag,dmmoab->range,PETSC_FALSE,PETSC_TRUE,gvec);CHKERRQ(ierr); 282*1d72bce8STim Tautges PetscFunctionReturn(0); 283*1d72bce8STim Tautges } 284*1d72bce8STim Tautges 285*1d72bce8STim Tautges 286*1d72bce8STim Tautges #undef __FUNCT__ 287*1d72bce8STim Tautges #define __FUNCT__ "DMCreateLocalVector_Moab" 288*1d72bce8STim Tautges PetscErrorCode DMCreateLocalVector_Moab(DM dm,Vec *gvec) 289*1d72bce8STim Tautges { 290*1d72bce8STim Tautges PetscErrorCode ierr; 291*1d72bce8STim Tautges DM_Moab *dmmoab = (DM_Moab*)dm->data; 292*1d72bce8STim Tautges 293*1d72bce8STim Tautges PetscFunctionBegin; 294*1d72bce8STim Tautges PetscInt bs = 1; 295*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 296*1d72bce8STim Tautges PetscValidPointer(gvec,2); 297*1d72bce8STim Tautges ierr = VecMoabCreate(dmmoab->mbint,dmmoab->pcomm,bs,dmmoab->ltog_tag,dmmoab->range,PETSC_TRUE,PETSC_TRUE,gvec);CHKERRQ(ierr); 298*1d72bce8STim Tautges PetscFunctionReturn(0); 299*1d72bce8STim Tautges } 300*1d72bce8STim Tautges 301*1d72bce8STim Tautges #undef __FUNCT__ 302*1d72bce8STim Tautges #define __FUNCT__ "DMMoabCreateVectorFromTag" 303*1d72bce8STim Tautges PetscErrorCode DMMoabCreateVectorFromTag(DM dm,moab::Tag tag,moab::Range range,PetscBool serial, PetscBool destroy_tag,Vec *X) 304*1d72bce8STim Tautges { 305*1d72bce8STim Tautges PetscErrorCode ierr; 306*1d72bce8STim Tautges PetscFunctionBegin; 307*1d72bce8STim Tautges DM_Moab *dmmoab = (DM_Moab*)dm->data; 308*1d72bce8STim Tautges ierr = VecMoabCreateFromTag(dmmoab->mbint, dmmoab->pcomm, tag, dmmoab->ltog_tag,range, serial, destroy_tag, X);CHKERRQ(ierr); 309*1d72bce8STim Tautges PetscFunctionReturn(0); 310*1d72bce8STim Tautges } 311*1d72bce8STim Tautges 312*1d72bce8STim Tautges #undef __FUNCT__ 313*1d72bce8STim Tautges #define __FUNCT__ "DMMoabCreateVector" 314*1d72bce8STim Tautges PetscErrorCode DMMoabCreateVector(DM dm,PetscInt tag_size,moab::Range range,PetscBool serial,PetscBool destroy_tag,Vec *vec) 315*1d72bce8STim Tautges { 316*1d72bce8STim Tautges PetscErrorCode ierr; 317*1d72bce8STim Tautges PetscFunctionBegin; 318*1d72bce8STim Tautges DM_Moab *dmmoab = (DM_Moab*)dm->data; 319*1d72bce8STim Tautges ierr = VecMoabCreate(dmmoab->mbint, dmmoab->pcomm, tag_size, dmmoab->ltog_tag, range, serial, destroy_tag, vec);CHKERRQ(ierr); 320*1d72bce8STim Tautges PetscFunctionReturn(0); 321*1d72bce8STim Tautges } 322*1d72bce8STim Tautges 323*1d72bce8STim Tautges 324*1d72bce8STim Tautges #undef __FUNCT__ 325*1d72bce8STim Tautges #define __FUNCT__ "DMMoabGetVecTag" 326*1d72bce8STim Tautges PetscErrorCode DMMoabGetVecTag(Vec vec,moab::Tag *tag) 327*1d72bce8STim Tautges { 328*1d72bce8STim Tautges PetscErrorCode ierr; 329*1d72bce8STim Tautges PetscFunctionBegin; 330*1d72bce8STim Tautges ierr = VecMoabGetTag(vec,tag);CHKERRQ(ierr); 331*1d72bce8STim Tautges PetscFunctionReturn(0); 332*1d72bce8STim Tautges } 333*1d72bce8STim Tautges 334*1d72bce8STim Tautges 335*1d72bce8STim Tautges #undef __FUNCT__ 336*1d72bce8STim Tautges #define __FUNCT__ "DMMoabGetVecRange" 337*1d72bce8STim Tautges PetscErrorCode DMMoabGetVecRange(Vec vec,moab::Range *range) 338*1d72bce8STim Tautges { 339*1d72bce8STim Tautges PetscErrorCode ierr; 340*1d72bce8STim Tautges PetscFunctionBegin; 341*1d72bce8STim Tautges ierr = VecMoabGetRange(vec,range); 342*1d72bce8STim Tautges PetscFunctionReturn(0); 343*1d72bce8STim Tautges } 344*1d72bce8STim Tautges 345*1d72bce8STim Tautges 346*1d72bce8STim Tautges #undef __FUNCT__ 347*1d72bce8STim Tautges #define __FUNCT__ "DMDestroy_Moab" 348*1d72bce8STim Tautges PetscErrorCode DMDestroy_Moab(DM dm) 349*1d72bce8STim Tautges { 350*1d72bce8STim Tautges PetscFunctionBegin; 351*1d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 352*1d72bce8STim Tautges 353*1d72bce8STim Tautges // Delete the DM_Moab: 354*1d72bce8STim Tautges if(dm->data) { 355*1d72bce8STim Tautges if (((DM_Moab*)dm->data)->icreatedinstance) { 356*1d72bce8STim Tautges delete ((DM_Moab*)dm->data)->mbint; 357*1d72bce8STim Tautges ((DM_Moab*)dm->data)->mbint = NULL; 358*1d72bce8STim Tautges ((DM_Moab*)dm->data)->pcomm = NULL; 359*1d72bce8STim Tautges } 360*1d72bce8STim Tautges delete (DM_Moab*)dm->data; 361*1d72bce8STim Tautges dm->data = NULL; 362*1d72bce8STim Tautges } 363*1d72bce8STim Tautges PetscFunctionReturn(0); 364*1d72bce8STim Tautges } 365*1d72bce8STim Tautges 366*1d72bce8STim Tautges 367