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 49*e5136372SVijay Mahadevan static void set_structured_coordinates(PetscInt i, PetscInt j, PetscInt k, PetscReal hxyz, PetscInt vcount, std::vector<double*>& vcoords) 50*e5136372SVijay Mahadevan { 51*e5136372SVijay Mahadevan vcoords[0][vcount] = i*hxyz; 52*e5136372SVijay Mahadevan vcoords[1][vcount] = j*hxyz; 53*e5136372SVijay Mahadevan vcoords[2][vcount] = k*hxyz; 54*e5136372SVijay Mahadevan } 55*e5136372SVijay Mahadevan 56*e5136372SVijay 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) 57*e5136372SVijay Mahadevan { 58*e5136372SVijay Mahadevan std::vector<int> subent_conn(pow(2,dim)); /* only linear edge, quad, hex supported now */ 59*e5136372SVijay Mahadevan 60*e5136372SVijay Mahadevan moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data()); 61*e5136372SVijay Mahadevan 62*e5136372SVijay Mahadevan switch(dim) { 63*e5136372SVijay Mahadevan case 1: 64*e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i; 65*e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1); 66*e5136372SVijay Mahadevan // PetscPrintf(PETSC_COMM_SELF, "[%D] ELEMENT[%D,%D,%D]: CONNECTIVITY = %D, %D\n", pcomm->rank(), i,j,k,connectivity[offset+subent_conn[0]], connectivity[offset+subent_conn[1]]); 67*e5136372SVijay Mahadevan break; 68*e5136372SVijay Mahadevan case 2: 69*e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1); 70*e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1); 71*e5136372SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1); 72*e5136372SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1); 73*e5136372SVijay Mahadevan // PetscPrintf(PETSC_COMM_SELF, "[%D] ELEMENT[%D,%D,%D]: CONNECTIVITY = %D, %D, %D, %D\n", pcomm->rank(), i,j,k,connectivity[offset+subent_conn[0]], connectivity[offset+subent_conn[1]], connectivity[offset+subent_conn[2]], connectivity[offset+subent_conn[3]]); 74*e5136372SVijay Mahadevan break; 75*e5136372SVijay Mahadevan case 3: 76*e5136372SVijay Mahadevan default: 77*e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 78*e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k); 79*e5136372SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 80*e5136372SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k); 81*e5136372SVijay Mahadevan connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1)); 82*e5136372SVijay Mahadevan connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 83*e5136372SVijay Mahadevan connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1)); 84*e5136372SVijay Mahadevan connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 85*e5136372SVijay Mahadevan break; 86*e5136372SVijay Mahadevan } 87*e5136372SVijay Mahadevan 88*e5136372SVijay Mahadevan } 8951d15aeeSVijay Mahadevan 9051d15aeeSVijay Mahadevan #undef __FUNCT__ 9151d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabCreateBoxMesh" 9251d15aeeSVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscInt nele, PetscInt nghost, DM *dm) 9351d15aeeSVijay Mahadevan { 9451d15aeeSVijay Mahadevan PetscErrorCode ierr; 9551d15aeeSVijay Mahadevan moab::ErrorCode merr; 96*e5136372SVijay Mahadevan PetscInt i,j,k,n,nprocs,rank; 9751d15aeeSVijay Mahadevan DM_Moab *dmmoab; 9851d15aeeSVijay Mahadevan moab::Interface *mbiface; 9951d15aeeSVijay Mahadevan moab::ParallelComm *pcomm; 100a4717116SVijay Mahadevan moab::ReadUtilIface* readMeshIface; 101a4717116SVijay Mahadevan 10251d15aeeSVijay Mahadevan moab::Tag id_tag=PETSC_NULL; 103a4717116SVijay Mahadevan moab::Range ownedvtx,ownedelms; 104a4717116SVijay Mahadevan moab::EntityHandle vfirst,efirst; 105a4717116SVijay Mahadevan std::vector<double*> vcoords; 106a4717116SVijay Mahadevan moab::EntityHandle *connectivity = 0; 10751d15aeeSVijay Mahadevan moab::EntityType etype; 108*e5136372SVijay Mahadevan PetscInt ise[6],imax,imin; 10951d15aeeSVijay Mahadevan PetscReal xse[6]; 11051d15aeeSVijay Mahadevan 111a4717116SVijay Mahadevan const PetscInt npts=nele+1; /* Number of points in every dimension */ 112*e5136372SVijay Mahadevan PetscInt vpere,locnele,locnpts,ghnele,ghnpts; /* Number of verts/element, vertices, elements owned by this process */ 11351d15aeeSVijay Mahadevan 11451d15aeeSVijay Mahadevan PetscFunctionBegin; 115*e5136372SVijay Mahadevan if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n"); 116*e5136372SVijay Mahadevan 117*e5136372SVijay Mahadevan /* total number of vertices in all dimensions */ 118*e5136372SVijay Mahadevan n=pow(npts,dim); 119*e5136372SVijay Mahadevan 120*e5136372SVijay Mahadevan /* do some error checking */ 121*e5136372SVijay Mahadevan if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n"); 122*e5136372SVijay 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"); 123*e5136372SVijay Mahadevan if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n"); 124*e5136372SVijay Mahadevan 125a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 12651d15aeeSVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 12751d15aeeSVijay Mahadevan 128a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 12951d15aeeSVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 13051d15aeeSVijay Mahadevan mbiface = dmmoab->mbiface; 13151d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm; 13251d15aeeSVijay Mahadevan id_tag = dmmoab->ltog_tag; 13351d15aeeSVijay Mahadevan nprocs = pcomm->size(); 134*e5136372SVijay Mahadevan rank = pcomm->rank(); 13551d15aeeSVijay Mahadevan 136a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */ 137a4717116SVijay Mahadevan merr = mbiface->query_interface(readMeshIface);MBERRNM(merr); 13851d15aeeSVijay Mahadevan 13951d15aeeSVijay Mahadevan ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr); 140a4717116SVijay Mahadevan 141a4717116SVijay Mahadevan /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */ 14251d15aeeSVijay Mahadevan ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr); 14351d15aeeSVijay Mahadevan 144a4717116SVijay Mahadevan /* set some variables based on dimension */ 14551d15aeeSVijay Mahadevan switch(dim) { 14651d15aeeSVijay Mahadevan case 1: 14751d15aeeSVijay Mahadevan vpere = 2; 148a4717116SVijay Mahadevan locnele = (ise[1]-ise[0]); 149a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1); 150*e5136372SVijay Mahadevan ghnele = (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0); 151*e5136372SVijay Mahadevan ghnpts = (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0); 15251d15aeeSVijay Mahadevan etype = moab::MBEDGE; 15351d15aeeSVijay Mahadevan break; 15451d15aeeSVijay Mahadevan case 2: 15551d15aeeSVijay Mahadevan vpere = 4; 156a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2]); 157a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1); 158*e5136372SVijay Mahadevan ghnele = (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0); 159*e5136372SVijay Mahadevan ghnpts = (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0); 16051d15aeeSVijay Mahadevan etype = moab::MBQUAD; 16151d15aeeSVijay Mahadevan break; 16251d15aeeSVijay Mahadevan case 3: 16351d15aeeSVijay Mahadevan vpere = 8; 164a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]); 165a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1); 166*e5136372SVijay Mahadevan ghnele = (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0); 167*e5136372SVijay Mahadevan ghnpts = (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0); 16851d15aeeSVijay Mahadevan etype = moab::MBHEX; 16951d15aeeSVijay Mahadevan break; 17051d15aeeSVijay Mahadevan } 171*e5136372SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "[%D] Element Partition Indices -[%D - %D], [%D - %D], [%D - %D]\n", rank, ise[0], ise[1], ise[2], ise[3], ise[4], ise[5]); 17251d15aeeSVijay Mahadevan 17351d15aeeSVijay Mahadevan /* we have a domain of size [1,1,1] - now compute local co-ordinate box */ 17451d15aeeSVijay Mahadevan ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr); 17551d15aeeSVijay Mahadevan for(int i=0; i<6; ++i) { 17651d15aeeSVijay Mahadevan xse[i]=(PetscReal)ise[i]/nele; 17751d15aeeSVijay Mahadevan } 17851d15aeeSVijay Mahadevan 179a4717116SVijay Mahadevan /* Create vertexes and set the coodinate of each vertex */ 180a4717116SVijay Mahadevan const int sequence_size = (locnele + vpere) + 1; 181*e5136372SVijay Mahadevan merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords/*,sequence_size*/);MBERRNM(merr); 18251d15aeeSVijay Mahadevan 183a4717116SVijay Mahadevan /* Compute the co-ordinates of vertices and global IDs */ 184a4717116SVijay Mahadevan std::vector<int> vgid(locnpts); 18551d15aeeSVijay Mahadevan int vcount=0; 18651d15aeeSVijay Mahadevan double hxyz=1.0/nele; 187*e5136372SVijay Mahadevan 188*e5136372SVijay Mahadevan /* create ghosted vertices requested by user - below the current plane */ 189*e5136372SVijay Mahadevan if (ise[2*dim-2] > 0) { 190*e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) { 191*e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) { 192*e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) { 193*e5136372SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "[%D] creating ghost vertex - %D, %D, %D\n", rank, i, j, k); 194*e5136372SVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 195*e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 196*e5136372SVijay Mahadevan } 197*e5136372SVijay Mahadevan } 198*e5136372SVijay Mahadevan } 199*e5136372SVijay Mahadevan } 200*e5136372SVijay Mahadevan 201*e5136372SVijay Mahadevan /* create all the owned vertices */ 202*e5136372SVijay Mahadevan for (k = ise[4]; k <= ise[5]; k++) { 203*e5136372SVijay Mahadevan for (j = ise[2]; j <= ise[3]; j++) { 204*e5136372SVijay Mahadevan for (i = ise[0]; i <= ise[1]; i++, vcount++) { 205*e5136372SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "[%D] creating owned vertex - %D, %D, %D\n", rank, i, j, k); 206*e5136372SVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 207*e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 208*e5136372SVijay Mahadevan } 209*e5136372SVijay Mahadevan } 210*e5136372SVijay Mahadevan } 211*e5136372SVijay Mahadevan 212*e5136372SVijay Mahadevan /* create ghosted vertices requested by user - above the current plane */ 213*e5136372SVijay Mahadevan if (ise[2*dim-1] < nele) { 214*e5136372SVijay Mahadevan for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) { 215*e5136372SVijay Mahadevan for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) { 216*e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) { 217*e5136372SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "[%D] creating ghost vertex - %D, %D, %D\n", rank, i, j, k); 218*e5136372SVijay Mahadevan set_structured_coordinates(i,j,k,hxyz,vcount,vcoords); 219*e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 220*e5136372SVijay Mahadevan } 22151d15aeeSVijay Mahadevan } 22251d15aeeSVijay Mahadevan } 22351d15aeeSVijay Mahadevan } 22451d15aeeSVijay Mahadevan 22551d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 226a4717116SVijay Mahadevan merr = mbiface->add_entities (dmmoab->fileset, ownedvtx);MBERRNM(merr); 227*e5136372SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,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 */ 237*e5136372SVijay 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 */ 240*e5136372SVijay Mahadevan vfirst-=vgid[0]-1; 24151d15aeeSVijay Mahadevan 242a4717116SVijay Mahadevan /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k): 243a4717116SVijay Mahadevan and then set the connectivity for each element appropriately */ 24451d15aeeSVijay Mahadevan int ecount=0; 24551d15aeeSVijay Mahadevan 246*e5136372SVijay Mahadevan /* create ghosted elements requested by user - below the current plane */ 247*e5136372SVijay Mahadevan if (ise[2*dim-2] >= nghost) { 248*e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) { 249*e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) { 250*e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++, ecount++) { 251*e5136372SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "[%D] creating ghost element - %D, %D, %D\n", rank, i, j, k); 252*e5136372SVijay Mahadevan const int offset = ecount*vpere; 253*e5136372SVijay Mahadevan set_element_connectivity(dim, etype, offset, nele, i, j, k, vfirst, connectivity); 25451d15aeeSVijay Mahadevan } 25551d15aeeSVijay Mahadevan } 25651d15aeeSVijay Mahadevan } 25751d15aeeSVijay Mahadevan } 258*e5136372SVijay Mahadevan 259*e5136372SVijay Mahadevan /* create owned elements requested by user */ 260*e5136372SVijay Mahadevan for (k = ise[4]; k < std::max(ise[5],1); k++) { 261*e5136372SVijay Mahadevan for (j = ise[2]; j < std::max(ise[3],1); j++) { 262*e5136372SVijay Mahadevan for (i = ise[0]; i < std::max(ise[1],1); i++,ecount++) { 263*e5136372SVijay Mahadevan const int offset = ecount*vpere; 264*e5136372SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "[%D] creating owned element - %D, %D, %D\n", rank, i, j, k); 265*e5136372SVijay Mahadevan set_element_connectivity(dim, etype, offset, nele, i, j, k, vfirst, connectivity); 266*e5136372SVijay Mahadevan } 267*e5136372SVijay Mahadevan } 268*e5136372SVijay Mahadevan } 269*e5136372SVijay Mahadevan 270*e5136372SVijay Mahadevan /* create ghosted elements requested by user - above the current plane */ 271*e5136372SVijay Mahadevan if (ise[2*dim-1] <= nele-nghost) { 272*e5136372SVijay Mahadevan for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) { 273*e5136372SVijay Mahadevan for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) { 274*e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++, ecount++) { 275*e5136372SVijay Mahadevan PetscPrintf(PETSC_COMM_SELF, "[%D] creating ghost element - %D, %D, %D\n", rank, i, j, k); 276*e5136372SVijay Mahadevan const int offset = ecount*vpere; 277*e5136372SVijay Mahadevan set_element_connectivity(dim, etype, offset, nele, i, j, k, vfirst, connectivity); 278*e5136372SVijay Mahadevan } 279*e5136372SVijay Mahadevan } 280*e5136372SVijay Mahadevan } 281*e5136372SVijay Mahadevan } 282*e5136372SVijay Mahadevan 283*e5136372SVijay Mahadevan merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr); 28451d15aeeSVijay Mahadevan 285a4717116SVijay 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. 286a4717116SVijay Mahadevan first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */ 28751d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr); 288a4717116SVijay Mahadevan merr = mbiface->add_entities (dmmoab->fileset, ownedelms);MBERRNM(merr); 289*e5136372SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr); 29051d15aeeSVijay Mahadevan 291*e5136372SVijay Mahadevan // merr = mbiface->unite_meshset(0, dmmoab->fileset);MBERRV(mbiface,merr); 292a4717116SVijay Mahadevan 293*e5136372SVijay Mahadevan if (locnele+ghnele != (int) ownedelms.size()) 294*e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size()); 295*e5136372SVijay Mahadevan else if(locnpts+ghnpts != (int) ownedvtx.size()) 296*e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size()); 29751d15aeeSVijay Mahadevan else 298a4717116SVijay Mahadevan PetscInfo2(*dm, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size()); 29951d15aeeSVijay Mahadevan 300a4717116SVijay Mahadevan /* check the handles */ 301a4717116SVijay Mahadevan merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr); 30251d15aeeSVijay Mahadevan 303*e5136372SVijay Mahadevan #if 1 304*e5136372SVijay Mahadevan std::stringstream sstr; 305*e5136372SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 306*e5136372SVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 307*e5136372SVijay Mahadevan #endif 308*e5136372SVijay Mahadevan 309a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */ 310a4717116SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr); 311a4717116SVijay Mahadevan merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,&id_tag);MBERRV(mbiface,merr); 31251d15aeeSVijay Mahadevan 313a4717116SVijay Mahadevan /* Reassign global IDs on all entities. */ 314*e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 31551d15aeeSVijay Mahadevan 316a4717116SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 317a4717116SVijay Mahadevan on all of the ghost vertexes */ 318a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr); 319a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr); 32051d15aeeSVijay Mahadevan 321*e5136372SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,nghost,dim,true,true,&dmmoab->fileset);MBERRV(mbiface,merr); 322*e5136372SVijay Mahadevan 323*e5136372SVijay Mahadevan if(rank) { 324*e5136372SVijay Mahadevan ownedvtx.print(0); 325*e5136372SVijay Mahadevan ownedelms.print(0); 326*e5136372SVijay Mahadevan } 32751d15aeeSVijay Mahadevan 32851d15aeeSVijay Mahadevan merr = mbiface->set_dimension(dim);MBERRNM(merr); 32951d15aeeSVijay Mahadevan 33051d15aeeSVijay Mahadevan PetscFunctionReturn(0); 33151d15aeeSVijay Mahadevan } 33251d15aeeSVijay Mahadevan 33351d15aeeSVijay Mahadevan 334a4717116SVijay Mahadevan #undef __FUNCT__ 335a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private" 3367023aa44SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* extra_opts, const char** read_opts) 33751d15aeeSVijay Mahadevan { 338a4717116SVijay Mahadevan std::ostringstream str; 33951d15aeeSVijay Mahadevan 340a4717116SVijay Mahadevan PetscFunctionBegin; 341a4717116SVijay Mahadevan // do parallel read unless only one processor 342a4717116SVijay Mahadevan if (numproc > 1) { 3437023aa44SVijay Mahadevan str << "PARALLEL=" << mode << ";PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"; 3447023aa44SVijay Mahadevan str << "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;"; 345a4717116SVijay Mahadevan if (by_rank) 346a4717116SVijay Mahadevan str << "PARTITION_BY_RANK;"; 347a4717116SVijay Mahadevan } 34851d15aeeSVijay Mahadevan 349a4717116SVijay Mahadevan if (extra_opts) 350a4717116SVijay Mahadevan str << extra_opts << ";"; 35151d15aeeSVijay Mahadevan 352a4717116SVijay Mahadevan if (dbglevel) 3537023aa44SVijay Mahadevan str << "DEBUG_IO=" << dbglevel << ";DEBUG_PIO=" << dbglevel << ";CPUTIME;"; 35451d15aeeSVijay Mahadevan 3557023aa44SVijay Mahadevan *read_opts = str.str().c_str(); 35651d15aeeSVijay Mahadevan PetscFunctionReturn(0); 35751d15aeeSVijay Mahadevan } 35851d15aeeSVijay Mahadevan 35951d15aeeSVijay Mahadevan 36051d15aeeSVijay Mahadevan #undef __FUNCT__ 361a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile" 362a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm) 36351d15aeeSVijay Mahadevan { 364a4717116SVijay Mahadevan moab::ErrorCode merr; 3657023aa44SVijay Mahadevan PetscInt nprocs,dbglevel=0; 3667023aa44SVijay Mahadevan PetscBool part_by_rank=PETSC_FALSE; 367a4717116SVijay Mahadevan DM_Moab *dmmoab; 368a4717116SVijay Mahadevan moab::Interface *mbiface; 369a4717116SVijay Mahadevan moab::ParallelComm *pcomm; 370a4717116SVijay Mahadevan moab::Range verts,elems; 3717023aa44SVijay Mahadevan const char *readopts; 372a4717116SVijay Mahadevan PetscErrorCode ierr; 37351d15aeeSVijay Mahadevan 37451d15aeeSVijay Mahadevan PetscFunctionBegin; 375a4717116SVijay Mahadevan PetscValidPointer(dm,5); 37651d15aeeSVijay Mahadevan 377a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 378a4717116SVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 37951d15aeeSVijay Mahadevan 380a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 381a4717116SVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 382a4717116SVijay Mahadevan mbiface = dmmoab->mbiface; 383a4717116SVijay Mahadevan pcomm = dmmoab->pcomm; 384a4717116SVijay Mahadevan nprocs = pcomm->size(); 385a4717116SVijay Mahadevan 3867023aa44SVijay Mahadevan /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */ 3877023aa44SVijay Mahadevan ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab"); 3887023aa44SVijay Mahadevan ierr = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr); 3897023aa44SVijay 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); 3907023aa44SVijay Mahadevan ierr = PetscOptionsEnd(); 3917023aa44SVijay Mahadevan 392a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */ 3937023aa44SVijay Mahadevan ierr = DMMoab_GetReadOptions_Private(part_by_rank, nprocs, dim, MOAB_PARROPTS_READ_PART, dbglevel, usrreadopts, &readopts);CHKERRQ(ierr); 3947023aa44SVijay Mahadevan 395*e5136372SVijay Mahadevan PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts); 396a4717116SVijay Mahadevan 397a4717116SVijay Mahadevan /* Load the mesh from a file. */ 3987023aa44SVijay Mahadevan merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr); 399a4717116SVijay Mahadevan 4006e40195eSVijay Mahadevan /* Reassign global IDs on all entities. */ 401*e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 402*e5136372SVijay Mahadevan 403*e5136372SVijay Mahadevan /* load the local vertices */ 404*e5136372SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 405*e5136372SVijay Mahadevan /* load the local elements */ 406*e5136372SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 407*e5136372SVijay Mahadevan 408*e5136372SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 409*e5136372SVijay Mahadevan on all of the ghost vertexes */ 410*e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr); 411*e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr); 412*e5136372SVijay Mahadevan 413*e5136372SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr); 4146e40195eSVijay Mahadevan 415a4717116SVijay Mahadevan merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 416a4717116SVijay Mahadevan 417*e5136372SVijay Mahadevan PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size()); 418*e5136372SVijay Mahadevan 419a4717116SVijay Mahadevan merr = mbiface->set_dimension(dim);MBERRNM(merr); 420a4717116SVijay Mahadevan 421*e5136372SVijay Mahadevan #if 1 422*e5136372SVijay Mahadevan std::stringstream sstr; 423*e5136372SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 424*e5136372SVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 425*e5136372SVijay Mahadevan #endif 426*e5136372SVijay Mahadevan 42751d15aeeSVijay Mahadevan PetscFunctionReturn(0); 42851d15aeeSVijay Mahadevan } 42951d15aeeSVijay Mahadevan 430