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 49e5136372SVijay Mahadevan static void set_structured_coordinates(PetscInt i, PetscInt j, PetscInt k, PetscReal hxyz, PetscInt vcount, std::vector<double*>& vcoords) 50e5136372SVijay Mahadevan { 51e5136372SVijay Mahadevan vcoords[0][vcount] = i*hxyz; 52e5136372SVijay Mahadevan vcoords[1][vcount] = j*hxyz; 53e5136372SVijay Mahadevan vcoords[2][vcount] = k*hxyz; 54e5136372SVijay Mahadevan } 55e5136372SVijay Mahadevan 56e5136372SVijay Mahadevan static void set_element_connectivity(PetscInt dim, moab::EntityType etype, PetscInt offset, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity) 57e5136372SVijay Mahadevan { 58e5136372SVijay Mahadevan std::vector<int> subent_conn(pow(2,dim)); /* only linear edge, quad, hex supported now */ 59e5136372SVijay Mahadevan 60e5136372SVijay Mahadevan moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data()); 61e5136372SVijay Mahadevan 62e5136372SVijay Mahadevan switch(dim) { 63e5136372SVijay Mahadevan case 1: 64e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i; 65e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1); 66e5136372SVijay Mahadevan break; 67e5136372SVijay Mahadevan case 2: 68e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1); 69e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1); 70e5136372SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1); 71e5136372SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1); 72e5136372SVijay Mahadevan break; 73e5136372SVijay Mahadevan case 3: 74e5136372SVijay Mahadevan default: 75e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 76e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k); 77e5136372SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 78e5136372SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k); 79e5136372SVijay Mahadevan connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1)); 80e5136372SVijay Mahadevan connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 81e5136372SVijay Mahadevan connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1)); 82e5136372SVijay Mahadevan connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 83e5136372SVijay Mahadevan break; 84e5136372SVijay Mahadevan } 85e5136372SVijay Mahadevan } 8651d15aeeSVijay Mahadevan 8751d15aeeSVijay Mahadevan #undef __FUNCT__ 8851d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabCreateBoxMesh" 8951d15aeeSVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt nele, PetscInt nghost, DM *dm) 9051d15aeeSVijay Mahadevan { 9151d15aeeSVijay Mahadevan PetscErrorCode ierr; 9251d15aeeSVijay Mahadevan moab::ErrorCode merr; 93e5136372SVijay Mahadevan PetscInt i,j,k,n,nprocs,rank; 9451d15aeeSVijay Mahadevan DM_Moab *dmmoab; 9551d15aeeSVijay Mahadevan moab::Interface *mbiface; 9651d15aeeSVijay Mahadevan moab::ParallelComm *pcomm; 97a4717116SVijay Mahadevan moab::ReadUtilIface* readMeshIface; 98a4717116SVijay Mahadevan 9951d15aeeSVijay Mahadevan moab::Tag id_tag=PETSC_NULL; 100a4717116SVijay Mahadevan moab::Range ownedvtx,ownedelms; 101a4717116SVijay Mahadevan moab::EntityHandle vfirst,efirst; 102a4717116SVijay Mahadevan std::vector<double*> vcoords; 103a4717116SVijay Mahadevan moab::EntityHandle *connectivity = 0; 10451d15aeeSVijay Mahadevan moab::EntityType etype; 105c8d5365dSVijay Mahadevan PetscInt ise[6]; 10651d15aeeSVijay Mahadevan PetscReal xse[6]; 10751d15aeeSVijay Mahadevan 108a4717116SVijay Mahadevan const PetscInt npts=nele+1; /* Number of points in every dimension */ 109e5136372SVijay Mahadevan PetscInt vpere,locnele,locnpts,ghnele,ghnpts; /* Number of verts/element, vertices, elements owned by this process */ 11051d15aeeSVijay Mahadevan 11151d15aeeSVijay Mahadevan PetscFunctionBegin; 112e5136372SVijay Mahadevan if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n"); 113e5136372SVijay Mahadevan 114c8d5365dSVijay Mahadevan ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr); 115e5136372SVijay Mahadevan /* total number of vertices in all dimensions */ 116e5136372SVijay Mahadevan n=pow(npts,dim); 117e5136372SVijay Mahadevan 118e5136372SVijay Mahadevan /* do some error checking */ 119e5136372SVijay Mahadevan if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n"); 120e5136372SVijay Mahadevan if(nprocs > n) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of processors must be less than or equal to number of elements.\n"); 121e5136372SVijay Mahadevan if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n"); 122e5136372SVijay Mahadevan 123a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 12451d15aeeSVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 12551d15aeeSVijay Mahadevan 126a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 12751d15aeeSVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 12851d15aeeSVijay Mahadevan mbiface = dmmoab->mbiface; 12951d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm; 13051d15aeeSVijay Mahadevan id_tag = dmmoab->ltog_tag; 13151d15aeeSVijay Mahadevan nprocs = pcomm->size(); 132e5136372SVijay Mahadevan rank = pcomm->rank(); 133c8d5365dSVijay Mahadevan dmmoab->dim = dim; 13451d15aeeSVijay Mahadevan 135a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */ 136a4717116SVijay Mahadevan merr = mbiface->query_interface(readMeshIface);MBERRNM(merr); 13751d15aeeSVijay Mahadevan 13851d15aeeSVijay Mahadevan ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr); 139a4717116SVijay Mahadevan 140a4717116SVijay Mahadevan /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */ 14151d15aeeSVijay Mahadevan ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr); 14251d15aeeSVijay Mahadevan 143a4717116SVijay Mahadevan /* set some variables based on dimension */ 14451d15aeeSVijay Mahadevan switch(dim) { 14551d15aeeSVijay Mahadevan case 1: 14651d15aeeSVijay Mahadevan vpere = 2; 147a4717116SVijay Mahadevan locnele = (ise[1]-ise[0]); 148a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1); 149c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 ); 150c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0); 15151d15aeeSVijay Mahadevan etype = moab::MBEDGE; 15251d15aeeSVijay Mahadevan break; 15351d15aeeSVijay Mahadevan case 2: 15451d15aeeSVijay Mahadevan vpere = 4; 155a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2]); 156a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1); 157c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0); 158c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0); 15951d15aeeSVijay Mahadevan etype = moab::MBQUAD; 16051d15aeeSVijay Mahadevan break; 16151d15aeeSVijay Mahadevan case 3: 16251d15aeeSVijay Mahadevan vpere = 8; 163a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]); 164a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1); 165c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0); 166c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0); 16751d15aeeSVijay Mahadevan etype = moab::MBHEX; 16851d15aeeSVijay Mahadevan break; 16951d15aeeSVijay Mahadevan } 17051d15aeeSVijay Mahadevan 17151d15aeeSVijay Mahadevan /* we have a domain of size [1,1,1] - now compute local co-ordinate box */ 17251d15aeeSVijay Mahadevan ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr); 173c8d5365dSVijay Mahadevan for(i=0; i<6; ++i) { 17451d15aeeSVijay Mahadevan xse[i]=(PetscReal)ise[i]/nele; 17551d15aeeSVijay Mahadevan } 17651d15aeeSVijay Mahadevan 177a4717116SVijay Mahadevan /* Create vertexes and set the coodinate of each vertex */ 178c8d5365dSVijay Mahadevan merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr); 17951d15aeeSVijay Mahadevan 180a4717116SVijay Mahadevan /* Compute the co-ordinates of vertices and global IDs */ 181c8d5365dSVijay Mahadevan std::vector<int> vgid(locnpts+ghnpts); 18251d15aeeSVijay Mahadevan int vcount=0; 18351d15aeeSVijay Mahadevan double hxyz=1.0/nele; 184e5136372SVijay Mahadevan 185c8d5365dSVijay Mahadevan /* create all the owned vertices */ 186c8d5365dSVijay Mahadevan for (k = ise[4]; k <= ise[5]; k++) { 187c8d5365dSVijay Mahadevan for (j = ise[2]; j <= ise[3]; j++) { 188c8d5365dSVijay Mahadevan for (i = ise[0]; i <= ise[1]; i++, vcount++) { 189c8d5365dSVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 190c8d5365dSVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 191c8d5365dSVijay Mahadevan } 192c8d5365dSVijay Mahadevan } 193c8d5365dSVijay Mahadevan } 194c8d5365dSVijay Mahadevan 195e5136372SVijay Mahadevan /* create ghosted vertices requested by user - below the current plane */ 196e5136372SVijay Mahadevan if (ise[2*dim-2] > 0) { 197e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) { 198e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) { 199e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) { 200e5136372SVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 201e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 202e5136372SVijay Mahadevan } 203e5136372SVijay Mahadevan } 204e5136372SVijay Mahadevan } 205e5136372SVijay Mahadevan } 206e5136372SVijay Mahadevan 207e5136372SVijay Mahadevan /* create ghosted vertices requested by user - above the current plane */ 208e5136372SVijay Mahadevan if (ise[2*dim-1] < nele) { 209e5136372SVijay Mahadevan for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) { 210e5136372SVijay Mahadevan for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) { 211e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) { 212e5136372SVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 213e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 214e5136372SVijay Mahadevan } 21551d15aeeSVijay Mahadevan } 21651d15aeeSVijay Mahadevan } 21751d15aeeSVijay Mahadevan } 21851d15aeeSVijay Mahadevan 219c8d5365dSVijay Mahadevan if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount); 220c8d5365dSVijay Mahadevan 22151d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 222a4717116SVijay Mahadevan merr = mbiface->add_entities (dmmoab->fileset, ownedvtx);MBERRNM(merr); 22351d15aeeSVijay Mahadevan 224a4717116SVijay Mahadevan /* The global ID tag is applied to each owned 225a4717116SVijay Mahadevan vertex. It acts as an global identifier which MOAB uses to 226a4717116SVijay Mahadevan assemble the individual pieces of the mesh: 227a4717116SVijay Mahadevan Set the global ID indices */ 22851d15aeeSVijay Mahadevan merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr); 22951d15aeeSVijay Mahadevan 230a4717116SVijay Mahadevan /* Create elements between mesh points using the ReadUtilInterface 231a4717116SVijay Mahadevan get the reference to element connectivities for all local elements from the ReadUtilInterface */ 232e5136372SVijay Mahadevan merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr); 23351d15aeeSVijay Mahadevan 234a4717116SVijay Mahadevan /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */ 235c8d5365dSVijay Mahadevan // PetscPrintf(PETSC_COMM_SELF, "[%D] first local handle %D\n", rank, vfirst); 236e5136372SVijay Mahadevan vfirst-=vgid[0]-1; 23751d15aeeSVijay Mahadevan 238a4717116SVijay Mahadevan /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k): 239a4717116SVijay Mahadevan and then set the connectivity for each element appropriately */ 24051d15aeeSVijay Mahadevan int ecount=0; 24151d15aeeSVijay Mahadevan 242e5136372SVijay Mahadevan /* create ghosted elements requested by user - below the current plane */ 243e5136372SVijay Mahadevan if (ise[2*dim-2] >= nghost) { 244e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) { 245e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) { 246e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++, ecount++) { 247c8d5365dSVijay Mahadevan set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 24851d15aeeSVijay Mahadevan } 24951d15aeeSVijay Mahadevan } 25051d15aeeSVijay Mahadevan } 25151d15aeeSVijay Mahadevan } 252e5136372SVijay Mahadevan 253e5136372SVijay Mahadevan /* create owned elements requested by user */ 254e5136372SVijay Mahadevan for (k = ise[4]; k < std::max(ise[5],1); k++) { 255e5136372SVijay Mahadevan for (j = ise[2]; j < std::max(ise[3],1); j++) { 256e5136372SVijay Mahadevan for (i = ise[0]; i < std::max(ise[1],1); i++,ecount++) { 257c8d5365dSVijay Mahadevan set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 258e5136372SVijay Mahadevan } 259e5136372SVijay Mahadevan } 260e5136372SVijay Mahadevan } 261e5136372SVijay Mahadevan 262e5136372SVijay Mahadevan /* create ghosted elements requested by user - above the current plane */ 263e5136372SVijay Mahadevan if (ise[2*dim-1] <= nele-nghost) { 264e5136372SVijay Mahadevan for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) { 265e5136372SVijay Mahadevan for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) { 266e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++, ecount++) { 267c8d5365dSVijay Mahadevan set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 268e5136372SVijay Mahadevan } 269e5136372SVijay Mahadevan } 270e5136372SVijay Mahadevan } 271e5136372SVijay Mahadevan } 272e5136372SVijay Mahadevan 273e5136372SVijay Mahadevan merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr); 27451d15aeeSVijay Mahadevan 275a4717116SVijay 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. 276a4717116SVijay Mahadevan first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */ 27751d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr); 278a4717116SVijay Mahadevan merr = mbiface->add_entities (dmmoab->fileset, ownedelms);MBERRNM(merr); 279a4717116SVijay Mahadevan 280e5136372SVijay Mahadevan if (locnele+ghnele != (int) ownedelms.size()) 281e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size()); 282e5136372SVijay Mahadevan else if(locnpts+ghnpts != (int) ownedvtx.size()) 283e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size()); 28451d15aeeSVijay Mahadevan else 285a4717116SVijay Mahadevan PetscInfo2(*dm, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size()); 28651d15aeeSVijay Mahadevan 287a4717116SVijay Mahadevan /* check the handles */ 288a4717116SVijay Mahadevan merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr); 28951d15aeeSVijay Mahadevan 290c8d5365dSVijay Mahadevan #if 0 291e5136372SVijay Mahadevan std::stringstream sstr; 292e5136372SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 293e5136372SVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 294e5136372SVijay Mahadevan #endif 295e5136372SVijay Mahadevan 296a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */ 297a4717116SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr); 298a4717116SVijay Mahadevan merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,&id_tag);MBERRV(mbiface,merr); 29951d15aeeSVijay Mahadevan 300a4717116SVijay Mahadevan /* Reassign global IDs on all entities. */ 301e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 30251d15aeeSVijay Mahadevan 303c8d5365dSVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,1/*nghost*/,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr); 304c8d5365dSVijay Mahadevan 305a4717116SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 306a4717116SVijay Mahadevan on all of the ghost vertexes */ 307c8d5365dSVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 308c8d5365dSVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr); 309c8d5365dSVijay Mahadevan 310a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr); 311a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr); 31251d15aeeSVijay Mahadevan PetscFunctionReturn(0); 31351d15aeeSVijay Mahadevan } 31451d15aeeSVijay Mahadevan 31551d15aeeSVijay Mahadevan 316a4717116SVijay Mahadevan #undef __FUNCT__ 317a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private" 3187023aa44SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* extra_opts, const char** read_opts) 31951d15aeeSVijay Mahadevan { 320a4717116SVijay Mahadevan std::ostringstream str; 32151d15aeeSVijay Mahadevan 322a4717116SVijay Mahadevan PetscFunctionBegin; 323*e23c60ebSVijay Mahadevan /* do parallel read unless using only one processor */ 324a4717116SVijay Mahadevan if (numproc > 1) { 3257023aa44SVijay Mahadevan str << "PARALLEL=" << mode << ";PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"; 3267023aa44SVijay Mahadevan str << "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;"; 327a4717116SVijay Mahadevan if (by_rank) 328a4717116SVijay Mahadevan str << "PARTITION_BY_RANK;"; 329a4717116SVijay Mahadevan } 33051d15aeeSVijay Mahadevan 331c8d5365dSVijay Mahadevan if (dbglevel) { 332c8d5365dSVijay Mahadevan str << "CPUTIME;DEBUG_IO=" << dbglevel << ";"; 333c8d5365dSVijay Mahadevan if (numproc>1) 334c8d5365dSVijay Mahadevan str << "DEBUG_PIO=" << dbglevel << ";"; 335c8d5365dSVijay Mahadevan } 33651d15aeeSVijay Mahadevan 337c8d5365dSVijay Mahadevan if (extra_opts) 338c8d5365dSVijay Mahadevan str << extra_opts; 33951d15aeeSVijay Mahadevan 3407023aa44SVijay Mahadevan *read_opts = str.str().c_str(); 34151d15aeeSVijay Mahadevan PetscFunctionReturn(0); 34251d15aeeSVijay Mahadevan } 34351d15aeeSVijay Mahadevan 34451d15aeeSVijay Mahadevan 34551d15aeeSVijay Mahadevan #undef __FUNCT__ 346a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile" 347a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm) 34851d15aeeSVijay Mahadevan { 349a4717116SVijay Mahadevan moab::ErrorCode merr; 3507023aa44SVijay Mahadevan PetscInt nprocs,dbglevel=0; 3517023aa44SVijay Mahadevan PetscBool part_by_rank=PETSC_FALSE; 352a4717116SVijay Mahadevan DM_Moab *dmmoab; 353a4717116SVijay Mahadevan moab::Interface *mbiface; 354a4717116SVijay Mahadevan moab::ParallelComm *pcomm; 355a4717116SVijay Mahadevan moab::Range verts,elems; 3567023aa44SVijay Mahadevan const char *readopts; 357a4717116SVijay Mahadevan PetscErrorCode ierr; 35851d15aeeSVijay Mahadevan 35951d15aeeSVijay Mahadevan PetscFunctionBegin; 360a4717116SVijay Mahadevan PetscValidPointer(dm,5); 36151d15aeeSVijay Mahadevan 362a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 363a4717116SVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 36451d15aeeSVijay Mahadevan 365a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 366a4717116SVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 367a4717116SVijay Mahadevan mbiface = dmmoab->mbiface; 368a4717116SVijay Mahadevan pcomm = dmmoab->pcomm; 369a4717116SVijay Mahadevan nprocs = pcomm->size(); 370c8d5365dSVijay Mahadevan dmmoab->dim = dim; 371a4717116SVijay Mahadevan 3727023aa44SVijay Mahadevan /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */ 3737023aa44SVijay Mahadevan ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab"); 3747023aa44SVijay Mahadevan ierr = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr); 3757023aa44SVijay Mahadevan ierr = PetscOptionsBool("-dmmb_part_by_rank", "Use partition by rank when reading MOAB meshes from file", "dmmbutil.cxx", part_by_rank, &part_by_rank, NULL);CHKERRQ(ierr); 3767023aa44SVijay Mahadevan ierr = PetscOptionsEnd(); 3777023aa44SVijay Mahadevan 378a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */ 3797023aa44SVijay Mahadevan ierr = DMMoab_GetReadOptions_Private(part_by_rank, nprocs, dim, MOAB_PARROPTS_READ_PART, dbglevel, usrreadopts, &readopts);CHKERRQ(ierr); 3807023aa44SVijay Mahadevan 381e5136372SVijay Mahadevan PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts); 382a4717116SVijay Mahadevan 383a4717116SVijay Mahadevan /* Load the mesh from a file. */ 3847023aa44SVijay Mahadevan merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr); 385a4717116SVijay Mahadevan 3866e40195eSVijay Mahadevan /* Reassign global IDs on all entities. */ 387e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 388e5136372SVijay Mahadevan 389e5136372SVijay Mahadevan /* load the local vertices */ 390e5136372SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 391e5136372SVijay Mahadevan /* load the local elements */ 392e5136372SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 393e5136372SVijay Mahadevan 394e5136372SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 395e5136372SVijay Mahadevan on all of the ghost vertexes */ 396e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr); 397e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr); 398e5136372SVijay Mahadevan 399e5136372SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr); 4006e40195eSVijay Mahadevan 401a4717116SVijay Mahadevan merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 402a4717116SVijay Mahadevan 403e5136372SVijay Mahadevan PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size()); 404e5136372SVijay Mahadevan 405e5136372SVijay Mahadevan #if 1 406e5136372SVijay Mahadevan std::stringstream sstr; 407e5136372SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 408e5136372SVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 409e5136372SVijay Mahadevan #endif 41051d15aeeSVijay Mahadevan PetscFunctionReturn(0); 41151d15aeeSVijay Mahadevan } 41251d15aeeSVijay Mahadevan 413