1032b8ab6SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I "petscdm.h" I*/ 2aa768e4cSTim Tautges #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/ 31d72bce8STim Tautges 41d72bce8STim Tautges #include <petscdmmoab.h> 588face26SJed Brown #include <MBTagConventions.hpp> 61cec0304SVijay Mahadevan #include <moab/Skinner.hpp> 7032b8ab6SVijay Mahadevan 8853cdec3SJed Brown #undef __FUNCT__ 9853cdec3SJed Brown #define __FUNCT__ "DMDestroy_Moab" 10853cdec3SJed Brown PetscErrorCode DMDestroy_Moab(DM dm) 11853cdec3SJed Brown { 12853cdec3SJed Brown PetscErrorCode ierr; 13032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 14212ad6d1SVijay Mahadevan PetscSection section; 15853cdec3SJed Brown 16853cdec3SJed Brown PetscFunctionBegin; 17853cdec3SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 18212ad6d1SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 19212ad6d1SVijay Mahadevan ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 20032b8ab6SVijay Mahadevan if (dmmoab->icreatedinstance) { 21032b8ab6SVijay Mahadevan delete dmmoab->mbiface; 22853cdec3SJed Brown } 23032b8ab6SVijay Mahadevan dmmoab->mbiface = NULL; 24032b8ab6SVijay Mahadevan dmmoab->pcomm = NULL; 25032b8ab6SVijay Mahadevan delete dmmoab->vlocal; 26032b8ab6SVijay Mahadevan delete dmmoab->vowned; 27032b8ab6SVijay Mahadevan delete dmmoab->vghost; 28032b8ab6SVijay Mahadevan delete dmmoab->elocal; 29032b8ab6SVijay Mahadevan delete dmmoab->eghost; 308d8d51c8SVijay Mahadevan 318d8d51c8SVijay Mahadevan ierr = PetscFree(dmmoab->isbndyvtx);CHKERRQ(ierr); 328d8d51c8SVijay Mahadevan ierr = PetscFree(dmmoab->isbndyfaces);CHKERRQ(ierr); 338d8d51c8SVijay Mahadevan ierr = PetscFree(dmmoab->isbndyelems);CHKERRQ(ierr); 34fc418013SVijay Mahadevan ierr = PetscFree(dmmoab->gsindices);CHKERRQ(ierr); 35032b8ab6SVijay Mahadevan ierr = VecScatterDestroy(&dmmoab->ltog_sendrecv);CHKERRQ(ierr); 36032b8ab6SVijay Mahadevan ierr = ISLocalToGlobalMappingDestroy(&dmmoab->ltog_map);CHKERRQ(ierr); 37853cdec3SJed Brown ierr = PetscFree(dm->data);CHKERRQ(ierr); 38853cdec3SJed Brown PetscFunctionReturn(0); 39853cdec3SJed Brown } 40853cdec3SJed Brown 41aa768e4cSTim Tautges #undef __FUNCT__ 42032b8ab6SVijay Mahadevan #define __FUNCT__ "DMSetUp_Moab" 43032b8ab6SVijay Mahadevan PetscErrorCode DMSetUp_Moab(DM dm) 44032b8ab6SVijay Mahadevan { 45032b8ab6SVijay Mahadevan PetscErrorCode ierr; 46032b8ab6SVijay Mahadevan moab::ErrorCode merr; 47032b8ab6SVijay Mahadevan Vec local, global; 48032b8ab6SVijay Mahadevan IS from; 49032b8ab6SVijay Mahadevan moab::Range::iterator iter; 50fc418013SVijay Mahadevan PetscInt i,j,bs,gsiz,lsiz; 51032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 52eb9d2429SVijay Mahadevan PetscInt totsize; 531cec0304SVijay Mahadevan PetscSection section; 54eb9d2429SVijay Mahadevan PetscInt gmin,lmin,lmax; 55032b8ab6SVijay Mahadevan 5669263071SVijay Mahadevan moab::Range adj; 5769263071SVijay Mahadevan 58032b8ab6SVijay Mahadevan PetscFunctionBegin; 59032b8ab6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 60032b8ab6SVijay Mahadevan /* Get the local and shared vertices and cache it */ 61032b8ab6SVijay Mahadevan if (dmmoab->mbiface == PETSC_NULL || dmmoab->pcomm == PETSC_NULL) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ORDER, "Set the MOAB Interface and ParallelComm objects before calling SetUp."); 62032b8ab6SVijay Mahadevan 631cec0304SVijay Mahadevan /* Get the entities recursively in the current part of the mesh, if user did not set the local vertices explicitly */ 64032b8ab6SVijay Mahadevan if (dmmoab->vlocal->empty()) { 651cec0304SVijay Mahadevan merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,*dmmoab->vlocal,true);MBERRNM(merr); 66032b8ab6SVijay Mahadevan 67032b8ab6SVijay Mahadevan /* filter based on parallel status */ 68fc418013SVijay Mahadevan merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,dmmoab->vowned);MBERRNM(merr); 69032b8ab6SVijay Mahadevan *dmmoab->vghost = moab::subtract(*dmmoab->vlocal, *dmmoab->vowned); 70032b8ab6SVijay Mahadevan 71032b8ab6SVijay Mahadevan dmmoab->nloc = dmmoab->vowned->size(); 72032b8ab6SVijay Mahadevan dmmoab->nghost = dmmoab->vghost->size(); 73032b8ab6SVijay Mahadevan ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 74fc418013SVijay Mahadevan 75e23c60ebSVijay Mahadevan #if 0 76fc418013SVijay Mahadevan if(dmmoab->pcomm->rank() || dmmoab->pcomm->size()==1) { 77fc418013SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "Vertices: global: %D, local: %D", dmmoab->n, dmmoab->nloc+dmmoab->nghost); 78fc418013SVijay Mahadevan dmmoab->vlocal->print(0); 79fc418013SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "Vertices: owned: %D", dmmoab->nloc); 80fc418013SVijay Mahadevan dmmoab->vowned->print(0); 81fc418013SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "Vertices: ghost: %D", dmmoab->nghost); 82fc418013SVijay Mahadevan dmmoab->vghost->print(0); 83fc418013SVijay Mahadevan } 84fc418013SVijay Mahadevan #endif 85032b8ab6SVijay Mahadevan } 86032b8ab6SVijay Mahadevan 87032b8ab6SVijay Mahadevan /* get the information about the local elements in the mesh */ 88032b8ab6SVijay Mahadevan { 89032b8ab6SVijay Mahadevan dmmoab->eghost->clear(); 90fc418013SVijay Mahadevan 91fc418013SVijay Mahadevan /* first decipher the leading dimension */ 92fc418013SVijay Mahadevan for (i=3;i>0;i--) { 93fc418013SVijay Mahadevan dmmoab->elocal->clear(); 94fc418013SVijay Mahadevan merr = dmmoab->mbiface->get_entities_by_dimension(dmmoab->fileset, i, *dmmoab->elocal, true);CHKERRQ(merr); 95fc418013SVijay Mahadevan 96fc418013SVijay Mahadevan /* store the current mesh dimension */ 97fc418013SVijay Mahadevan if (dmmoab->elocal->size()) { 98fc418013SVijay Mahadevan dmmoab->dim=i; 99fc418013SVijay Mahadevan break; 100fc418013SVijay Mahadevan } 101fc418013SVijay Mahadevan } 102fc418013SVijay Mahadevan 103032b8ab6SVijay Mahadevan *dmmoab->eghost = *dmmoab->elocal; 104032b8ab6SVijay Mahadevan merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 105032b8ab6SVijay Mahadevan *dmmoab->eghost = moab::subtract(*dmmoab->eghost, *dmmoab->elocal); 106032b8ab6SVijay Mahadevan 107032b8ab6SVijay Mahadevan dmmoab->neleloc = dmmoab->elocal->size(); 108032b8ab6SVijay Mahadevan ierr = MPI_Allreduce(&dmmoab->neleloc, &dmmoab->nele, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 109032b8ab6SVijay Mahadevan } 110032b8ab6SVijay Mahadevan 111032b8ab6SVijay Mahadevan bs = dmmoab->bs; 112032b8ab6SVijay Mahadevan if (!dmmoab->ltog_tag) { 113db66d124SVijay Mahadevan /* Get the global ID tag. The global ID tag is applied to each 114db66d124SVijay Mahadevan vertex. It acts as an global identifier which MOAB uses to 115db66d124SVijay Mahadevan assemble the individual pieces of the mesh */ 116032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr); 117032b8ab6SVijay Mahadevan } 118032b8ab6SVijay Mahadevan 119032b8ab6SVijay Mahadevan totsize=dmmoab->vlocal->size(); 120fc418013SVijay Mahadevan ierr = PetscMalloc(totsize*sizeof(PetscInt), &dmmoab->gsindices);CHKERRQ(ierr); 1211cec0304SVijay Mahadevan { 122032b8ab6SVijay Mahadevan /* first get the local indices */ 123fc418013SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vowned,&dmmoab->gsindices[0]);MBERRNM(merr); 1244a40b570SVijay Mahadevan /* next get the ghosted indices */ 125fc418013SVijay Mahadevan if (dmmoab->nghost) { 126fc418013SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vghost,&dmmoab->gsindices[dmmoab->nloc]);MBERRNM(merr); 1271cec0304SVijay Mahadevan } 1286e40195eSVijay Mahadevan 1296e40195eSVijay Mahadevan /* find out the local and global minima of GLOBAL_ID */ 130eb9d2429SVijay Mahadevan lmin=lmax=dmmoab->gsindices[0]; 13169263071SVijay Mahadevan for (i=0; i<totsize; ++i) { 132eb9d2429SVijay Mahadevan if(lmin>dmmoab->gsindices[i]) lmin=dmmoab->gsindices[i]; 133eb9d2429SVijay Mahadevan if(lmax<dmmoab->gsindices[i]) lmax=dmmoab->gsindices[i]; 134fc418013SVijay Mahadevan } 1356e40195eSVijay Mahadevan 136eb9d2429SVijay Mahadevan ierr = MPI_Allreduce(&lmin, &gmin, 1, MPI_INT, MPI_MIN, ((PetscObject)dm)->comm);CHKERRQ(ierr); 137eb9d2429SVijay Mahadevan PetscInfo3(dm, "GLOBAL_ID: Local minima - %D, Local maxima - %D, Global minima - %D.\n", lmin, lmax, gmin); 1381cec0304SVijay Mahadevan } 139032b8ab6SVijay Mahadevan 1401cec0304SVijay Mahadevan { 1411cec0304SVijay Mahadevan ierr = PetscSectionCreate(((PetscObject)dm)->comm, §ion);CHKERRQ(ierr); 1421cec0304SVijay Mahadevan ierr = PetscSectionSetNumFields(section, dmmoab->nfields);CHKERRQ(ierr); 143eb9d2429SVijay Mahadevan ierr = PetscSectionSetChart(section, lmin, lmax+1);CHKERRQ(ierr); 144fc418013SVijay Mahadevan for (j=0; j<totsize; ++j) { 145fc418013SVijay Mahadevan PetscInt locgid = dmmoab->gsindices[j]; 1461cec0304SVijay Mahadevan for (i=0; i < dmmoab->nfields; ++i) { 1471cec0304SVijay Mahadevan ierr = PetscSectionSetFieldName(section, i, dmmoab->fields[i]);CHKERRQ(ierr); 1481cec0304SVijay Mahadevan if (bs>1) { 149eb9d2429SVijay Mahadevan ierr = PetscSectionSetFieldDof(section, locgid, i, (locgid-gmin)*dmmoab->nfields+i);CHKERRQ(ierr); 150eb9d2429SVijay Mahadevan ierr = PetscSectionSetFieldOffset(section, locgid, i, (locgid-gmin)*dmmoab->nfields); 1511cec0304SVijay Mahadevan } 1521cec0304SVijay Mahadevan else { 153eb9d2429SVijay Mahadevan ierr = PetscSectionSetFieldDof(section, locgid, i, dmmoab->n*i+locgid-gmin);CHKERRQ(ierr); 15469263071SVijay Mahadevan ierr = PetscSectionSetFieldOffset(section, locgid, i, i*dmmoab->n); 1551cec0304SVijay Mahadevan } 1561cec0304SVijay Mahadevan } 157fc418013SVijay Mahadevan ierr = PetscSectionSetDof(section, locgid, dmmoab->nfields);CHKERRQ(ierr); 1581cec0304SVijay Mahadevan } 1591cec0304SVijay Mahadevan ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 1601cec0304SVijay Mahadevan ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr); 1611cec0304SVijay Mahadevan } 1621cec0304SVijay Mahadevan 1631cec0304SVijay Mahadevan { 164fc418013SVijay Mahadevan for (i=0; i<totsize; ++i) { 165eb9d2429SVijay Mahadevan dmmoab->gsindices[i]-=gmin; /* zero based index needed for IS */ 166fc418013SVijay Mahadevan } 167fc418013SVijay Mahadevan 168032b8ab6SVijay Mahadevan /* Create Global to Local Vector Scatter Context */ 169032b8ab6SVijay Mahadevan ierr = DMCreateGlobalVector_Moab(dm, &global);CHKERRQ(ierr); 170032b8ab6SVijay Mahadevan ierr = DMCreateLocalVector_Moab(dm, &local);CHKERRQ(ierr); 171032b8ab6SVijay Mahadevan 172032b8ab6SVijay Mahadevan /* global to local must retrieve ghost points */ 173fc418013SVijay Mahadevan ierr = ISCreateBlock(((PetscObject)dm)->comm,bs,totsize,&dmmoab->gsindices[0],PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 174032b8ab6SVijay Mahadevan 175db66d124SVijay Mahadevan ierr = VecGetLocalSize(global,&gsiz);CHKERRQ(ierr); 176db66d124SVijay Mahadevan ierr = VecGetLocalSize(local,&lsiz);CHKERRQ(ierr); 177032b8ab6SVijay Mahadevan 178032b8ab6SVijay Mahadevan ierr = VecScatterCreate(local,from,global,from,&dmmoab->ltog_sendrecv);CHKERRQ(ierr); 179032b8ab6SVijay Mahadevan ierr = ISDestroy(&from);CHKERRQ(ierr); 180032b8ab6SVijay Mahadevan ierr = VecDestroy(&local);CHKERRQ(ierr); 181032b8ab6SVijay Mahadevan ierr = VecDestroy(&global);CHKERRQ(ierr); 182032b8ab6SVijay Mahadevan } 183032b8ab6SVijay Mahadevan 1841cec0304SVijay Mahadevan /* skin the boundary and store nodes */ 1851cec0304SVijay Mahadevan { 186*0c8a2322SVijay Mahadevan moab::Range bndyfaces, bndyvtx, bndyelems; 187eb9d2429SVijay Mahadevan /* get the skin vertices of boundary faces for the current partition and then filter 188eb9d2429SVijay Mahadevan the local, boundary faces, vertices and elements alone via PSTATUS flags; 189eb9d2429SVijay Mahadevan this should not give us any ghosted boundary, but if user needs such a functionality 190eb9d2429SVijay Mahadevan it would be easy to add it based on the find_skin query below */ 1911cec0304SVijay Mahadevan moab::Skinner skinner(dmmoab->mbiface); 192eb9d2429SVijay Mahadevan 193eb9d2429SVijay Mahadevan /* get the entities on the skin - only the faces */ 194*0c8a2322SVijay Mahadevan merr = skinner.find_skin(dmmoab->fileset, *dmmoab->elocal, false, bndyfaces, NULL, false, true, false, false);MBERRNM(merr); // 'false' param indicates we want faces back, not vertices 195eb9d2429SVijay Mahadevan 196eb9d2429SVijay Mahadevan /* filter all the non-owned and shared entities out of the list */ 197*0c8a2322SVijay Mahadevan merr = dmmoab->pcomm->filter_pstatus(bndyfaces,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 198*0c8a2322SVijay Mahadevan merr = dmmoab->pcomm->filter_pstatus(bndyfaces,PSTATUS_SHARED,PSTATUS_NOT);MBERRNM(merr); 19969263071SVijay Mahadevan 200eb9d2429SVijay Mahadevan /* get all the nodes via connectivity and the parent elements via adjacency information */ 201*0c8a2322SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity(bndyfaces, bndyvtx, false);MBERRNM(ierr); 202*0c8a2322SVijay Mahadevan merr = dmmoab->mbiface->get_adjacencies(bndyfaces, dmmoab->dim, false, bndyelems, moab::Interface::UNION);MBERRNM(ierr); 203*0c8a2322SVijay Mahadevan PetscInfo3(dm, "Found %D boundary vertices, %D boundary faces and %D boundary elements.\n", bndyvtx.size(), bndyvtx.size(), bndyelems.size()); 2048d8d51c8SVijay Mahadevan 2058d8d51c8SVijay Mahadevan /* cache a bit-vector for easy query */ 206*0c8a2322SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscBool)*((PetscInt)(*bndyvtx.rbegin())+1),&dmmoab->isbndyvtx);CHKERRQ(ierr); 207*0c8a2322SVijay Mahadevan ierr = PetscMemzero(dmmoab->isbndyvtx,sizeof(PetscBool)*((PetscInt)(*bndyvtx.rbegin())+1));CHKERRQ(ierr); 208*0c8a2322SVijay Mahadevan for(moab::Range::iterator iter = bndyvtx.begin(); iter != bndyvtx.end(); iter++) { 209*0c8a2322SVijay Mahadevan dmmoab->isbndyvtx[(PetscInt)*iter]=PETSC_TRUE; 2108d8d51c8SVijay Mahadevan } 2118d8d51c8SVijay Mahadevan 212*0c8a2322SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscBool)*((PetscInt)(*bndyelems.rbegin())+1),&dmmoab->isbndyelems);CHKERRQ(ierr); 213*0c8a2322SVijay Mahadevan ierr = PetscMemzero(dmmoab->isbndyelems,sizeof(PetscBool)*((PetscInt)(*bndyelems.rbegin())+1));CHKERRQ(ierr); 214*0c8a2322SVijay Mahadevan for(moab::Range::iterator iter = bndyelems.begin(); iter != bndyelems.end(); iter++) { 215*0c8a2322SVijay Mahadevan dmmoab->isbndyelems[(PetscInt)*iter]=PETSC_TRUE; 2168d8d51c8SVijay Mahadevan } 2178d8d51c8SVijay Mahadevan 218*0c8a2322SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscBool)*((PetscInt)(*bndyfaces.rbegin())+1),&dmmoab->isbndyfaces);CHKERRQ(ierr); 219*0c8a2322SVijay Mahadevan ierr = PetscMemzero(dmmoab->isbndyfaces,sizeof(PetscBool)*((PetscInt)(*bndyfaces.rbegin())+1));CHKERRQ(ierr); 220*0c8a2322SVijay Mahadevan for(moab::Range::iterator iter = bndyfaces.begin(); iter != bndyfaces.end(); iter++) { 221*0c8a2322SVijay Mahadevan dmmoab->isbndyfaces[(PetscInt)*iter]=PETSC_TRUE; 2228d8d51c8SVijay Mahadevan } 2231cec0304SVijay Mahadevan } 224032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 225032b8ab6SVijay Mahadevan } 226032b8ab6SVijay Mahadevan 2271cec0304SVijay Mahadevan 228032b8ab6SVijay Mahadevan #undef __FUNCT__ 229aa768e4cSTim Tautges #define __FUNCT__ "DMCreate_Moab" 230853cdec3SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Moab(DM dm) 231aa768e4cSTim Tautges { 232aa768e4cSTim Tautges PetscErrorCode ierr; 233aa768e4cSTim Tautges 234aa768e4cSTim Tautges PetscFunctionBegin; 235aa768e4cSTim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 236032b8ab6SVijay Mahadevan ierr = PetscNewLog(dm,&dm->data);CHKERRQ(ierr); 237032b8ab6SVijay Mahadevan 238032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->bs = 1; 239324f1edfSVijay Mahadevan ((DM_Moab*)dm->data)->nfields = 1; 240032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->n = 0; 241032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->nloc = 0; 242032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->nele = 0; 243032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->neleloc = 0; 244032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->nghost = 0; 245032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->ltog_map = PETSC_NULL; 246032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->ltog_sendrecv = PETSC_NULL; 247032b8ab6SVijay Mahadevan 248032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->vlocal = new moab::Range(); 249032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->vowned = new moab::Range(); 250032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->vghost = new moab::Range(); 251032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->elocal = new moab::Range(); 252032b8ab6SVijay Mahadevan ((DM_Moab*)dm->data)->eghost = new moab::Range(); 253aa768e4cSTim Tautges 25497ea90e6SJed Brown dm->ops->createglobalvector = DMCreateGlobalVector_Moab; 25597ea90e6SJed Brown dm->ops->createlocalvector = DMCreateLocalVector_Moab; 256032b8ab6SVijay Mahadevan dm->ops->creatematrix = DMCreateMatrix_Moab; 257032b8ab6SVijay Mahadevan dm->ops->setup = DMSetUp_Moab; 25897ea90e6SJed Brown dm->ops->destroy = DMDestroy_Moab; 259032b8ab6SVijay Mahadevan dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Moab; 260032b8ab6SVijay Mahadevan dm->ops->globaltolocalend = DMGlobalToLocalEnd_Moab; 261032b8ab6SVijay Mahadevan dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Moab; 262032b8ab6SVijay Mahadevan dm->ops->localtoglobalend = DMLocalToGlobalEnd_Moab; 263aa768e4cSTim Tautges PetscFunctionReturn(0); 264aa768e4cSTim Tautges } 265fd349b41STim Tautges 266fd349b41STim Tautges #undef __FUNCT__ 2671d72bce8STim Tautges #define __FUNCT__ "DMMoabCreate" 2681d72bce8STim Tautges /*@ 2691d72bce8STim Tautges DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance 2701d72bce8STim Tautges 2711d72bce8STim Tautges Collective on MPI_Comm 2721d72bce8STim Tautges 2731d72bce8STim Tautges Input Parameter: 2741d72bce8STim Tautges . comm - The communicator for the DMMoab object 2751d72bce8STim Tautges 2761d72bce8STim Tautges Output Parameter: 277032b8ab6SVijay Mahadevan . dmb - The DMMoab object 2781d72bce8STim Tautges 2791d72bce8STim Tautges Level: beginner 2801d72bce8STim Tautges 2811d72bce8STim Tautges .keywords: DMMoab, create 2821d72bce8STim Tautges @*/ 283032b8ab6SVijay Mahadevan PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *dmb) 2841d72bce8STim Tautges { 2851d72bce8STim Tautges PetscErrorCode ierr; 2861d72bce8STim Tautges 2871d72bce8STim Tautges PetscFunctionBegin; 288032b8ab6SVijay Mahadevan PetscValidPointer(dmb,2); 289032b8ab6SVijay Mahadevan ierr = DMCreate(comm, dmb);CHKERRQ(ierr); 290032b8ab6SVijay Mahadevan ierr = DMSetType(*dmb, DMMOAB);CHKERRQ(ierr); 2911d72bce8STim Tautges PetscFunctionReturn(0); 2921d72bce8STim Tautges } 2931d72bce8STim Tautges 2941d72bce8STim Tautges #undef __FUNCT__ 295aa768e4cSTim Tautges #define __FUNCT__ "DMMoabCreateMoab" 2961d72bce8STim Tautges /*@ 297a4d2169cSTim Tautges DMMoabCreate - Creates a DMMoab object, optionally from an instance and other data 2981d72bce8STim Tautges 2991d72bce8STim Tautges Collective on MPI_Comm 3001d72bce8STim Tautges 3011d72bce8STim Tautges Input Parameter: 3021d72bce8STim Tautges . comm - The communicator for the DMMoab object 303032b8ab6SVijay Mahadevan . mbiface - (ptr to) the MOAB Instance; if passed in NULL, MOAB instance is created inside PETSc, and destroyed 304a4d2169cSTim Tautges along with the DMMoab 305a4d2169cSTim Tautges . pcomm - (ptr to) a ParallelComm; if NULL, creates one internally for the whole communicator 3061d72bce8STim Tautges . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag 3071d72bce8STim Tautges . range - If non-NULL, contains range of entities to which DOFs will be assigned 3081d72bce8STim Tautges 3091d72bce8STim Tautges Output Parameter: 310032b8ab6SVijay Mahadevan . dmb - The DMMoab object 3111d72bce8STim Tautges 312032b8ab6SVijay Mahadevan Level: intermediate 3131d72bce8STim Tautges 3141d72bce8STim Tautges .keywords: DMMoab, create 3151d72bce8STim Tautges @*/ 316032b8ab6SVijay Mahadevan PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *dmb) 3171d72bce8STim Tautges { 3181d72bce8STim Tautges PetscErrorCode ierr; 319032b8ab6SVijay Mahadevan moab::ErrorCode merr; 3201cec0304SVijay Mahadevan moab::EntityHandle partnset; 3211cec0304SVijay Mahadevan PetscInt rank, nprocs; 322853cdec3SJed Brown DM_Moab *dmmoab; 3231d72bce8STim Tautges 3241d72bce8STim Tautges PetscFunctionBegin; 325032b8ab6SVijay Mahadevan PetscValidPointer(dmb,6); 326032b8ab6SVijay Mahadevan ierr = DMMoabCreate(comm, dmb);CHKERRQ(ierr); 327032b8ab6SVijay Mahadevan dmmoab = (DM_Moab*)(*dmb)->data; 328a4d2169cSTim Tautges 329a4d2169cSTim Tautges if (!mbiface) { 33072ff976dSVijay Mahadevan dmmoab->mbiface = new moab::Core(); 3317d89fc02STim Tautges dmmoab->icreatedinstance = PETSC_TRUE; 3321d72bce8STim Tautges } 3331cec0304SVijay Mahadevan else { 3341cec0304SVijay Mahadevan dmmoab->mbiface = mbiface; 3357d89fc02STim Tautges dmmoab->icreatedinstance = PETSC_FALSE; 3361cec0304SVijay Mahadevan } 3371cec0304SVijay Mahadevan 3381cec0304SVijay Mahadevan /* create a fileset to store the hierarchy of entities belonging to current DM */ 3391cec0304SVijay Mahadevan merr = dmmoab->mbiface->create_meshset(moab::MESHSET_ORDERED, dmmoab->fileset);MBERR("Creating file set failed", merr); 3407d89fc02STim Tautges 341a4d2169cSTim Tautges if (!pcomm) { 342032b8ab6SVijay Mahadevan ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 343032b8ab6SVijay Mahadevan ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr); 344032b8ab6SVijay Mahadevan 345db66d124SVijay Mahadevan /* Create root sets for each mesh. Then pass these 346db66d124SVijay Mahadevan to the load_file functions to be populated. */ 347*0c8a2322SVijay Mahadevan merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, partnset);MBERR("Creating partition set failed", merr); 348032b8ab6SVijay Mahadevan 349db66d124SVijay Mahadevan /* Create the parallel communicator object with the partition handle associated with MOAB */ 35072ff976dSVijay Mahadevan dmmoab->pcomm = moab::ParallelComm::get_pcomm(dmmoab->mbiface, partnset, &comm); 35172ff976dSVijay Mahadevan } 35272ff976dSVijay Mahadevan else { 35372ff976dSVijay Mahadevan ierr = DMMoabSetParallelComm(*dmb, pcomm);CHKERRQ(ierr); 354032b8ab6SVijay Mahadevan } 355032b8ab6SVijay Mahadevan 3564973de03SVijay Mahadevan /* do the remaining initializations for DMMoab */ 3574973de03SVijay Mahadevan dmmoab->bs = 1; 358324f1edfSVijay Mahadevan dmmoab->nfields = 1; 3594973de03SVijay Mahadevan 3604973de03SVijay Mahadevan /* set global ID tag handle */ 361032b8ab6SVijay Mahadevan if (!ltog_tag) { 3624973de03SVijay Mahadevan merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr); 363032b8ab6SVijay Mahadevan } 364032b8ab6SVijay Mahadevan else { 365032b8ab6SVijay Mahadevan ierr = DMMoabSetLocalToGlobalTag(*dmb, *ltog_tag);CHKERRQ(ierr); 366a4d2169cSTim Tautges } 367a4d2169cSTim Tautges 3684973de03SVijay Mahadevan /* set the local range of entities (vertices) of interest */ 369a4d2169cSTim Tautges if (range) { 3705eb88e9dSVijay Mahadevan ierr = DMMoabSetLocalVertices(*dmb, range);CHKERRQ(ierr); 371a4d2169cSTim Tautges } 3721d72bce8STim Tautges PetscFunctionReturn(0); 3731d72bce8STim Tautges } 3741d72bce8STim Tautges 3751d72bce8STim Tautges #undef __FUNCT__ 3761d72bce8STim Tautges #define __FUNCT__ "DMMoabSetParallelComm" 377aa768e4cSTim Tautges /*@ 378aa768e4cSTim Tautges DMMoabSetParallelComm - Set the ParallelComm used with this DMMoab 379aa768e4cSTim Tautges 380aa768e4cSTim Tautges Collective on MPI_Comm 381aa768e4cSTim Tautges 382aa768e4cSTim Tautges Input Parameter: 383aa768e4cSTim Tautges . dm - The DMMoab object being set 384aa768e4cSTim Tautges . pcomm - The ParallelComm being set on the DMMoab 385aa768e4cSTim Tautges 386aa768e4cSTim Tautges Level: beginner 387aa768e4cSTim Tautges 388aa768e4cSTim Tautges .keywords: DMMoab, create 389aa768e4cSTim Tautges @*/ 3901d72bce8STim Tautges PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm) 3911d72bce8STim Tautges { 392032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 393032b8ab6SVijay Mahadevan 3941d72bce8STim Tautges PetscFunctionBegin; 3951d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3961cec0304SVijay Mahadevan PetscValidPointer(pcomm,2); 397032b8ab6SVijay Mahadevan dmmoab->pcomm = pcomm; 398032b8ab6SVijay Mahadevan dmmoab->mbiface = pcomm->get_moab(); 399032b8ab6SVijay Mahadevan dmmoab->icreatedinstance = PETSC_FALSE; 4001d72bce8STim Tautges PetscFunctionReturn(0); 4011d72bce8STim Tautges } 4021d72bce8STim Tautges 4031d72bce8STim Tautges 4041d72bce8STim Tautges #undef __FUNCT__ 4051d72bce8STim Tautges #define __FUNCT__ "DMMoabGetParallelComm" 406aa768e4cSTim Tautges /*@ 407aa768e4cSTim Tautges DMMoabGetParallelComm - Get the ParallelComm used with this DMMoab 408aa768e4cSTim Tautges 409aa768e4cSTim Tautges Collective on MPI_Comm 410aa768e4cSTim Tautges 411aa768e4cSTim Tautges Input Parameter: 412aa768e4cSTim Tautges . dm - The DMMoab object being set 413aa768e4cSTim Tautges 414aa768e4cSTim Tautges Output Parameter: 415aa768e4cSTim Tautges . pcomm - The ParallelComm for the DMMoab 416aa768e4cSTim Tautges 417aa768e4cSTim Tautges Level: beginner 418aa768e4cSTim Tautges 419aa768e4cSTim Tautges .keywords: DMMoab, create 420aa768e4cSTim Tautges @*/ 4211d72bce8STim Tautges PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm) 4221d72bce8STim Tautges { 4231d72bce8STim Tautges PetscFunctionBegin; 4241d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 425032b8ab6SVijay Mahadevan *pcomm = ((DM_Moab*)(dm)->data)->pcomm; 4261d72bce8STim Tautges PetscFunctionReturn(0); 4271d72bce8STim Tautges } 4281d72bce8STim Tautges 4291d72bce8STim Tautges 4301d72bce8STim Tautges #undef __FUNCT__ 4311d72bce8STim Tautges #define __FUNCT__ "DMMoabSetInterface" 432aa768e4cSTim Tautges /*@ 433aa768e4cSTim Tautges DMMoabSetInterface - Set the MOAB instance used with this DMMoab 434aa768e4cSTim Tautges 435aa768e4cSTim Tautges Collective on MPI_Comm 436aa768e4cSTim Tautges 437aa768e4cSTim Tautges Input Parameter: 438aa768e4cSTim Tautges . dm - The DMMoab object being set 439aa768e4cSTim Tautges . mbiface - The MOAB instance being set on this DMMoab 440aa768e4cSTim Tautges 441aa768e4cSTim Tautges Level: beginner 442aa768e4cSTim Tautges 443aa768e4cSTim Tautges .keywords: DMMoab, create 444aa768e4cSTim Tautges @*/ 445a4d2169cSTim Tautges PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *mbiface) 4461d72bce8STim Tautges { 447032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 448032b8ab6SVijay Mahadevan 4491d72bce8STim Tautges PetscFunctionBegin; 4501d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4511cec0304SVijay Mahadevan PetscValidPointer(mbiface,2); 452032b8ab6SVijay Mahadevan dmmoab->pcomm = NULL; 453032b8ab6SVijay Mahadevan dmmoab->mbiface = mbiface; 454032b8ab6SVijay Mahadevan dmmoab->icreatedinstance = PETSC_FALSE; 4551d72bce8STim Tautges PetscFunctionReturn(0); 4561d72bce8STim Tautges } 4571d72bce8STim Tautges 4581d72bce8STim Tautges 4591d72bce8STim Tautges #undef __FUNCT__ 4601d72bce8STim Tautges #define __FUNCT__ "DMMoabGetInterface" 461aa768e4cSTim Tautges /*@ 462aa768e4cSTim Tautges DMMoabGetInterface - Get the MOAB instance used with this DMMoab 463aa768e4cSTim Tautges 464aa768e4cSTim Tautges Collective on MPI_Comm 465aa768e4cSTim Tautges 466aa768e4cSTim Tautges Input Parameter: 467aa768e4cSTim Tautges . dm - The DMMoab object being set 468aa768e4cSTim Tautges 469aa768e4cSTim Tautges Output Parameter: 470aa768e4cSTim Tautges . mbiface - The MOAB instance set on this DMMoab 471aa768e4cSTim Tautges 472aa768e4cSTim Tautges Level: beginner 473aa768e4cSTim Tautges 474aa768e4cSTim Tautges .keywords: DMMoab, create 475aa768e4cSTim Tautges @*/ 476a4d2169cSTim Tautges PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbiface) 4771d72bce8STim Tautges { 4789426e041SSatish Balay PetscErrorCode ierr; 479cabb514dSBarry Smith static PetscBool cite = PETSC_FALSE; 480cabb514dSBarry Smith 4811d72bce8STim Tautges PetscFunctionBegin; 4821d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 483cabb514dSBarry Smith ierr = PetscCitationsRegister("@techreport{tautges_moab:_2004,\n type = {{SAND2004-1592}},\n title = {{MOAB:} A Mesh-Oriented Database}, institution = {Sandia National Laboratories},\n author = {Tautges, T. J. and Meyers, R. and Merkley, K. and Stimpson, C. and Ernst, C.},\n year = {2004}, note = {Report}\n}\n",&cite);CHKERRQ(ierr); 484a4d2169cSTim Tautges *mbiface = ((DM_Moab*)dm->data)->mbiface; 4851d72bce8STim Tautges PetscFunctionReturn(0); 4861d72bce8STim Tautges } 4871d72bce8STim Tautges 4881d72bce8STim Tautges 4891d72bce8STim Tautges #undef __FUNCT__ 4905eb88e9dSVijay Mahadevan #define __FUNCT__ "DMMoabSetLocalVertices" 491aa768e4cSTim Tautges /*@ 4925eb88e9dSVijay Mahadevan DMMoabSetLocalVertices - Set the entities having DOFs on this DMMoab 493aa768e4cSTim Tautges 494aa768e4cSTim Tautges Collective on MPI_Comm 495aa768e4cSTim Tautges 496aa768e4cSTim Tautges Input Parameter: 497aa768e4cSTim Tautges . dm - The DMMoab object being set 498aa768e4cSTim Tautges . range - The entities treated by this DMMoab 499aa768e4cSTim Tautges 500aa768e4cSTim Tautges Level: beginner 501aa768e4cSTim Tautges 502aa768e4cSTim Tautges .keywords: DMMoab, create 503aa768e4cSTim Tautges @*/ 5045eb88e9dSVijay Mahadevan PetscErrorCode DMMoabSetLocalVertices(DM dm,moab::Range *range) 5051d72bce8STim Tautges { 506032b8ab6SVijay Mahadevan moab::ErrorCode merr; 507032b8ab6SVijay Mahadevan PetscErrorCode ierr; 508032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 509032b8ab6SVijay Mahadevan 5101d72bce8STim Tautges PetscFunctionBegin; 5111d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 512032b8ab6SVijay Mahadevan dmmoab->vlocal->clear(); 513032b8ab6SVijay Mahadevan dmmoab->vowned->clear(); 514032b8ab6SVijay Mahadevan dmmoab->vlocal->insert(range->begin(), range->end()); 515032b8ab6SVijay Mahadevan *dmmoab->vowned = *dmmoab->vlocal; 516032b8ab6SVijay Mahadevan merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 517032b8ab6SVijay Mahadevan *dmmoab->vghost = moab::subtract(*range, *dmmoab->vowned); 518032b8ab6SVijay Mahadevan dmmoab->nloc=dmmoab->vowned->size(); 519032b8ab6SVijay Mahadevan dmmoab->nghost=dmmoab->vghost->size(); 520032b8ab6SVijay Mahadevan ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 5211d72bce8STim Tautges PetscFunctionReturn(0); 5221d72bce8STim Tautges } 5231d72bce8STim Tautges 5241d72bce8STim Tautges 5251d72bce8STim Tautges #undef __FUNCT__ 5268d8d51c8SVijay Mahadevan #define __FUNCT__ "DMMoabGetAllVertices" 5278d8d51c8SVijay Mahadevan /*@ 5288d8d51c8SVijay Mahadevan DMMoabGetAllVertices - Get the entities having DOFs on this DMMoab 5298d8d51c8SVijay Mahadevan 5308d8d51c8SVijay Mahadevan Collective on MPI_Comm 5318d8d51c8SVijay Mahadevan 5328d8d51c8SVijay Mahadevan Input Parameter: 5338d8d51c8SVijay Mahadevan . dm - The DMMoab object being set 5348d8d51c8SVijay Mahadevan 5358d8d51c8SVijay Mahadevan Output Parameter: 5368d8d51c8SVijay Mahadevan . owned - The local vertex entities in this DMMoab = (owned+ghosted) 5378d8d51c8SVijay Mahadevan 5388d8d51c8SVijay Mahadevan Level: beginner 5398d8d51c8SVijay Mahadevan 5408d8d51c8SVijay Mahadevan .keywords: DMMoab, create 5418d8d51c8SVijay Mahadevan @*/ 5428d8d51c8SVijay Mahadevan PetscErrorCode DMMoabGetAllVertices(DM dm,moab::Range *local) 5438d8d51c8SVijay Mahadevan { 5448d8d51c8SVijay Mahadevan PetscFunctionBegin; 5458d8d51c8SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5468d8d51c8SVijay Mahadevan if (local) *local = *((DM_Moab*)dm->data)->vlocal; 5478d8d51c8SVijay Mahadevan PetscFunctionReturn(0); 5488d8d51c8SVijay Mahadevan } 5498d8d51c8SVijay Mahadevan 5508d8d51c8SVijay Mahadevan 5518d8d51c8SVijay Mahadevan 5528d8d51c8SVijay Mahadevan #undef __FUNCT__ 5535eb88e9dSVijay Mahadevan #define __FUNCT__ "DMMoabGetLocalVertices" 554aa768e4cSTim Tautges /*@ 5555eb88e9dSVijay Mahadevan DMMoabGetLocalVertices - Get the entities having DOFs on this DMMoab 556aa768e4cSTim Tautges 557aa768e4cSTim Tautges Collective on MPI_Comm 558aa768e4cSTim Tautges 559aa768e4cSTim Tautges Input Parameter: 560aa768e4cSTim Tautges . dm - The DMMoab object being set 561aa768e4cSTim Tautges 562aa768e4cSTim Tautges Output Parameter: 5635eb88e9dSVijay Mahadevan . owned - The owned vertex entities in this DMMoab 5645eb88e9dSVijay Mahadevan . ghost - The ghosted entities (non-owned) stored locally in this partition 565aa768e4cSTim Tautges 566aa768e4cSTim Tautges Level: beginner 567aa768e4cSTim Tautges 568aa768e4cSTim Tautges .keywords: DMMoab, create 569aa768e4cSTim Tautges @*/ 5701cec0304SVijay Mahadevan PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost) 5711d72bce8STim Tautges { 5721d72bce8STim Tautges PetscFunctionBegin; 5731d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5741cec0304SVijay Mahadevan if (owned) *owned = *((DM_Moab*)dm->data)->vowned; 5751cec0304SVijay Mahadevan if (ghost) *ghost = *((DM_Moab*)dm->data)->vghost; 5761d72bce8STim Tautges PetscFunctionReturn(0); 5771d72bce8STim Tautges } 5781d72bce8STim Tautges 5791d72bce8STim Tautges #undef __FUNCT__ 5805eb88e9dSVijay Mahadevan #define __FUNCT__ "DMMoabGetLocalElements" 5815eb88e9dSVijay Mahadevan /*@ 5825eb88e9dSVijay Mahadevan DMMoabGetLocalElements - Get the higher-dimensional entities that are locally owned 5835eb88e9dSVijay Mahadevan 5845eb88e9dSVijay Mahadevan Collective on MPI_Comm 5855eb88e9dSVijay Mahadevan 5865eb88e9dSVijay Mahadevan Input Parameter: 5875eb88e9dSVijay Mahadevan . dm - The DMMoab object being set 5885eb88e9dSVijay Mahadevan 5895eb88e9dSVijay Mahadevan Output Parameter: 5905eb88e9dSVijay Mahadevan . range - The entities owned locally 5915eb88e9dSVijay Mahadevan 5925eb88e9dSVijay Mahadevan Level: beginner 5935eb88e9dSVijay Mahadevan 5945eb88e9dSVijay Mahadevan .keywords: DMMoab, create 5955eb88e9dSVijay Mahadevan @*/ 5961cec0304SVijay Mahadevan PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range) 5975eb88e9dSVijay Mahadevan { 5985eb88e9dSVijay Mahadevan PetscFunctionBegin; 5995eb88e9dSVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6001cec0304SVijay Mahadevan if (range) *range = *((DM_Moab*)dm->data)->elocal; 6011cec0304SVijay Mahadevan PetscFunctionReturn(0); 6021cec0304SVijay Mahadevan } 6031cec0304SVijay Mahadevan 6041cec0304SVijay Mahadevan 6051cec0304SVijay Mahadevan #undef __FUNCT__ 6061cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabSetLocalElements" 6071cec0304SVijay Mahadevan /*@ 6081cec0304SVijay Mahadevan DMMoabSetLocalElements - Set the entities having DOFs on this DMMoab 6091cec0304SVijay Mahadevan 6101cec0304SVijay Mahadevan Collective on MPI_Comm 6111cec0304SVijay Mahadevan 6121cec0304SVijay Mahadevan Input Parameter: 6131cec0304SVijay Mahadevan . dm - The DMMoab object being set 6141cec0304SVijay Mahadevan . range - The entities treated by this DMMoab 6151cec0304SVijay Mahadevan 6161cec0304SVijay Mahadevan Level: beginner 6171cec0304SVijay Mahadevan 6181cec0304SVijay Mahadevan .keywords: DMMoab, create 6191cec0304SVijay Mahadevan @*/ 6201cec0304SVijay Mahadevan PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range) 6211cec0304SVijay Mahadevan { 6221cec0304SVijay Mahadevan moab::ErrorCode merr; 6231cec0304SVijay Mahadevan PetscErrorCode ierr; 6241cec0304SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)(dm)->data; 6251cec0304SVijay Mahadevan 6261cec0304SVijay Mahadevan PetscFunctionBegin; 6271cec0304SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6281cec0304SVijay Mahadevan dmmoab->elocal->clear(); 6291cec0304SVijay Mahadevan dmmoab->eghost->clear(); 6301cec0304SVijay Mahadevan dmmoab->elocal->insert(range->begin(), range->end()); 6311cec0304SVijay Mahadevan PetscInfo2(dm, "Range size = %D; elocal size = %D.\n", range->size(), dmmoab->elocal->size()); 6321cec0304SVijay Mahadevan merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr); 6331cec0304SVijay Mahadevan *dmmoab->eghost = moab::subtract(*range, *dmmoab->elocal); 6341cec0304SVijay Mahadevan dmmoab->neleloc=dmmoab->elocal->size(); 6351cec0304SVijay Mahadevan ierr = MPI_Allreduce(&dmmoab->nele, &dmmoab->neleloc, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr); 6361cec0304SVijay Mahadevan PetscInfo2(dm, "Created %D local and %D glocal elements.\n", dmmoab->neleloc, dmmoab->nele); 6375eb88e9dSVijay Mahadevan PetscFunctionReturn(0); 6385eb88e9dSVijay Mahadevan } 6395eb88e9dSVijay Mahadevan 6405eb88e9dSVijay Mahadevan 6415eb88e9dSVijay Mahadevan #undef __FUNCT__ 6421d72bce8STim Tautges #define __FUNCT__ "DMMoabSetLocalToGlobalTag" 643aa768e4cSTim Tautges /*@ 644aa768e4cSTim Tautges DMMoabSetLocalToGlobalTag - Set the tag used for local to global numbering 645aa768e4cSTim Tautges 646aa768e4cSTim Tautges Collective on MPI_Comm 647aa768e4cSTim Tautges 648aa768e4cSTim Tautges Input Parameter: 649aa768e4cSTim Tautges . dm - The DMMoab object being set 650aa768e4cSTim Tautges . ltogtag - The MOAB tag used for local to global ids 651aa768e4cSTim Tautges 652aa768e4cSTim Tautges Level: beginner 653aa768e4cSTim Tautges 654aa768e4cSTim Tautges .keywords: DMMoab, create 655aa768e4cSTim Tautges @*/ 6561d72bce8STim Tautges PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag) 6571d72bce8STim Tautges { 6581d72bce8STim Tautges PetscFunctionBegin; 6591d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6601d72bce8STim Tautges ((DM_Moab*)dm->data)->ltog_tag = ltogtag; 6611d72bce8STim Tautges PetscFunctionReturn(0); 6621d72bce8STim Tautges } 6631d72bce8STim Tautges 6641d72bce8STim Tautges 6651d72bce8STim Tautges #undef __FUNCT__ 6661d72bce8STim Tautges #define __FUNCT__ "DMMoabGetLocalToGlobalTag" 667aa768e4cSTim Tautges /*@ 668aa768e4cSTim Tautges DMMoabGetLocalToGlobalTag - Get the tag used for local to global numbering 669aa768e4cSTim Tautges 670aa768e4cSTim Tautges Collective on MPI_Comm 671aa768e4cSTim Tautges 672aa768e4cSTim Tautges Input Parameter: 673aa768e4cSTim Tautges . dm - The DMMoab object being set 674aa768e4cSTim Tautges 675aa768e4cSTim Tautges Output Parameter: 676aa768e4cSTim Tautges . ltogtag - The MOAB tag used for local to global ids 677aa768e4cSTim Tautges 678aa768e4cSTim Tautges Level: beginner 679aa768e4cSTim Tautges 680aa768e4cSTim Tautges .keywords: DMMoab, create 681aa768e4cSTim Tautges @*/ 6821d72bce8STim Tautges PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag) 6831d72bce8STim Tautges { 6841d72bce8STim Tautges PetscFunctionBegin; 6851d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6861d72bce8STim Tautges *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag; 6871d72bce8STim Tautges PetscFunctionReturn(0); 6881d72bce8STim Tautges } 6891d72bce8STim Tautges 6901d72bce8STim Tautges 6911d72bce8STim Tautges #undef __FUNCT__ 6921d72bce8STim Tautges #define __FUNCT__ "DMMoabSetBlockSize" 693aa768e4cSTim Tautges /*@ 694aa768e4cSTim Tautges DMMoabSetBlockSize - Set the block size used with this DMMoab 695aa768e4cSTim Tautges 696aa768e4cSTim Tautges Collective on MPI_Comm 697aa768e4cSTim Tautges 698aa768e4cSTim Tautges Input Parameter: 699aa768e4cSTim Tautges . dm - The DMMoab object being set 700aa768e4cSTim Tautges . bs - The block size used with this DMMoab 701aa768e4cSTim Tautges 702aa768e4cSTim Tautges Level: beginner 703aa768e4cSTim Tautges 704aa768e4cSTim Tautges .keywords: DMMoab, create 705aa768e4cSTim Tautges @*/ 7061d72bce8STim Tautges PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs) 7071d72bce8STim Tautges { 7081d72bce8STim Tautges PetscFunctionBegin; 7091d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 7101d72bce8STim Tautges ((DM_Moab*)dm->data)->bs = bs; 7111d72bce8STim Tautges PetscFunctionReturn(0); 7121d72bce8STim Tautges } 7131d72bce8STim Tautges 7141d72bce8STim Tautges 7151d72bce8STim Tautges #undef __FUNCT__ 7161d72bce8STim Tautges #define __FUNCT__ "DMMoabGetBlockSize" 717aa768e4cSTim Tautges /*@ 718aa768e4cSTim Tautges DMMoabGetBlockSize - Get the block size used with this DMMoab 719aa768e4cSTim Tautges 720aa768e4cSTim Tautges Collective on MPI_Comm 721aa768e4cSTim Tautges 722aa768e4cSTim Tautges Input Parameter: 723aa768e4cSTim Tautges . dm - The DMMoab object being set 724aa768e4cSTim Tautges 725aa768e4cSTim Tautges Output Parameter: 726aa768e4cSTim Tautges . bs - The block size used with this DMMoab 727aa768e4cSTim Tautges 728aa768e4cSTim Tautges Level: beginner 729aa768e4cSTim Tautges 730aa768e4cSTim Tautges .keywords: DMMoab, create 731aa768e4cSTim Tautges @*/ 7321d72bce8STim Tautges PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs) 7331d72bce8STim Tautges { 7341d72bce8STim Tautges PetscFunctionBegin; 7351d72bce8STim Tautges PetscValidHeaderSpecific(dm,DM_CLASSID,1); 7361d72bce8STim Tautges *bs = ((DM_Moab*)dm->data)->bs; 7371d72bce8STim Tautges PetscFunctionReturn(0); 7381d72bce8STim Tautges } 7391d72bce8STim Tautges 7401cec0304SVijay Mahadevan 7411cec0304SVijay Mahadevan #undef __FUNCT__ 742212ad6d1SVijay Mahadevan #define __FUNCT__ "DMMoabGetSize" 743212ad6d1SVijay Mahadevan /*@ 744212ad6d1SVijay Mahadevan DMMoabGetSize - Get the global vertex size used with this DMMoab 745212ad6d1SVijay Mahadevan 746212ad6d1SVijay Mahadevan Collective on MPI_Comm 747212ad6d1SVijay Mahadevan 748212ad6d1SVijay Mahadevan Input Parameter: 749212ad6d1SVijay Mahadevan . dm - The DMMoab object being set 750212ad6d1SVijay Mahadevan 751212ad6d1SVijay Mahadevan Output Parameter: 752212ad6d1SVijay Mahadevan . ng - The global size of the DMMoab instance 753212ad6d1SVijay Mahadevan 754212ad6d1SVijay Mahadevan Level: beginner 755212ad6d1SVijay Mahadevan 756212ad6d1SVijay Mahadevan .keywords: DMMoab, create 757212ad6d1SVijay Mahadevan @*/ 758212ad6d1SVijay Mahadevan PetscErrorCode DMMoabGetSize(DM dm,PetscInt *ng) 759212ad6d1SVijay Mahadevan { 760212ad6d1SVijay Mahadevan PetscFunctionBegin; 761212ad6d1SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 762212ad6d1SVijay Mahadevan if(ng) *ng = ((DM_Moab*)dm->data)->n; 763212ad6d1SVijay Mahadevan PetscFunctionReturn(0); 764212ad6d1SVijay Mahadevan } 765212ad6d1SVijay Mahadevan 766212ad6d1SVijay Mahadevan 767212ad6d1SVijay Mahadevan #undef __FUNCT__ 768212ad6d1SVijay Mahadevan #define __FUNCT__ "DMMoabGetLocalSize" 769212ad6d1SVijay Mahadevan /*@ 770212ad6d1SVijay Mahadevan DMMoabGetLocalSize - Get the local and ghosted vertex size used with this DMMoab 771212ad6d1SVijay Mahadevan 772212ad6d1SVijay Mahadevan Collective on MPI_Comm 773212ad6d1SVijay Mahadevan 774212ad6d1SVijay Mahadevan Input Parameter: 775212ad6d1SVijay Mahadevan . dm - The DMMoab object being set 776212ad6d1SVijay Mahadevan 777212ad6d1SVijay Mahadevan Output Parameter: 778212ad6d1SVijay Mahadevan . nl - The local size of the DMMoab instance 779212ad6d1SVijay Mahadevan . ng - The ghosted size of the DMMoab instance 780212ad6d1SVijay Mahadevan 781212ad6d1SVijay Mahadevan Level: beginner 782212ad6d1SVijay Mahadevan 783212ad6d1SVijay Mahadevan .keywords: DMMoab, create 784212ad6d1SVijay Mahadevan @*/ 785212ad6d1SVijay Mahadevan PetscErrorCode DMMoabGetLocalSize(DM dm,PetscInt *nl,PetscInt *ng) 786212ad6d1SVijay Mahadevan { 787212ad6d1SVijay Mahadevan PetscFunctionBegin; 788212ad6d1SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 789212ad6d1SVijay Mahadevan if(nl) *nl = ((DM_Moab*)dm->data)->nloc; 790212ad6d1SVijay Mahadevan if(ng) *ng = ((DM_Moab*)dm->data)->nghost; 791212ad6d1SVijay Mahadevan PetscFunctionReturn(0); 792212ad6d1SVijay Mahadevan } 793212ad6d1SVijay Mahadevan 794212ad6d1SVijay Mahadevan 795212ad6d1SVijay Mahadevan #undef __FUNCT__ 7964920ab11SVijay Mahadevan #define __FUNCT__ "DMMoabGetDimension" 7974920ab11SVijay Mahadevan /*@ 7984920ab11SVijay Mahadevan DMMoabGetDimension - Get the dimension of the DM Mesh 7994920ab11SVijay Mahadevan 8004920ab11SVijay Mahadevan Collective on MPI_Comm 8014920ab11SVijay Mahadevan 8024920ab11SVijay Mahadevan Input Parameter: 8034920ab11SVijay Mahadevan . dm - The DMMoab object being set 8044920ab11SVijay Mahadevan 8054920ab11SVijay Mahadevan Output Parameter: 8064920ab11SVijay Mahadevan . dim - The dimension of DM 8074920ab11SVijay Mahadevan 8084920ab11SVijay Mahadevan Level: beginner 8094920ab11SVijay Mahadevan 8104920ab11SVijay Mahadevan .keywords: DMMoab, create 8114920ab11SVijay Mahadevan @*/ 8124920ab11SVijay Mahadevan PetscErrorCode DMMoabGetDimension(DM dm,PetscInt *dim) 8134920ab11SVijay Mahadevan { 8144920ab11SVijay Mahadevan PetscFunctionBegin; 8154920ab11SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 8164920ab11SVijay Mahadevan *dim = ((DM_Moab*)dm->data)->dim; 8174920ab11SVijay Mahadevan PetscFunctionReturn(0); 8184920ab11SVijay Mahadevan } 8194920ab11SVijay Mahadevan 8204920ab11SVijay Mahadevan 8214920ab11SVijay Mahadevan 8224920ab11SVijay Mahadevan #undef __FUNCT__ 8231cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabSetFieldVector" 8241cec0304SVijay Mahadevan PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec) 8251cec0304SVijay Mahadevan { 8261cec0304SVijay Mahadevan DM_Moab *dmmoab; 8271cec0304SVijay Mahadevan moab::Tag vtag,ntag; 828212ad6d1SVijay Mahadevan const PetscScalar *varray; 829212ad6d1SVijay Mahadevan PetscScalar *farray; 8301cec0304SVijay Mahadevan moab::ErrorCode merr; 8311cec0304SVijay Mahadevan PetscErrorCode ierr; 8321cec0304SVijay Mahadevan PetscInt doff; 8331cec0304SVijay Mahadevan std::string tag_name; 8341cec0304SVijay Mahadevan moab::Range::iterator iter; 8351cec0304SVijay Mahadevan 8361cec0304SVijay Mahadevan PetscFunctionBegin; 8371cec0304SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 8381cec0304SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 8391cec0304SVijay Mahadevan 8401cec0304SVijay Mahadevan /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 8411cec0304SVijay Mahadevan merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 8421cec0304SVijay Mahadevan moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 8431cec0304SVijay Mahadevan 8441cec0304SVijay Mahadevan ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 8451cec0304SVijay Mahadevan 8461cec0304SVijay Mahadevan merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 8471cec0304SVijay Mahadevan if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 848212ad6d1SVijay Mahadevan ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 849212ad6d1SVijay Mahadevan 8501cec0304SVijay Mahadevan for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) { 8511cec0304SVijay Mahadevan moab::EntityHandle vtx = (*iter); 8521cec0304SVijay Mahadevan 8531cec0304SVijay Mahadevan /* get field dof index */ 854212ad6d1SVijay Mahadevan ierr = DMMoabGetFieldDof(dm, vtx, ifield, &doff); 8551cec0304SVijay Mahadevan 8561cec0304SVijay Mahadevan /* use the entity handle and the Dof index to set the right value */ 8571cec0304SVijay Mahadevan merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr); 8581cec0304SVijay Mahadevan } 859212ad6d1SVijay Mahadevan ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 8601cec0304SVijay Mahadevan } 8611cec0304SVijay Mahadevan else { 862212ad6d1SVijay Mahadevan ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr); 8631cec0304SVijay Mahadevan /* we are using a MOAB Vec - directly copy the tag data to new one */ 864212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)farray);MBERRNM(merr); 865212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 8661cec0304SVijay Mahadevan /* make sure the parallel exchange for ghosts are done appropriately */ 8671cec0304SVijay Mahadevan merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr); 868212ad6d1SVijay Mahadevan ierr = PetscFree(farray);CHKERRQ(ierr); 869212ad6d1SVijay Mahadevan } 870212ad6d1SVijay Mahadevan PetscFunctionReturn(0); 871212ad6d1SVijay Mahadevan } 872212ad6d1SVijay Mahadevan 873212ad6d1SVijay Mahadevan 874212ad6d1SVijay Mahadevan #undef __FUNCT__ 875212ad6d1SVijay Mahadevan #define __FUNCT__ "DMMoabSetGlobalFieldVector" 876212ad6d1SVijay Mahadevan PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec) 877212ad6d1SVijay Mahadevan { 878212ad6d1SVijay Mahadevan DM_Moab *dmmoab; 879212ad6d1SVijay Mahadevan moab::Tag vtag,ntag; 880212ad6d1SVijay Mahadevan const PetscScalar *varray; 881212ad6d1SVijay Mahadevan PetscScalar *farray; 882212ad6d1SVijay Mahadevan moab::ErrorCode merr; 883212ad6d1SVijay Mahadevan PetscErrorCode ierr; 884212ad6d1SVijay Mahadevan PetscSection section; 885212ad6d1SVijay Mahadevan PetscInt i,doff,ifield; 886212ad6d1SVijay Mahadevan std::string tag_name; 887212ad6d1SVijay Mahadevan moab::Range::iterator iter; 888212ad6d1SVijay Mahadevan 889212ad6d1SVijay Mahadevan PetscFunctionBegin; 890212ad6d1SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 891212ad6d1SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 892212ad6d1SVijay Mahadevan 893212ad6d1SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 894212ad6d1SVijay Mahadevan 895212ad6d1SVijay Mahadevan /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */ 896212ad6d1SVijay Mahadevan ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 897212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 898212ad6d1SVijay Mahadevan if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 899212ad6d1SVijay Mahadevan /* not a MOAB vector - use VecGetSubVector to get the parts as needed */ 900212ad6d1SVijay Mahadevan 901212ad6d1SVijay Mahadevan ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 902212ad6d1SVijay Mahadevan for (ifield=0; ifield<dmmoab->nfields; ++ifield) { 903212ad6d1SVijay Mahadevan 904212ad6d1SVijay Mahadevan /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 905212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 906212ad6d1SVijay Mahadevan moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 907212ad6d1SVijay Mahadevan 908212ad6d1SVijay Mahadevan for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) { 909212ad6d1SVijay Mahadevan moab::EntityHandle vtx = (*iter); 910212ad6d1SVijay Mahadevan 911212ad6d1SVijay Mahadevan /* get field dof index */ 912212ad6d1SVijay Mahadevan ierr = DMMoabGetFieldDof(dm, vtx, ifield, &doff); 913212ad6d1SVijay Mahadevan 914212ad6d1SVijay Mahadevan /* use the entity handle and the Dof index to set the right value */ 915212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr); 916212ad6d1SVijay Mahadevan } 917212ad6d1SVijay Mahadevan } 918212ad6d1SVijay Mahadevan ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 919212ad6d1SVijay Mahadevan } 920212ad6d1SVijay Mahadevan else { 921212ad6d1SVijay Mahadevan ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr); 922212ad6d1SVijay Mahadevan ierr = PetscMalloc(dmmoab->nloc*dmmoab->bs*sizeof(PetscScalar),&varray);CHKERRQ(ierr); 923212ad6d1SVijay Mahadevan 924212ad6d1SVijay Mahadevan /* we are using a MOAB Vec - directly copy the tag data to new one */ 925212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr); 926212ad6d1SVijay Mahadevan for (ifield=0; ifield<dmmoab->nfields; ++ifield) { 927212ad6d1SVijay Mahadevan 928212ad6d1SVijay Mahadevan /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 929212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 930212ad6d1SVijay Mahadevan moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr); 931212ad6d1SVijay Mahadevan 932212ad6d1SVijay Mahadevan /* we are using a MOAB Vec - directly copy the tag data to new one */ 933212ad6d1SVijay Mahadevan for(i=0; i < dmmoab->nloc; i++) { 934212ad6d1SVijay Mahadevan farray[i] = varray[i*dmmoab->bs+ifield]; 935212ad6d1SVijay Mahadevan } 936212ad6d1SVijay Mahadevan 937212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 938212ad6d1SVijay Mahadevan /* make sure the parallel exchange for ghosts are done appropriately */ 939212ad6d1SVijay Mahadevan merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr); 940212ad6d1SVijay Mahadevan } 941212ad6d1SVijay Mahadevan ierr = PetscFree(farray);CHKERRQ(ierr); 9421cec0304SVijay Mahadevan ierr = PetscFree(varray);CHKERRQ(ierr); 9431cec0304SVijay Mahadevan } 9441cec0304SVijay Mahadevan PetscFunctionReturn(0); 9451cec0304SVijay Mahadevan } 9461cec0304SVijay Mahadevan 9471cec0304SVijay Mahadevan 948212ad6d1SVijay Mahadevan 9491cec0304SVijay Mahadevan #undef __FUNCT__ 9507023aa44SVijay Mahadevan #define __FUNCT__ "DMMoabGetVertexCoordinates" 9517023aa44SVijay Mahadevan PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos) 9527023aa44SVijay Mahadevan { 9537023aa44SVijay Mahadevan DM_Moab *dmmoab; 9547023aa44SVijay Mahadevan PetscErrorCode ierr; 9557023aa44SVijay Mahadevan moab::ErrorCode merr; 9567023aa44SVijay Mahadevan 9577023aa44SVijay Mahadevan PetscFunctionBegin; 9587023aa44SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 9597023aa44SVijay Mahadevan PetscValidPointer(conn,3); 9607023aa44SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 9617023aa44SVijay Mahadevan 9627023aa44SVijay Mahadevan if (!vpos) { 9637023aa44SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscScalar)*nconn*3, &vpos);CHKERRQ(ierr); 9647023aa44SVijay Mahadevan } 9657023aa44SVijay Mahadevan 9667023aa44SVijay Mahadevan /* Get connectivity information in MOAB canonical ordering */ 9677023aa44SVijay Mahadevan merr = dmmoab->mbiface->get_coords(conn, nconn, vpos);MBERRNM(merr); 9687023aa44SVijay Mahadevan PetscFunctionReturn(0); 9697023aa44SVijay Mahadevan } 9707023aa44SVijay Mahadevan 9717023aa44SVijay Mahadevan 9727023aa44SVijay Mahadevan #undef __FUNCT__ 9738d8d51c8SVijay Mahadevan #define __FUNCT__ "DMMoabGetVertexConnectivity" 9748d8d51c8SVijay Mahadevan PetscErrorCode DMMoabGetVertexConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn, moab::EntityHandle **conn) 9758d8d51c8SVijay Mahadevan { 9768d8d51c8SVijay Mahadevan DM_Moab *dmmoab; 9778d8d51c8SVijay Mahadevan std::vector<moab::EntityHandle> adj_entities,connect; 9788d8d51c8SVijay Mahadevan PetscErrorCode ierr; 9798d8d51c8SVijay Mahadevan moab::ErrorCode merr; 9808d8d51c8SVijay Mahadevan 9818d8d51c8SVijay Mahadevan PetscFunctionBegin; 9828d8d51c8SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 9838d8d51c8SVijay Mahadevan PetscValidPointer(conn,4); 9848d8d51c8SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 9858d8d51c8SVijay Mahadevan 9868d8d51c8SVijay Mahadevan /* Get connectivity information in MOAB canonical ordering */ 9878d8d51c8SVijay Mahadevan merr = dmmoab->mbiface->get_adjacencies(&ehandle, 1, 1, true, adj_entities, moab::Interface::UNION);MBERRNM(merr); 9888d8d51c8SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity(&adj_entities[0],adj_entities.size(),connect);MBERRNM(merr); 9898d8d51c8SVijay Mahadevan 9908d8d51c8SVijay Mahadevan #if 0 9918d8d51c8SVijay Mahadevan for(unsigned int jter = 0; jter < connect.size(); jter++) { 9928d8d51c8SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF,"Handle=%D\tAdj_Size=%D\tAdj_Entity=%D\n",ehandle,connect.size(),connect[jter]); 9938d8d51c8SVijay Mahadevan } 9948d8d51c8SVijay Mahadevan #endif 9958d8d51c8SVijay Mahadevan 9968d8d51c8SVijay Mahadevan if (conn) { 9978d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(moab::EntityHandle)*connect.size(), conn);CHKERRQ(ierr); 9988d8d51c8SVijay Mahadevan ierr = PetscMemcpy(*conn, &connect[0], sizeof(moab::EntityHandle)*connect.size());CHKERRQ(ierr); 9998d8d51c8SVijay Mahadevan } 10008d8d51c8SVijay Mahadevan if (nconn) *nconn=connect.size(); 10018d8d51c8SVijay Mahadevan PetscFunctionReturn(0); 10028d8d51c8SVijay Mahadevan } 10038d8d51c8SVijay Mahadevan 10048d8d51c8SVijay Mahadevan 10058d8d51c8SVijay Mahadevan #undef __FUNCT__ 10068d8d51c8SVijay Mahadevan #define __FUNCT__ "DMMoabRestoreVertexConnectivity" 10078d8d51c8SVijay Mahadevan PetscErrorCode DMMoabRestoreVertexConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn, moab::EntityHandle **conn) 10088d8d51c8SVijay Mahadevan { 10098d8d51c8SVijay Mahadevan PetscErrorCode ierr; 10108d8d51c8SVijay Mahadevan 10118d8d51c8SVijay Mahadevan PetscFunctionBegin; 10128d8d51c8SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 10138d8d51c8SVijay Mahadevan PetscValidPointer(conn,4); 10148d8d51c8SVijay Mahadevan 10158d8d51c8SVijay Mahadevan if (conn) { 10168d8d51c8SVijay Mahadevan ierr = PetscFree(*conn);CHKERRQ(ierr); 10178d8d51c8SVijay Mahadevan } 10188d8d51c8SVijay Mahadevan if (nconn) *nconn=0; 10198d8d51c8SVijay Mahadevan PetscFunctionReturn(0); 10208d8d51c8SVijay Mahadevan } 10218d8d51c8SVijay Mahadevan 10228d8d51c8SVijay Mahadevan 10238d8d51c8SVijay Mahadevan 10248d8d51c8SVijay Mahadevan #undef __FUNCT__ 10257023aa44SVijay Mahadevan #define __FUNCT__ "DMMoabGetElementConnectivity" 10267023aa44SVijay Mahadevan PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn) 10277023aa44SVijay Mahadevan { 10287023aa44SVijay Mahadevan DM_Moab *dmmoab; 10297023aa44SVijay Mahadevan const moab::EntityHandle *connect; 10307023aa44SVijay Mahadevan moab::ErrorCode merr; 10317023aa44SVijay Mahadevan PetscInt nnodes; 10327023aa44SVijay Mahadevan 10337023aa44SVijay Mahadevan PetscFunctionBegin; 10347023aa44SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 10357023aa44SVijay Mahadevan PetscValidPointer(conn,4); 10367023aa44SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 10377023aa44SVijay Mahadevan 10387023aa44SVijay Mahadevan /* Get connectivity information in MOAB canonical ordering */ 10397023aa44SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity(ehandle, connect, nnodes);MBERRNM(merr); 10407023aa44SVijay Mahadevan if (conn) *conn=connect; 10417023aa44SVijay Mahadevan if (nconn) *nconn=nnodes; 10427023aa44SVijay Mahadevan PetscFunctionReturn(0); 10437023aa44SVijay Mahadevan } 10447023aa44SVijay Mahadevan 10457023aa44SVijay Mahadevan 10467023aa44SVijay Mahadevan #undef __FUNCT__ 104769263071SVijay Mahadevan #define __FUNCT__ "DMMoabIsEntityOnBoundary" 104869263071SVijay Mahadevan PetscErrorCode DMMoabIsEntityOnBoundary(DM dm,const moab::EntityHandle ent,PetscBool* ent_on_boundary) 104969263071SVijay Mahadevan { 105069263071SVijay Mahadevan moab::EntityType etype; 105169263071SVijay Mahadevan DM_Moab *dmmoab; 105269263071SVijay Mahadevan PetscInt edim; 105369263071SVijay Mahadevan 105469263071SVijay Mahadevan PetscFunctionBegin; 105569263071SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 105669263071SVijay Mahadevan PetscValidPointer(ent_on_boundary,3); 105769263071SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 105869263071SVijay Mahadevan 105969263071SVijay Mahadevan /* get the entity type and handle accordingly */ 106069263071SVijay Mahadevan etype=dmmoab->mbiface->type_from_handle(ent); 106169263071SVijay Mahadevan if(etype >= moab::MBPOLYHEDRON) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Entity type on the boundary skin is invalid. EntityType = %D\n",etype); 106269263071SVijay Mahadevan 106369263071SVijay Mahadevan /* get the entity dimension */ 106469263071SVijay Mahadevan edim=dmmoab->mbiface->dimension_from_handle(ent); 106569263071SVijay Mahadevan 106669263071SVijay Mahadevan *ent_on_boundary=PETSC_FALSE; 106769263071SVijay Mahadevan if(etype == moab::MBVERTEX && edim == 0) { 10688d8d51c8SVijay Mahadevan if (ent < (*dmmoab->vlocal)[0] || ent > (*dmmoab->vlocal)[dmmoab->nloc-1]) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid boundary vertex entity handle: %D\n",ent); 1069*0c8a2322SVijay Mahadevan *ent_on_boundary=dmmoab->isbndyvtx[(PetscInt)ent]; 107069263071SVijay Mahadevan } 107169263071SVijay Mahadevan else { 107269263071SVijay Mahadevan if (edim == dmmoab->dim) { /* check the higher-dimensional elements first */ 10738d8d51c8SVijay Mahadevan if (ent < (*dmmoab->elocal)[0] || ent > (*dmmoab->elocal)[dmmoab->neleloc-1]) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid boundary element entity handle: %D\n",ent); 1074*0c8a2322SVijay Mahadevan *ent_on_boundary=dmmoab->isbndyelems[(PetscInt)ent]; 107569263071SVijay Mahadevan } 107669263071SVijay Mahadevan else { /* next check the lower-dimensional faces */ 1077*0c8a2322SVijay Mahadevan /* how do we check the bounds before accessing ? will segfault for non-boundary faces */ 1078*0c8a2322SVijay Mahadevan *ent_on_boundary=dmmoab->isbndyfaces[(PetscInt)ent]; 107969263071SVijay Mahadevan } 108069263071SVijay Mahadevan } 108169263071SVijay Mahadevan PetscFunctionReturn(0); 108269263071SVijay Mahadevan } 108369263071SVijay Mahadevan 108469263071SVijay Mahadevan 108569263071SVijay Mahadevan #undef __FUNCT__ 10867023aa44SVijay Mahadevan #define __FUNCT__ "DMMoabCheckBoundaryVertices" 108769263071SVijay Mahadevan PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx) 10887023aa44SVijay Mahadevan { 10897023aa44SVijay Mahadevan DM_Moab *dmmoab; 10907023aa44SVijay Mahadevan PetscInt i; 10917023aa44SVijay Mahadevan 10927023aa44SVijay Mahadevan PetscFunctionBegin; 10937023aa44SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 10947023aa44SVijay Mahadevan PetscValidPointer(cnt,3); 10957023aa44SVijay Mahadevan PetscValidPointer(isbdvtx,4); 10967023aa44SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 10977023aa44SVijay Mahadevan 10987023aa44SVijay Mahadevan for (i=0; i < nconn; ++i) { 1099*0c8a2322SVijay Mahadevan isbdvtx[i]=dmmoab->isbndyvtx[(PetscInt)cnt[i]]; 11007023aa44SVijay Mahadevan } 11017023aa44SVijay Mahadevan PetscFunctionReturn(0); 11027023aa44SVijay Mahadevan } 11037023aa44SVijay Mahadevan 11047023aa44SVijay Mahadevan 11057023aa44SVijay Mahadevan #undef __FUNCT__ 1106*0c8a2322SVijay Mahadevan #define __FUNCT__ "DMMoabGetBoundaryMarkers" 1107*0c8a2322SVijay Mahadevan PetscErrorCode DMMoabGetBoundaryMarkers(DM dm,PetscBool **bdvtx,PetscBool** bdelems,PetscBool** bdfaces) 11081cec0304SVijay Mahadevan { 11091cec0304SVijay Mahadevan DM_Moab *dmmoab; 11101cec0304SVijay Mahadevan 11111cec0304SVijay Mahadevan PetscFunctionBegin; 11121cec0304SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 11131cec0304SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 11141cec0304SVijay Mahadevan 1115*0c8a2322SVijay Mahadevan if (bdvtx) *bdvtx = dmmoab->isbndyvtx; 1116*0c8a2322SVijay Mahadevan if (bdfaces) *bdfaces = dmmoab->isbndyfaces; 1117*0c8a2322SVijay Mahadevan if (bdelems) *bdfaces = dmmoab->isbndyelems; 11181cec0304SVijay Mahadevan PetscFunctionReturn(0); 11191cec0304SVijay Mahadevan } 11201cec0304SVijay Mahadevan 11211cec0304SVijay Mahadevan 11221cec0304SVijay Mahadevan #undef __FUNCT__ 11231cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabSetFields" 11241cec0304SVijay Mahadevan PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields) 11251cec0304SVijay Mahadevan { 11261cec0304SVijay Mahadevan DM_Moab *dmmoab; 11271cec0304SVijay Mahadevan 11281cec0304SVijay Mahadevan PetscFunctionBegin; 11291cec0304SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 11301cec0304SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 11311cec0304SVijay Mahadevan 11321cec0304SVijay Mahadevan dmmoab->fields = fields; 11331cec0304SVijay Mahadevan dmmoab->nfields = nfields; 11341cec0304SVijay Mahadevan PetscFunctionReturn(0); 11351cec0304SVijay Mahadevan } 11361cec0304SVijay Mahadevan 11371cec0304SVijay Mahadevan 11381cec0304SVijay Mahadevan #undef __FUNCT__ 11391cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldDof" 11401cec0304SVijay Mahadevan PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof) 11411cec0304SVijay Mahadevan { 11421cec0304SVijay Mahadevan PetscSection section; 1143fc418013SVijay Mahadevan PetscInt gid; 11441cec0304SVijay Mahadevan PetscErrorCode ierr; 1145fc418013SVijay Mahadevan moab::ErrorCode merr; 1146fc418013SVijay Mahadevan DM_Moab *dmmoab; 11471cec0304SVijay Mahadevan 11481cec0304SVijay Mahadevan PetscFunctionBegin; 11491cec0304SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1150fc418013SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1151fc418013SVijay Mahadevan 11521cec0304SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1153fc418013SVijay Mahadevan 1154fc418013SVijay Mahadevan /* first get the global ID for the point */ 1155fc418013SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&point,1,&gid);MBERRNM(merr); 1156fc418013SVijay Mahadevan 1157fc418013SVijay Mahadevan /* get the dof value for the field */ 1158fc418013SVijay Mahadevan ierr = PetscSectionGetFieldDof(section, gid, field, dof);CHKERRQ(ierr); 11591cec0304SVijay Mahadevan PetscFunctionReturn(0); 11601cec0304SVijay Mahadevan } 11611cec0304SVijay Mahadevan 11621cec0304SVijay Mahadevan 11631cec0304SVijay Mahadevan #undef __FUNCT__ 11641cec0304SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldDofs" 11651cec0304SVijay Mahadevan PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 11661cec0304SVijay Mahadevan { 1167fc418013SVijay Mahadevan PetscInt i,gid; 11681cec0304SVijay Mahadevan PetscSection section; 11691cec0304SVijay Mahadevan PetscErrorCode ierr; 1170fc418013SVijay Mahadevan moab::ErrorCode merr; 1171fc418013SVijay Mahadevan DM_Moab *dmmoab; 11721cec0304SVijay Mahadevan 11731cec0304SVijay Mahadevan PetscFunctionBegin; 11741cec0304SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1175324f1edfSVijay Mahadevan PetscValidPointer(points,2); 1176fc418013SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1177fc418013SVijay Mahadevan 11781cec0304SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 11791cec0304SVijay Mahadevan if (!dof) { 11808d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr); 11811cec0304SVijay Mahadevan } 1182fc418013SVijay Mahadevan 1183fc418013SVijay Mahadevan /* first get the local indices */ 1184fc418013SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1185fc418013SVijay Mahadevan 11861cec0304SVijay Mahadevan for (i=0; i<npoints; ++i) { 1187fc418013SVijay Mahadevan gid=dof[i]; 1188fc418013SVijay Mahadevan ierr = PetscSectionGetFieldDof(section, gid, field, &dof[i]);CHKERRQ(ierr); 11891cec0304SVijay Mahadevan } 11901cec0304SVijay Mahadevan PetscFunctionReturn(0); 11911cec0304SVijay Mahadevan } 11921cec0304SVijay Mahadevan 11931cec0304SVijay Mahadevan 1194fc418013SVijay Mahadevan #undef __FUNCT__ 1195212ad6d1SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldDofsLocal" 1196212ad6d1SVijay Mahadevan PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 1197eb9d2429SVijay Mahadevan { 1198eb9d2429SVijay Mahadevan PetscInt i,offset; 1199eb9d2429SVijay Mahadevan PetscErrorCode ierr; 1200eb9d2429SVijay Mahadevan DM_Moab *dmmoab; 1201eb9d2429SVijay Mahadevan 1202eb9d2429SVijay Mahadevan PetscFunctionBegin; 1203eb9d2429SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1204324f1edfSVijay Mahadevan PetscValidPointer(points,2); 1205eb9d2429SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1206eb9d2429SVijay Mahadevan 1207eb9d2429SVijay Mahadevan if (!dof) { 12088d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr); 1209eb9d2429SVijay Mahadevan } 1210eb9d2429SVijay Mahadevan 1211eb9d2429SVijay Mahadevan if (dmmoab->bs > 1) { 1212eb9d2429SVijay Mahadevan for (i=0; i<npoints; ++i) 1213eb9d2429SVijay Mahadevan dof[i] = (points[i]-1)*dmmoab->bs+field; 1214eb9d2429SVijay Mahadevan } 1215eb9d2429SVijay Mahadevan else { 1216eb9d2429SVijay Mahadevan offset = field*dmmoab->n; /* assume all fields have equal distribution */ 1217eb9d2429SVijay Mahadevan for (i=0; i<npoints; ++i) 1218eb9d2429SVijay Mahadevan dof[i] = offset+points[i]-1; 1219eb9d2429SVijay Mahadevan } 1220eb9d2429SVijay Mahadevan PetscFunctionReturn(0); 1221eb9d2429SVijay Mahadevan } 1222eb9d2429SVijay Mahadevan 1223eb9d2429SVijay Mahadevan 1224eb9d2429SVijay Mahadevan #undef __FUNCT__ 1225212ad6d1SVijay Mahadevan #define __FUNCT__ "DMMoabGetDofs" 1226212ad6d1SVijay Mahadevan PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1227212ad6d1SVijay Mahadevan { 1228212ad6d1SVijay Mahadevan PetscInt i,f,gid; 1229212ad6d1SVijay Mahadevan PetscSection section; 1230212ad6d1SVijay Mahadevan PetscErrorCode ierr; 1231212ad6d1SVijay Mahadevan moab::ErrorCode merr; 1232212ad6d1SVijay Mahadevan DM_Moab *dmmoab; 1233212ad6d1SVijay Mahadevan 1234212ad6d1SVijay Mahadevan PetscFunctionBegin; 1235212ad6d1SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1236324f1edfSVijay Mahadevan PetscValidPointer(points,2); 1237212ad6d1SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1238212ad6d1SVijay Mahadevan 1239212ad6d1SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1240212ad6d1SVijay Mahadevan if (!dof) { 12418d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr); 1242212ad6d1SVijay Mahadevan } 1243212ad6d1SVijay Mahadevan 1244212ad6d1SVijay Mahadevan /* first get the local indices */ 1245212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1246212ad6d1SVijay Mahadevan 1247212ad6d1SVijay Mahadevan for (i=0; i<npoints; ++i) { 1248212ad6d1SVijay Mahadevan gid=dof[i]; 1249212ad6d1SVijay Mahadevan for (f=0; f<dmmoab->nfields; ++f) { 1250212ad6d1SVijay Mahadevan ierr = PetscSectionGetFieldDof(section, gid, f, &dof[i*dmmoab->nfields+f]);CHKERRQ(ierr); 1251212ad6d1SVijay Mahadevan } 1252212ad6d1SVijay Mahadevan } 1253212ad6d1SVijay Mahadevan PetscFunctionReturn(0); 1254212ad6d1SVijay Mahadevan } 1255212ad6d1SVijay Mahadevan 1256212ad6d1SVijay Mahadevan 1257212ad6d1SVijay Mahadevan #undef __FUNCT__ 1258212ad6d1SVijay Mahadevan #define __FUNCT__ "DMMoabGetDofsLocal" 1259212ad6d1SVijay Mahadevan PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1260212ad6d1SVijay Mahadevan { 1261212ad6d1SVijay Mahadevan PetscInt i,f,offset; 1262212ad6d1SVijay Mahadevan PetscErrorCode ierr; 1263212ad6d1SVijay Mahadevan DM_Moab *dmmoab; 1264212ad6d1SVijay Mahadevan 1265212ad6d1SVijay Mahadevan PetscFunctionBegin; 1266212ad6d1SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1267324f1edfSVijay Mahadevan PetscValidPointer(points,2); 1268212ad6d1SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1269212ad6d1SVijay Mahadevan 1270212ad6d1SVijay Mahadevan if (!dof) { 12718d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr); 1272212ad6d1SVijay Mahadevan } 1273212ad6d1SVijay Mahadevan 1274212ad6d1SVijay Mahadevan if (dmmoab->bs > 1) { 1275212ad6d1SVijay Mahadevan for (f=0; f<dmmoab->nfields; ++f) 1276212ad6d1SVijay Mahadevan for (i=0; i<npoints; ++i) 1277212ad6d1SVijay Mahadevan dof[i*dmmoab->nfields+f] = (points[i]-1)*dmmoab->bs+f; 1278212ad6d1SVijay Mahadevan } 1279212ad6d1SVijay Mahadevan else { 1280212ad6d1SVijay Mahadevan for (f=0; f<dmmoab->nfields; ++f) { 1281212ad6d1SVijay Mahadevan offset = f*dmmoab->n; /* assume all fields have equal distribution - say all vertex based */ 1282212ad6d1SVijay Mahadevan for (i=0; i<npoints; ++i) 1283212ad6d1SVijay Mahadevan dof[i*dmmoab->nfields+f] = offset+points[i]-1; 1284212ad6d1SVijay Mahadevan } 1285212ad6d1SVijay Mahadevan } 1286212ad6d1SVijay Mahadevan PetscFunctionReturn(0); 1287212ad6d1SVijay Mahadevan } 1288212ad6d1SVijay Mahadevan 1289212ad6d1SVijay Mahadevan 1290212ad6d1SVijay Mahadevan #undef __FUNCT__ 1291212ad6d1SVijay Mahadevan #define __FUNCT__ "DMMoabGetDofsBlocked" 1292212ad6d1SVijay Mahadevan PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1293212ad6d1SVijay Mahadevan { 1294212ad6d1SVijay Mahadevan PetscInt i,gid,dofindx; 1295212ad6d1SVijay Mahadevan PetscSection section; 1296212ad6d1SVijay Mahadevan PetscErrorCode ierr; 1297212ad6d1SVijay Mahadevan moab::ErrorCode merr; 1298212ad6d1SVijay Mahadevan DM_Moab *dmmoab; 1299212ad6d1SVijay Mahadevan 1300212ad6d1SVijay Mahadevan PetscFunctionBegin; 1301212ad6d1SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1302324f1edfSVijay Mahadevan PetscValidPointer(points,2); 1303212ad6d1SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1304212ad6d1SVijay Mahadevan 1305212ad6d1SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1306212ad6d1SVijay Mahadevan if (!dof) { 13078d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr); 1308212ad6d1SVijay Mahadevan } 1309212ad6d1SVijay Mahadevan 1310212ad6d1SVijay Mahadevan /* first get the local indices */ 1311212ad6d1SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr); 1312212ad6d1SVijay Mahadevan 1313212ad6d1SVijay Mahadevan for (i=0; i<npoints; ++i) { 1314212ad6d1SVijay Mahadevan gid=dof[i]; 1315212ad6d1SVijay Mahadevan ierr = PetscSectionGetFieldDof(section, gid, 0, &dofindx);CHKERRQ(ierr); 1316212ad6d1SVijay Mahadevan if (dmmoab->bs > 1) dof[i]=dofindx/dmmoab->bs; 1317212ad6d1SVijay Mahadevan else dof[i]=dofindx; 1318212ad6d1SVijay Mahadevan } 1319212ad6d1SVijay Mahadevan PetscFunctionReturn(0); 1320212ad6d1SVijay Mahadevan } 1321212ad6d1SVijay Mahadevan 1322212ad6d1SVijay Mahadevan 1323212ad6d1SVijay Mahadevan #undef __FUNCT__ 1324212ad6d1SVijay Mahadevan #define __FUNCT__ "DMMoabGetDofsBlockedLocal" 1325212ad6d1SVijay Mahadevan PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 1326212ad6d1SVijay Mahadevan { 1327212ad6d1SVijay Mahadevan PetscInt i; 1328212ad6d1SVijay Mahadevan PetscErrorCode ierr; 1329212ad6d1SVijay Mahadevan 1330212ad6d1SVijay Mahadevan PetscFunctionBegin; 1331212ad6d1SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1332324f1edfSVijay Mahadevan PetscValidPointer(points,2); 1333212ad6d1SVijay Mahadevan 1334212ad6d1SVijay Mahadevan if (!dof) { 13358d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr); 1336212ad6d1SVijay Mahadevan } 1337212ad6d1SVijay Mahadevan 1338212ad6d1SVijay Mahadevan for (i=0; i<npoints; ++i) 1339212ad6d1SVijay Mahadevan dof[i] = points[i]-1; 1340212ad6d1SVijay Mahadevan PetscFunctionReturn(0); 1341212ad6d1SVijay Mahadevan } 1342212ad6d1SVijay Mahadevan 13438d8d51c8SVijay Mahadevan 13448d8d51c8SVijay Mahadevan #undef __FUNCT__ 13458d8d51c8SVijay Mahadevan #define __FUNCT__ "DMMoabGetVertexDofsBlocked" 13468d8d51c8SVijay Mahadevan PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm,PetscInt** dof) 13478d8d51c8SVijay Mahadevan { 13488d8d51c8SVijay Mahadevan PetscInt i,gid; 13498d8d51c8SVijay Mahadevan DM_Moab *dmmoab; 13508d8d51c8SVijay Mahadevan PetscSection section; 13518d8d51c8SVijay Mahadevan PetscErrorCode ierr; 13528d8d51c8SVijay Mahadevan moab::ErrorCode merr; 13538d8d51c8SVijay Mahadevan 13548d8d51c8SVijay Mahadevan PetscFunctionBegin; 13558d8d51c8SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 13568d8d51c8SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 13578d8d51c8SVijay Mahadevan 13588d8d51c8SVijay Mahadevan *dof = dmmoab->gsindices; 13598d8d51c8SVijay Mahadevan 13608d8d51c8SVijay Mahadevan if (false) { 13618d8d51c8SVijay Mahadevan if (!dof) { 13628d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*(dmmoab->nloc+dmmoab->nghost), dof);CHKERRQ(ierr); 13638d8d51c8SVijay Mahadevan } 13648d8d51c8SVijay Mahadevan 13658d8d51c8SVijay Mahadevan /* first get the local indices */ 13668d8d51c8SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vlocal,*dof);MBERRNM(merr); 13678d8d51c8SVijay Mahadevan 13688d8d51c8SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 13698d8d51c8SVijay Mahadevan 13708d8d51c8SVijay Mahadevan /* Compute function over the locally owned part of the grid */ 13718d8d51c8SVijay Mahadevan for(i=0; i<dmmoab->nloc+dmmoab->nghost; i++) { 13728d8d51c8SVijay Mahadevan gid=(*dof)[i]; 13738d8d51c8SVijay Mahadevan ierr = PetscSectionGetFieldDof(section, gid, 0, &(*dof)[i]);CHKERRQ(ierr); 13748d8d51c8SVijay Mahadevan } 13758d8d51c8SVijay Mahadevan } 13768d8d51c8SVijay Mahadevan PetscFunctionReturn(0); 13778d8d51c8SVijay Mahadevan } 13788d8d51c8SVijay Mahadevan 13798d8d51c8SVijay Mahadevan 13808d8d51c8SVijay Mahadevan #undef __FUNCT__ 13818d8d51c8SVijay Mahadevan #define __FUNCT__ "DMMoabGetVertexDofsBlockedLocal" 13828d8d51c8SVijay Mahadevan PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm,PetscInt** dof) 13838d8d51c8SVijay Mahadevan { 13848d8d51c8SVijay Mahadevan PetscInt i; 13858d8d51c8SVijay Mahadevan DM_Moab *dmmoab; 13868d8d51c8SVijay Mahadevan PetscErrorCode ierr; 13878d8d51c8SVijay Mahadevan 13888d8d51c8SVijay Mahadevan PetscFunctionBegin; 13898d8d51c8SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 13908d8d51c8SVijay Mahadevan PetscValidPointer(dof,2); 13918d8d51c8SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 13928d8d51c8SVijay Mahadevan 13938d8d51c8SVijay Mahadevan if (!(*dof)) { 13948d8d51c8SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*(dmmoab->nloc+dmmoab->nghost), dof);CHKERRQ(ierr); 13958d8d51c8SVijay Mahadevan } 13968d8d51c8SVijay Mahadevan 13978d8d51c8SVijay Mahadevan i=0; 13988d8d51c8SVijay Mahadevan /* Compute function over the locally owned part of the grid */ 13998d8d51c8SVijay Mahadevan for(moab::Range::iterator iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++,i++) { 14008d8d51c8SVijay Mahadevan (*dof)[i] = (*iter)-1; 14018d8d51c8SVijay Mahadevan } 14028d8d51c8SVijay Mahadevan for(moab::Range::iterator iter = dmmoab->vghost->begin(); iter != dmmoab->vghost->end(); iter++,i++) { 14038d8d51c8SVijay Mahadevan (*dof)[i] = (*iter)-1; 14048d8d51c8SVijay Mahadevan } 14058d8d51c8SVijay Mahadevan PetscFunctionReturn(0); 14068d8d51c8SVijay Mahadevan } 14078d8d51c8SVijay Mahadevan 14088d8d51c8SVijay Mahadevan 1409212ad6d1SVijay Mahadevan #undef __FUNCT__ 1410fc418013SVijay Mahadevan #define __FUNCT__ "DMMoab_GetWriteOptions_Private" 1411fc418013SVijay Mahadevan PetscErrorCode DMMoab_GetWriteOptions_Private(PetscInt fsetid, PetscInt numproc, PetscInt dim, MoabWriteMode mode, PetscInt dbglevel, const char* extra_opts, const char** write_opts) 1412fc418013SVijay Mahadevan { 1413fc418013SVijay Mahadevan std::ostringstream str; 1414fc418013SVijay Mahadevan 1415fc418013SVijay Mahadevan PetscFunctionBegin; 1416fc418013SVijay Mahadevan 1417fc418013SVijay Mahadevan // do parallel read unless only one processor 1418fc418013SVijay Mahadevan if (numproc > 1) { 1419fc418013SVijay Mahadevan str << "PARALLEL=" << mode << ";"; 1420fc418013SVijay Mahadevan if (fsetid>=0) str << "PARALLEL_COMM=" << fsetid << ";"; 1421fc418013SVijay Mahadevan } 1422fc418013SVijay Mahadevan 1423fc418013SVijay Mahadevan if (dbglevel) 1424fc418013SVijay Mahadevan str << "CPUTIME;DEBUG_IO=" << dbglevel << ";"; 1425fc418013SVijay Mahadevan 1426fc418013SVijay Mahadevan if (extra_opts) 1427fc418013SVijay Mahadevan str << extra_opts ; 1428fc418013SVijay Mahadevan 1429fc418013SVijay Mahadevan *write_opts = str.str().c_str(); 1430fc418013SVijay Mahadevan PetscFunctionReturn(0); 1431fc418013SVijay Mahadevan } 1432fc418013SVijay Mahadevan 1433fc418013SVijay Mahadevan 1434fc418013SVijay Mahadevan #undef __FUNCT__ 1435fc418013SVijay Mahadevan #define __FUNCT__ "DMMoabOutput" 1436fc418013SVijay Mahadevan PetscErrorCode DMMoabOutput(DM dm,const char* filename,const char* usrwriteopts) 1437fc418013SVijay Mahadevan { 1438fc418013SVijay Mahadevan DM_Moab *dmmoab; 1439fc418013SVijay Mahadevan PetscInt dbglevel=0; 1440fc418013SVijay Mahadevan const char *writeopts; 1441fc418013SVijay Mahadevan 1442fc418013SVijay Mahadevan PetscErrorCode ierr; 1443fc418013SVijay Mahadevan moab::ErrorCode merr; 1444fc418013SVijay Mahadevan 1445fc418013SVijay Mahadevan PetscFunctionBegin; 1446fc418013SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1447fc418013SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1448fc418013SVijay Mahadevan 1449fc418013SVijay Mahadevan PetscBarrier((PetscObject)dm); 1450fc418013SVijay Mahadevan 1451fc418013SVijay Mahadevan /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */ 1452fc418013SVijay Mahadevan ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab"); 1453fc418013SVijay Mahadevan ierr = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr); 1454fc418013SVijay Mahadevan ierr = PetscOptionsEnd(); 1455fc418013SVijay Mahadevan 1456fc418013SVijay Mahadevan /* add mesh loading options specific to the DM */ 1457fc418013SVijay Mahadevan ierr = DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, MOAB_PARWOPTS_WRITE_PART, dbglevel, usrwriteopts, &writeopts);CHKERRQ(ierr); 1458fc418013SVijay Mahadevan PetscInfo2(dm, "Writing file %s with options: %s\n",filename,writeopts); 1459fc418013SVijay Mahadevan 1460fc418013SVijay Mahadevan /* output file, using parallel write */ 1461fc418013SVijay Mahadevan merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1);MBERRVM(dmmoab->mbiface,"Writing output of DMMoab failed.",merr); 1462fc418013SVijay Mahadevan PetscFunctionReturn(0); 1463fc418013SVijay Mahadevan } 1464fc418013SVijay Mahadevan 1465