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" 89*3a4aeca1SVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt nele, PetscInt user_nghost, DM *dm) 9051d15aeeSVijay Mahadevan { 9151d15aeeSVijay Mahadevan PetscErrorCode ierr; 9251d15aeeSVijay Mahadevan moab::ErrorCode merr; 93*3a4aeca1SVijay Mahadevan PetscInt i,j,k,n,nprocs; 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; 101*3a4aeca1SVijay Mahadevan moab::EntityHandle vfirst,efirst,regionset,faceset,edgeset,vtxset; 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]; 107*3a4aeca1SVijay Mahadevan /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */ 108*3a4aeca1SVijay Mahadevan const PetscInt nghost=0; 109*3a4aeca1SVijay Mahadevan 110*3a4aeca1SVijay Mahadevan moab::Tag geom_tag; 111*3a4aeca1SVijay Mahadevan 112*3a4aeca1SVijay Mahadevan moab::Range adj,dim3,dim2; 113*3a4aeca1SVijay Mahadevan bool build_adjacencies=false; 11451d15aeeSVijay Mahadevan 115a4717116SVijay Mahadevan const PetscInt npts=nele+1; /* Number of points in every dimension */ 116e5136372SVijay Mahadevan PetscInt vpere,locnele,locnpts,ghnele,ghnpts; /* Number of verts/element, vertices, elements owned by this process */ 11751d15aeeSVijay Mahadevan 11851d15aeeSVijay Mahadevan PetscFunctionBegin; 119e5136372SVijay Mahadevan if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n"); 120e5136372SVijay Mahadevan 121c8d5365dSVijay Mahadevan ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr); 122e5136372SVijay Mahadevan /* total number of vertices in all dimensions */ 123e5136372SVijay Mahadevan n=pow(npts,dim); 124e5136372SVijay Mahadevan 125e5136372SVijay Mahadevan /* do some error checking */ 126e5136372SVijay Mahadevan if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n"); 127e5136372SVijay 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"); 128e5136372SVijay Mahadevan if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n"); 129e5136372SVijay Mahadevan 130a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 13151d15aeeSVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 13251d15aeeSVijay Mahadevan 133a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 13451d15aeeSVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 13551d15aeeSVijay Mahadevan mbiface = dmmoab->mbiface; 13651d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm; 13751d15aeeSVijay Mahadevan id_tag = dmmoab->ltog_tag; 13851d15aeeSVijay Mahadevan nprocs = pcomm->size(); 139c8d5365dSVijay Mahadevan dmmoab->dim = dim; 14051d15aeeSVijay Mahadevan 141a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */ 142a4717116SVijay Mahadevan merr = mbiface->query_interface(readMeshIface);MBERRNM(merr); 14351d15aeeSVijay Mahadevan 14451d15aeeSVijay Mahadevan ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr); 145a4717116SVijay Mahadevan 146a4717116SVijay Mahadevan /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */ 14751d15aeeSVijay Mahadevan ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr); 14851d15aeeSVijay Mahadevan 149a4717116SVijay Mahadevan /* set some variables based on dimension */ 15051d15aeeSVijay Mahadevan switch(dim) { 15151d15aeeSVijay Mahadevan case 1: 15251d15aeeSVijay Mahadevan vpere = 2; 153a4717116SVijay Mahadevan locnele = (ise[1]-ise[0]); 154a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1); 155c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 ); 156c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0); 15751d15aeeSVijay Mahadevan etype = moab::MBEDGE; 15851d15aeeSVijay Mahadevan break; 15951d15aeeSVijay Mahadevan case 2: 16051d15aeeSVijay Mahadevan vpere = 4; 161a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2]); 162a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1); 163c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0); 164c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0); 16551d15aeeSVijay Mahadevan etype = moab::MBQUAD; 16651d15aeeSVijay Mahadevan break; 16751d15aeeSVijay Mahadevan case 3: 16851d15aeeSVijay Mahadevan vpere = 8; 169a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]); 170a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1); 171c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0); 172c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0); 17351d15aeeSVijay Mahadevan etype = moab::MBHEX; 17451d15aeeSVijay Mahadevan break; 17551d15aeeSVijay Mahadevan } 17651d15aeeSVijay Mahadevan 17751d15aeeSVijay Mahadevan /* we have a domain of size [1,1,1] - now compute local co-ordinate box */ 17851d15aeeSVijay Mahadevan ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr); 179c8d5365dSVijay Mahadevan for(i=0; i<6; ++i) { 18051d15aeeSVijay Mahadevan xse[i]=(PetscReal)ise[i]/nele; 18151d15aeeSVijay Mahadevan } 18251d15aeeSVijay Mahadevan 183a4717116SVijay Mahadevan /* Create vertexes and set the coodinate of each vertex */ 184c8d5365dSVijay Mahadevan merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr); 18551d15aeeSVijay Mahadevan 186a4717116SVijay Mahadevan /* Compute the co-ordinates of vertices and global IDs */ 187c8d5365dSVijay Mahadevan std::vector<int> vgid(locnpts+ghnpts); 18851d15aeeSVijay Mahadevan int vcount=0; 18951d15aeeSVijay Mahadevan double hxyz=1.0/nele; 190e5136372SVijay Mahadevan 191c8d5365dSVijay Mahadevan /* create all the owned vertices */ 192c8d5365dSVijay Mahadevan for (k = ise[4]; k <= ise[5]; k++) { 193c8d5365dSVijay Mahadevan for (j = ise[2]; j <= ise[3]; j++) { 194c8d5365dSVijay Mahadevan for (i = ise[0]; i <= ise[1]; i++, vcount++) { 195c8d5365dSVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 196c8d5365dSVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 197c8d5365dSVijay Mahadevan } 198c8d5365dSVijay Mahadevan } 199c8d5365dSVijay Mahadevan } 200c8d5365dSVijay Mahadevan 201e5136372SVijay Mahadevan /* create ghosted vertices requested by user - below the current plane */ 202e5136372SVijay Mahadevan if (ise[2*dim-2] > 0) { 203e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) { 204e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) { 205e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) { 206e5136372SVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 207e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 208e5136372SVijay Mahadevan } 209e5136372SVijay Mahadevan } 210e5136372SVijay Mahadevan } 211e5136372SVijay Mahadevan } 212e5136372SVijay Mahadevan 213e5136372SVijay Mahadevan /* create ghosted vertices requested by user - above the current plane */ 214e5136372SVijay Mahadevan if (ise[2*dim-1] < nele) { 215e5136372SVijay Mahadevan for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) { 216e5136372SVijay Mahadevan for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) { 217e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) { 218e5136372SVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 219e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 220e5136372SVijay Mahadevan } 22151d15aeeSVijay Mahadevan } 22251d15aeeSVijay Mahadevan } 22351d15aeeSVijay Mahadevan } 22451d15aeeSVijay Mahadevan 225c8d5365dSVijay Mahadevan if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount); 226c8d5365dSVijay Mahadevan 22751d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 22851d15aeeSVijay Mahadevan 229a4717116SVijay Mahadevan /* The global ID tag is applied to each owned 230a4717116SVijay Mahadevan vertex. It acts as an global identifier which MOAB uses to 231a4717116SVijay Mahadevan assemble the individual pieces of the mesh: 232a4717116SVijay Mahadevan Set the global ID indices */ 23351d15aeeSVijay Mahadevan merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr); 23451d15aeeSVijay Mahadevan 235a4717116SVijay Mahadevan /* Create elements between mesh points using the ReadUtilInterface 236a4717116SVijay Mahadevan get the reference to element connectivities for all local elements from the ReadUtilInterface */ 237e5136372SVijay Mahadevan merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr); 23851d15aeeSVijay Mahadevan 239a4717116SVijay Mahadevan /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */ 240c8d5365dSVijay Mahadevan // PetscPrintf(PETSC_COMM_SELF, "[%D] first local handle %D\n", rank, vfirst); 241e5136372SVijay Mahadevan vfirst-=vgid[0]-1; 24251d15aeeSVijay Mahadevan 243a4717116SVijay Mahadevan /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k): 244a4717116SVijay Mahadevan and then set the connectivity for each element appropriately */ 24551d15aeeSVijay Mahadevan int ecount=0; 24651d15aeeSVijay Mahadevan 247e5136372SVijay Mahadevan /* create ghosted elements requested by user - below the current plane */ 248e5136372SVijay Mahadevan if (ise[2*dim-2] >= nghost) { 249e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) { 250e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) { 251e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++, ecount++) { 252c8d5365dSVijay Mahadevan set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 25351d15aeeSVijay Mahadevan } 25451d15aeeSVijay Mahadevan } 25551d15aeeSVijay Mahadevan } 25651d15aeeSVijay Mahadevan } 257e5136372SVijay Mahadevan 258e5136372SVijay Mahadevan /* create owned elements requested by user */ 259e5136372SVijay Mahadevan for (k = ise[4]; k < std::max(ise[5],1); k++) { 260e5136372SVijay Mahadevan for (j = ise[2]; j < std::max(ise[3],1); j++) { 261e5136372SVijay Mahadevan for (i = ise[0]; i < std::max(ise[1],1); i++,ecount++) { 262c8d5365dSVijay Mahadevan set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 263e5136372SVijay Mahadevan } 264e5136372SVijay Mahadevan } 265e5136372SVijay Mahadevan } 266e5136372SVijay Mahadevan 267e5136372SVijay Mahadevan /* create ghosted elements requested by user - above the current plane */ 268e5136372SVijay Mahadevan if (ise[2*dim-1] <= nele-nghost) { 269e5136372SVijay Mahadevan for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) { 270e5136372SVijay Mahadevan for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) { 271e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++, ecount++) { 272c8d5365dSVijay Mahadevan set_element_connectivity(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 273e5136372SVijay Mahadevan } 274e5136372SVijay Mahadevan } 275e5136372SVijay Mahadevan } 276e5136372SVijay Mahadevan } 277e5136372SVijay Mahadevan 278e5136372SVijay Mahadevan merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr); 27951d15aeeSVijay Mahadevan 280a4717116SVijay 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. 281a4717116SVijay Mahadevan first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */ 28251d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr); 283a4717116SVijay Mahadevan 284e5136372SVijay Mahadevan if (locnele+ghnele != (int) ownedelms.size()) 285e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size()); 286e5136372SVijay Mahadevan else if(locnpts+ghnpts != (int) ownedvtx.size()) 287e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size()); 28851d15aeeSVijay Mahadevan else 289a4717116SVijay Mahadevan PetscInfo2(*dm, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size()); 29051d15aeeSVijay Mahadevan 291*3a4aeca1SVijay Mahadevan /* lets create some sets */ 292*3a4aeca1SVijay Mahadevan merr = mbiface->tag_get_handle(GEOM_DIMENSION_TAG_NAME, 1, moab::MB_TYPE_INTEGER, geom_tag, moab::MB_TAG_SPARSE|moab::MB_TAG_CREAT);MBERRNM(merr); 293*3a4aeca1SVijay Mahadevan 294*3a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr); 295*3a4aeca1SVijay Mahadevan merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr); 296*3a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr); 297*3a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, ®ionset, 1, &dmmoab->dim);MBERRNM(merr); 298*3a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr); 299*3a4aeca1SVijay Mahadevan 300*3a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr); 301*3a4aeca1SVijay Mahadevan merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr); 302*3a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr); 303*3a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr); 304*3a4aeca1SVijay Mahadevan 305*3a4aeca1SVijay Mahadevan if (build_adjacencies) { 306*3a4aeca1SVijay Mahadevan // generate all lower dimensional adjacencies 307*3a4aeca1SVijay Mahadevan merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr); 308*3a4aeca1SVijay Mahadevan merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr); 309*3a4aeca1SVijay Mahadevan adj.merge(dim2); 310*3a4aeca1SVijay Mahadevan 311*3a4aeca1SVijay Mahadevan /* create face sets */ 312*3a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr); 313*3a4aeca1SVijay Mahadevan merr = mbiface->add_entities(faceset, adj);MBERRNM(merr); 314*3a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr); 315*3a4aeca1SVijay Mahadevan i=dim-1; 316*3a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr); 317*3a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr); 318*3a4aeca1SVijay Mahadevan PetscInfo2(dm, "Found %d %d-Dim quantities.\n", adj.size(), dim-1); 319*3a4aeca1SVijay Mahadevan 320*3a4aeca1SVijay Mahadevan if (dim > 2) { 321*3a4aeca1SVijay Mahadevan dim2.clear(); 322*3a4aeca1SVijay Mahadevan /* create edge sets, if appropriate i.e., if dim=3 */ 323*3a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr); 324*3a4aeca1SVijay Mahadevan merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr); 325*3a4aeca1SVijay Mahadevan merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr); 326*3a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr); 327*3a4aeca1SVijay Mahadevan i=dim-2; 328*3a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr); 329*3a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr); 330*3a4aeca1SVijay Mahadevan PetscInfo2(dm, "Found %d %d-Dim quantities.\n", adj.size(), dim-2); 331*3a4aeca1SVijay Mahadevan } 332*3a4aeca1SVijay Mahadevan } 333*3a4aeca1SVijay Mahadevan 334a4717116SVijay Mahadevan /* check the handles */ 335a4717116SVijay Mahadevan merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr); 33651d15aeeSVijay Mahadevan 337c8d5365dSVijay Mahadevan #if 0 338e5136372SVijay Mahadevan std::stringstream sstr; 339e5136372SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 340e5136372SVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 341e5136372SVijay Mahadevan #endif 342e5136372SVijay Mahadevan 343a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */ 344a4717116SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr); 345a4717116SVijay Mahadevan merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,&id_tag);MBERRV(mbiface,merr); 34651d15aeeSVijay Mahadevan 347a4717116SVijay Mahadevan /* Reassign global IDs on all entities. */ 348e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 34951d15aeeSVijay Mahadevan 350*3a4aeca1SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr); 351c8d5365dSVijay Mahadevan 352a4717116SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 353a4717116SVijay Mahadevan on all of the ghost vertexes */ 354c8d5365dSVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 355c8d5365dSVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr); 356c8d5365dSVijay Mahadevan 357a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr); 358a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr); 35951d15aeeSVijay Mahadevan PetscFunctionReturn(0); 36051d15aeeSVijay Mahadevan } 36151d15aeeSVijay Mahadevan 36251d15aeeSVijay Mahadevan 363a4717116SVijay Mahadevan #undef __FUNCT__ 364a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private" 3657023aa44SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* extra_opts, const char** read_opts) 36651d15aeeSVijay Mahadevan { 367a4717116SVijay Mahadevan std::ostringstream str; 36851d15aeeSVijay Mahadevan 369a4717116SVijay Mahadevan PetscFunctionBegin; 370e23c60ebSVijay Mahadevan /* do parallel read unless using only one processor */ 371a4717116SVijay Mahadevan if (numproc > 1) { 3727023aa44SVijay Mahadevan str << "PARALLEL=" << mode << ";PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"; 3737023aa44SVijay Mahadevan str << "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;"; 374a4717116SVijay Mahadevan if (by_rank) 375a4717116SVijay Mahadevan str << "PARTITION_BY_RANK;"; 376a4717116SVijay Mahadevan } 37751d15aeeSVijay Mahadevan 378c8d5365dSVijay Mahadevan if (dbglevel) { 379c8d5365dSVijay Mahadevan str << "CPUTIME;DEBUG_IO=" << dbglevel << ";"; 380c8d5365dSVijay Mahadevan if (numproc>1) 381c8d5365dSVijay Mahadevan str << "DEBUG_PIO=" << dbglevel << ";"; 382c8d5365dSVijay Mahadevan } 38351d15aeeSVijay Mahadevan 384c8d5365dSVijay Mahadevan if (extra_opts) 385c8d5365dSVijay Mahadevan str << extra_opts; 38651d15aeeSVijay Mahadevan 3877023aa44SVijay Mahadevan *read_opts = str.str().c_str(); 38851d15aeeSVijay Mahadevan PetscFunctionReturn(0); 38951d15aeeSVijay Mahadevan } 39051d15aeeSVijay Mahadevan 39151d15aeeSVijay Mahadevan 39251d15aeeSVijay Mahadevan #undef __FUNCT__ 393a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile" 394a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm) 39551d15aeeSVijay Mahadevan { 396a4717116SVijay Mahadevan moab::ErrorCode merr; 3977023aa44SVijay Mahadevan PetscInt nprocs,dbglevel=0; 3987023aa44SVijay Mahadevan PetscBool part_by_rank=PETSC_FALSE; 399a4717116SVijay Mahadevan DM_Moab *dmmoab; 400a4717116SVijay Mahadevan moab::Interface *mbiface; 401a4717116SVijay Mahadevan moab::ParallelComm *pcomm; 402a4717116SVijay Mahadevan moab::Range verts,elems; 4037023aa44SVijay Mahadevan const char *readopts; 404a4717116SVijay Mahadevan PetscErrorCode ierr; 40551d15aeeSVijay Mahadevan 40651d15aeeSVijay Mahadevan PetscFunctionBegin; 407a4717116SVijay Mahadevan PetscValidPointer(dm,5); 40851d15aeeSVijay Mahadevan 409a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 410a4717116SVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 41151d15aeeSVijay Mahadevan 412a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 413a4717116SVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 414a4717116SVijay Mahadevan mbiface = dmmoab->mbiface; 415a4717116SVijay Mahadevan pcomm = dmmoab->pcomm; 416a4717116SVijay Mahadevan nprocs = pcomm->size(); 417c8d5365dSVijay Mahadevan dmmoab->dim = dim; 418a4717116SVijay Mahadevan 4197023aa44SVijay Mahadevan /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */ 4207023aa44SVijay Mahadevan ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab"); 4217023aa44SVijay Mahadevan ierr = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr); 4227023aa44SVijay 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); 4237023aa44SVijay Mahadevan ierr = PetscOptionsEnd(); 4247023aa44SVijay Mahadevan 425a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */ 4267023aa44SVijay Mahadevan ierr = DMMoab_GetReadOptions_Private(part_by_rank, nprocs, dim, MOAB_PARROPTS_READ_PART, dbglevel, usrreadopts, &readopts);CHKERRQ(ierr); 4277023aa44SVijay Mahadevan 428e5136372SVijay Mahadevan PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts); 429a4717116SVijay Mahadevan 430a4717116SVijay Mahadevan /* Load the mesh from a file. */ 4317023aa44SVijay Mahadevan merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr); 432a4717116SVijay Mahadevan 4336e40195eSVijay Mahadevan /* Reassign global IDs on all entities. */ 434e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 435e5136372SVijay Mahadevan 436e5136372SVijay Mahadevan /* load the local vertices */ 437e5136372SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 438e5136372SVijay Mahadevan /* load the local elements */ 439e5136372SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 440e5136372SVijay Mahadevan 441e5136372SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 442e5136372SVijay Mahadevan on all of the ghost vertexes */ 443e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr); 444e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr); 445e5136372SVijay Mahadevan 446e5136372SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr); 4476e40195eSVijay Mahadevan 448a4717116SVijay Mahadevan merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 449a4717116SVijay Mahadevan 450e5136372SVijay Mahadevan PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size()); 451e5136372SVijay Mahadevan 452e5136372SVijay Mahadevan #if 1 453e5136372SVijay Mahadevan std::stringstream sstr; 454e5136372SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 455e5136372SVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 456e5136372SVijay Mahadevan #endif 45751d15aeeSVijay Mahadevan PetscFunctionReturn(0); 45851d15aeeSVijay Mahadevan } 45951d15aeeSVijay Mahadevan 460