151d15aeeSVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I "petscdm.h" I*/ 251d15aeeSVijay Mahadevan #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/ 351d15aeeSVijay Mahadevan 451d15aeeSVijay Mahadevan #include <petscdmmoab.h> 551d15aeeSVijay Mahadevan #include <MBTagConventions.hpp> 651d15aeeSVijay Mahadevan #include <moab/ReadUtilIface.hpp> 751d15aeeSVijay Mahadevan #include <moab/ScdInterface.hpp> 851d15aeeSVijay Mahadevan #include <moab/CN.hpp> 951d15aeeSVijay Mahadevan 1051d15aeeSVijay Mahadevan 1151d15aeeSVijay Mahadevan #undef __FUNCT__ 1251d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabComputeDomainBounds_Private" 1351d15aeeSVijay Mahadevan PetscErrorCode DMMoabComputeDomainBounds_Private(moab::ParallelComm* pcomm, PetscInt dim, PetscInt neleglob, PetscInt *ise) 1451d15aeeSVijay Mahadevan { 1551d15aeeSVijay Mahadevan PetscInt size,rank; 1651d15aeeSVijay Mahadevan PetscInt fraction,remainder; 1751d15aeeSVijay Mahadevan PetscInt starts[3],sizes[3]; 1851d15aeeSVijay Mahadevan 1951d15aeeSVijay Mahadevan PetscFunctionBegin; 2051d15aeeSVijay Mahadevan if(dim<1 && dim>3) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The problem dimension is invalid: %D",dim); 2151d15aeeSVijay Mahadevan 2251d15aeeSVijay Mahadevan size=pcomm->size(); 2351d15aeeSVijay Mahadevan rank=pcomm->rank(); 2451d15aeeSVijay Mahadevan fraction=neleglob/size; /* partition only by the largest dimension */ 2551d15aeeSVijay Mahadevan remainder=neleglob%size; /* remainder after partition which gets evenly distributed by round-robin */ 2651d15aeeSVijay Mahadevan 2751d15aeeSVijay Mahadevan if(fraction==0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The leading dimension size should be greater than number of processors: %D < %D",neleglob,size); 2851d15aeeSVijay Mahadevan 2951d15aeeSVijay Mahadevan starts[0]=starts[1]=starts[2]=0; /* default dimensional element = 1 */ 3051d15aeeSVijay Mahadevan sizes[0]=sizes[1]=sizes[2]=neleglob; /* default dimensional element = 1 */ 3151d15aeeSVijay Mahadevan 3251d15aeeSVijay Mahadevan if(rank < remainder) { 3351d15aeeSVijay Mahadevan /* This process gets "fraction+1" elements */ 3451d15aeeSVijay Mahadevan sizes[dim-1] = (fraction + 1); 3551d15aeeSVijay Mahadevan starts[dim-1] = rank * (fraction+1); 3651d15aeeSVijay Mahadevan } else { 3751d15aeeSVijay Mahadevan /* This process gets "fraction" elements */ 3851d15aeeSVijay Mahadevan sizes[dim-1] = fraction; 3951d15aeeSVijay Mahadevan starts[dim-1] = (remainder*(fraction+1) + fraction*(rank-remainder)); 4051d15aeeSVijay Mahadevan } 4151d15aeeSVijay Mahadevan 4251d15aeeSVijay Mahadevan for(int i=dim-1; i>=0; --i) { 4351d15aeeSVijay Mahadevan ise[2*i]=starts[i];ise[2*i+1]=starts[i]+sizes[i]; 4451d15aeeSVijay Mahadevan } 4551d15aeeSVijay Mahadevan 4651d15aeeSVijay Mahadevan PetscFunctionReturn(0); 4751d15aeeSVijay Mahadevan } 4851d15aeeSVijay Mahadevan 4951d15aeeSVijay Mahadevan 5051d15aeeSVijay Mahadevan #undef __FUNCT__ 5151d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabCreateBoxMesh" 5251d15aeeSVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt nele, PetscInt nghost, DM *dm) 5351d15aeeSVijay Mahadevan { 5451d15aeeSVijay Mahadevan PetscErrorCode ierr; 5551d15aeeSVijay Mahadevan moab::ErrorCode merr; 56*a4717116SVijay Mahadevan PetscInt nprocs; 5751d15aeeSVijay Mahadevan DM_Moab *dmmoab; 5851d15aeeSVijay Mahadevan moab::Interface *mbiface; 5951d15aeeSVijay Mahadevan moab::ParallelComm *pcomm; 60*a4717116SVijay Mahadevan moab::ReadUtilIface* readMeshIface; 61*a4717116SVijay Mahadevan 6251d15aeeSVijay Mahadevan moab::Tag id_tag=PETSC_NULL; 63*a4717116SVijay Mahadevan moab::Range ownedvtx,ownedelms; 64*a4717116SVijay Mahadevan moab::EntityHandle vfirst,efirst; 65*a4717116SVijay Mahadevan std::vector<double*> vcoords; 66*a4717116SVijay Mahadevan moab::EntityHandle *connectivity = 0; 67*a4717116SVijay Mahadevan std::vector<int> subent_conn; 6851d15aeeSVijay Mahadevan moab::EntityType etype; 6951d15aeeSVijay Mahadevan PetscInt ise[6]; 7051d15aeeSVijay Mahadevan PetscReal xse[6]; 7151d15aeeSVijay Mahadevan 72*a4717116SVijay Mahadevan const PetscInt npts=nele+1; /* Number of points in every dimension */ 73*a4717116SVijay Mahadevan PetscInt vpere,locnele,locnpts; /* Number of verts/element, vertices, elements owned by this process */ 7451d15aeeSVijay Mahadevan 7551d15aeeSVijay Mahadevan PetscFunctionBegin; 76*a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 7751d15aeeSVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 7851d15aeeSVijay Mahadevan 79*a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 8051d15aeeSVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 8151d15aeeSVijay Mahadevan mbiface = dmmoab->mbiface; 8251d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm; 8351d15aeeSVijay Mahadevan id_tag = dmmoab->ltog_tag; 8451d15aeeSVijay Mahadevan nprocs = pcomm->size(); 8551d15aeeSVijay Mahadevan 86*a4717116SVijay Mahadevan /* do some error checking */ 87*a4717116SVijay Mahadevan if(pow(npts,dim) < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2"); 88*a4717116SVijay Mahadevan if(nprocs > pow(nele,dim)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of processors must be less than or equal to number of elements"); 8951d15aeeSVijay Mahadevan 90*a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */ 91*a4717116SVijay Mahadevan merr = mbiface->query_interface(readMeshIface);MBERRNM(merr); 9251d15aeeSVijay Mahadevan 9351d15aeeSVijay Mahadevan ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr); 94*a4717116SVijay Mahadevan 95*a4717116SVijay Mahadevan /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */ 9651d15aeeSVijay Mahadevan ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr); 9751d15aeeSVijay Mahadevan 98*a4717116SVijay Mahadevan /* set some variables based on dimension */ 9951d15aeeSVijay Mahadevan switch(dim) { 10051d15aeeSVijay Mahadevan case 1: 10151d15aeeSVijay Mahadevan vpere = 2; 102*a4717116SVijay Mahadevan locnele = (ise[1]-ise[0]); 103*a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1); 10451d15aeeSVijay Mahadevan etype = moab::MBEDGE; 10551d15aeeSVijay Mahadevan break; 10651d15aeeSVijay Mahadevan case 2: 10751d15aeeSVijay Mahadevan vpere = 4; 108*a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2]); 109*a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1); 11051d15aeeSVijay Mahadevan etype = moab::MBQUAD; 11151d15aeeSVijay Mahadevan break; 11251d15aeeSVijay Mahadevan case 3: 11351d15aeeSVijay Mahadevan vpere = 8; 114*a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]); 115*a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1); 11651d15aeeSVijay Mahadevan etype = moab::MBHEX; 11751d15aeeSVijay Mahadevan break; 11851d15aeeSVijay Mahadevan } 11951d15aeeSVijay Mahadevan 12051d15aeeSVijay Mahadevan /* we have a domain of size [1,1,1] - now compute local co-ordinate box */ 12151d15aeeSVijay Mahadevan ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr); 12251d15aeeSVijay Mahadevan for(int i=0; i<6; ++i) { 12351d15aeeSVijay Mahadevan xse[i]=(PetscReal)ise[i]/nele; 12451d15aeeSVijay Mahadevan } 12551d15aeeSVijay Mahadevan 126*a4717116SVijay Mahadevan /* Create vertexes and set the coodinate of each vertex */ 127*a4717116SVijay Mahadevan const int sequence_size = (locnele + vpere) + 1; 128*a4717116SVijay Mahadevan merr = readMeshIface->get_node_coords(3,locnpts,0,vfirst,vcoords,sequence_size);MBERRNM(merr); 12951d15aeeSVijay Mahadevan 130*a4717116SVijay Mahadevan /* Compute the co-ordinates of vertices and global IDs */ 131*a4717116SVijay Mahadevan std::vector<int> vgid(locnpts); 13251d15aeeSVijay Mahadevan int vcount=0; 13351d15aeeSVijay Mahadevan double hxyz=1.0/nele; 13451d15aeeSVijay Mahadevan for (int k = ise[4]; k <= ise[5]; k++) { 13551d15aeeSVijay Mahadevan for (int j = ise[2]; j <= ise[3]; j++) { 13651d15aeeSVijay Mahadevan for (int i = ise[0]; i <= ise[1]; i++, vcount++) { 13751d15aeeSVijay Mahadevan vcoords[0][vcount] = i*hxyz; 13851d15aeeSVijay Mahadevan vcoords[1][vcount] = j*hxyz; 13951d15aeeSVijay Mahadevan vcoords[2][vcount] = k*hxyz; 14051d15aeeSVijay Mahadevan vgid[vcount] = (k*nele+j)*nele+i; 14151d15aeeSVijay Mahadevan } 14251d15aeeSVijay Mahadevan } 14351d15aeeSVijay Mahadevan } 14451d15aeeSVijay Mahadevan 14551d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 146*a4717116SVijay Mahadevan merr = mbiface->add_entities (dmmoab->fileset, ownedvtx);MBERRNM(merr); 14751d15aeeSVijay Mahadevan 148*a4717116SVijay Mahadevan /* The global ID tag is applied to each owned 149*a4717116SVijay Mahadevan vertex. It acts as an global identifier which MOAB uses to 150*a4717116SVijay Mahadevan assemble the individual pieces of the mesh: 151*a4717116SVijay Mahadevan Set the global ID indices */ 15251d15aeeSVijay Mahadevan merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr); 15351d15aeeSVijay Mahadevan 154*a4717116SVijay Mahadevan /* Create elements between mesh points using the ReadUtilInterface 155*a4717116SVijay Mahadevan get the reference to element connectivities for all local elements from the ReadUtilInterface */ 156*a4717116SVijay Mahadevan merr = readMeshIface->get_element_connect (locnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr); 15751d15aeeSVijay Mahadevan 158*a4717116SVijay Mahadevan /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */ 159*a4717116SVijay Mahadevan vfirst-=vgid[0]; 16051d15aeeSVijay Mahadevan 161*a4717116SVijay Mahadevan /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k): 162*a4717116SVijay Mahadevan and then set the connectivity for each element appropriately */ 16351d15aeeSVijay Mahadevan int ecount=0; 164*a4717116SVijay Mahadevan subent_conn.resize(vpere); 16551d15aeeSVijay Mahadevan for (int k = ise[4]; k < std::max(ise[5],1); k++) { 16651d15aeeSVijay Mahadevan for (int j = ise[2]; j < std::max(ise[3],1); j++) { 16751d15aeeSVijay Mahadevan for (int i = ise[0]; i < std::max(ise[1],1); i++,ecount++) { 16851d15aeeSVijay Mahadevan const int offset = ecount*vpere; 16951d15aeeSVijay Mahadevan moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data()); 17051d15aeeSVijay Mahadevan 17151d15aeeSVijay Mahadevan switch(dim) { 17251d15aeeSVijay Mahadevan case 1: 17351d15aeeSVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i; 17451d15aeeSVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1); 175*a4717116SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "ELEMENT[%D,%D,%D]: CONNECTIVITY = %D, %D\n", i,j,k,connectivity[offset+subent_conn[0]], connectivity[offset+subent_conn[1]]); 17651d15aeeSVijay Mahadevan break; 17751d15aeeSVijay Mahadevan case 2: 17851d15aeeSVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1); 17951d15aeeSVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1); 18051d15aeeSVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1); 18151d15aeeSVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1); 182*a4717116SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "ELEMENT[%D,%D,%D]: CONNECTIVITY = %D, %D, %D, %D\n", i,j,k,connectivity[offset+subent_conn[0]], connectivity[offset+subent_conn[1]], connectivity[offset+subent_conn[2]], connectivity[offset+subent_conn[3]]); 18351d15aeeSVijay Mahadevan break; 18451d15aeeSVijay Mahadevan case 3: 18551d15aeeSVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 18651d15aeeSVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k); 18751d15aeeSVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 18851d15aeeSVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k); 18951d15aeeSVijay Mahadevan connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1)); 19051d15aeeSVijay Mahadevan connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 19151d15aeeSVijay Mahadevan connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1)); 19251d15aeeSVijay Mahadevan connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 19351d15aeeSVijay Mahadevan break; 19451d15aeeSVijay Mahadevan } 19551d15aeeSVijay Mahadevan } 19651d15aeeSVijay Mahadevan } 19751d15aeeSVijay Mahadevan } 198*a4717116SVijay Mahadevan merr = readMeshIface->update_adjacencies(efirst,locnele,vpere,connectivity);MBERRNM(merr); 19951d15aeeSVijay Mahadevan 200*a4717116SVijay Mahadevan /* 2. Get the vertices and hexes from moab and check their numbers against I*J*K and (I-1)*(J-1)*(K-1), resp. 201*a4717116SVijay Mahadevan first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */ 20251d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr); 203*a4717116SVijay Mahadevan merr = mbiface->add_entities (dmmoab->fileset, ownedelms);MBERRNM(merr); 20451d15aeeSVijay Mahadevan 205*a4717116SVijay Mahadevan // merr = mbiface->unite_meshset(dmmoab->fileset, 0);MBERRV(mbiface,merr); 206*a4717116SVijay Mahadevan 207*a4717116SVijay Mahadevan if (locnele != (int) ownedelms.size()) 208*a4717116SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele,ownedelms.size()); 209*a4717116SVijay Mahadevan else if(locnpts != (int) ownedvtx.size()) 210*a4717116SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts,ownedvtx.size()); 21151d15aeeSVijay Mahadevan else 212*a4717116SVijay Mahadevan PetscInfo2(*dm, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size()); 21351d15aeeSVijay Mahadevan 214*a4717116SVijay Mahadevan /* check the handles */ 215*a4717116SVijay Mahadevan merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr); 21651d15aeeSVijay Mahadevan 217*a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */ 218*a4717116SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr); 219*a4717116SVijay Mahadevan merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,&id_tag);MBERRV(mbiface,merr); 22051d15aeeSVijay Mahadevan 221*a4717116SVijay Mahadevan /* Reassign global IDs on all entities. */ 222*a4717116SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true);MBERRNM(merr); 223*a4717116SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,nghost,0,true);MBERRV(mbiface,merr); 22451d15aeeSVijay Mahadevan 225*a4717116SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 226*a4717116SVijay Mahadevan on all of the ghost vertexes */ 227*a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr); 228*a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr); 22951d15aeeSVijay Mahadevan 230*a4717116SVijay Mahadevan // ierr = DMMoabSetLocalVertices(*dm, &ownedvtx);CHKERRQ(ierr); 231*a4717116SVijay Mahadevan // ierr = DMMoabSetLocalElements(*dm, &ownedelms);CHKERRQ(ierr); 23251d15aeeSVijay Mahadevan 23351d15aeeSVijay Mahadevan merr = mbiface->set_dimension(dim);MBERRNM(merr); 23451d15aeeSVijay Mahadevan 23551d15aeeSVijay Mahadevan std::stringstream sstr; 236*a4717116SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 23751d15aeeSVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 23851d15aeeSVijay Mahadevan PetscFunctionReturn(0); 23951d15aeeSVijay Mahadevan } 24051d15aeeSVijay Mahadevan 24151d15aeeSVijay Mahadevan 242*a4717116SVijay Mahadevan #undef __FUNCT__ 243*a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private" 244*a4717116SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* extra_opts, const char* read_opts) 24551d15aeeSVijay Mahadevan { 246*a4717116SVijay Mahadevan std::ostringstream str; 24751d15aeeSVijay Mahadevan 248*a4717116SVijay Mahadevan PetscFunctionBegin; 249*a4717116SVijay Mahadevan // do parallel read unless only one processor 250*a4717116SVijay Mahadevan if (numproc > 1) { 251*a4717116SVijay Mahadevan str << "PARALLEL=" << mode << ";PARTITION=PARALLEL_PARTITION;"; 252*a4717116SVijay Mahadevan str << "PARTITION_DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;"; 253*a4717116SVijay Mahadevan if (by_rank) 254*a4717116SVijay Mahadevan str << "PARTITION_BY_RANK;"; 255*a4717116SVijay Mahadevan } 25651d15aeeSVijay Mahadevan 257*a4717116SVijay Mahadevan if (extra_opts) 258*a4717116SVijay Mahadevan str << extra_opts << ";"; 25951d15aeeSVijay Mahadevan 260*a4717116SVijay Mahadevan if (dbglevel) 261*a4717116SVijay Mahadevan str << "DEBUG_IO=" << dbglevel << ";CPUTIME;"; 26251d15aeeSVijay Mahadevan 263*a4717116SVijay Mahadevan read_opts = str.str().c_str(); 26451d15aeeSVijay Mahadevan PetscFunctionReturn(0); 26551d15aeeSVijay Mahadevan } 26651d15aeeSVijay Mahadevan 26751d15aeeSVijay Mahadevan 26851d15aeeSVijay Mahadevan #undef __FUNCT__ 269*a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile" 270*a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm) 27151d15aeeSVijay Mahadevan { 272*a4717116SVijay Mahadevan moab::ErrorCode merr; 273*a4717116SVijay Mahadevan PetscInt nprocs; 274*a4717116SVijay Mahadevan DM_Moab *dmmoab; 275*a4717116SVijay Mahadevan moab::Interface *mbiface; 276*a4717116SVijay Mahadevan moab::ParallelComm *pcomm; 277*a4717116SVijay Mahadevan moab::Range verts,elems; 278*a4717116SVijay Mahadevan const char* readopts=0; 279*a4717116SVijay Mahadevan PetscErrorCode ierr; 28051d15aeeSVijay Mahadevan 28151d15aeeSVijay Mahadevan PetscFunctionBegin; 282*a4717116SVijay Mahadevan PetscValidPointer(dm,5); 28351d15aeeSVijay Mahadevan 284*a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 285*a4717116SVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 28651d15aeeSVijay Mahadevan 287*a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 288*a4717116SVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 289*a4717116SVijay Mahadevan mbiface = dmmoab->mbiface; 290*a4717116SVijay Mahadevan pcomm = dmmoab->pcomm; 291*a4717116SVijay Mahadevan nprocs = pcomm->size(); 292*a4717116SVijay Mahadevan 293*a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */ 294*a4717116SVijay Mahadevan ierr = DMMoab_GetReadOptions_Private(PETSC_TRUE, nprocs, dim, MOAB_PARROPTS_READ_PART, 0, usrreadopts, readopts);CHKERRQ(ierr); 295*a4717116SVijay Mahadevan 296*a4717116SVijay Mahadevan /* Load the mesh from a file. */ 297*a4717116SVijay Mahadevan merr = mbiface->load_file(filename, &dmmoab->fileset, ""/*readopts*/);MBERRVM(mbiface,"Reading MOAB file failed.", merr); 298*a4717116SVijay Mahadevan 299*a4717116SVijay Mahadevan merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 300*a4717116SVijay Mahadevan 301*a4717116SVijay Mahadevan /* load the local vertices */ 302*a4717116SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 303*a4717116SVijay Mahadevan ierr = DMMoabSetLocalVertices(*dm, &verts);CHKERRQ(ierr); 304*a4717116SVijay Mahadevan 305*a4717116SVijay Mahadevan /* load the local elements */ 306*a4717116SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 307*a4717116SVijay Mahadevan // ierr = DMMoabSetLocalElements(*dm, &elems);CHKERRQ(ierr); 308*a4717116SVijay Mahadevan 309*a4717116SVijay Mahadevan merr = mbiface->set_dimension(dim);MBERRNM(merr); 310*a4717116SVijay Mahadevan 311*a4717116SVijay Mahadevan PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size()); 31251d15aeeSVijay Mahadevan PetscFunctionReturn(0); 31351d15aeeSVijay Mahadevan } 31451d15aeeSVijay Mahadevan 315