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" 13*aa859218SVijay Mahadevan static 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*aa859218SVijay Mahadevan #undef __FUNCT__ 50*aa859218SVijay Mahadevan #define __FUNCT__ "DMMoab_SetStructuredCoords_Private" 51*aa859218SVijay Mahadevan static void DMMoab_SetStructuredCoords_Private(PetscInt i, PetscInt j, PetscInt k, PetscReal hx, PetscReal hy, PetscReal hz, PetscInt vcount, std::vector<double*>& vcoords) 52e5136372SVijay Mahadevan { 5366f88a78SVijay Mahadevan vcoords[0][vcount] = i*hx; 5466f88a78SVijay Mahadevan vcoords[1][vcount] = j*hy; 5566f88a78SVijay Mahadevan vcoords[2][vcount] = k*hz; 56e5136372SVijay Mahadevan } 57e5136372SVijay Mahadevan 58*aa859218SVijay Mahadevan #undef __FUNCT__ 59*aa859218SVijay Mahadevan #define __FUNCT__ "DMMoab_SetElementConnectivity_Private" 60*aa859218SVijay Mahadevan static void DMMoab_SetElementConnectivity_Private(PetscInt dim, moab::EntityType etype, PetscInt offset, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity) 61e5136372SVijay Mahadevan { 62e5136372SVijay Mahadevan std::vector<int> subent_conn(pow(2,dim)); /* only linear edge, quad, hex supported now */ 63e5136372SVijay Mahadevan 64e5136372SVijay Mahadevan moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data()); 65e5136372SVijay Mahadevan 66e5136372SVijay Mahadevan switch(dim) { 67e5136372SVijay Mahadevan case 1: 68e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i; 69e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1); 70e5136372SVijay Mahadevan break; 71e5136372SVijay Mahadevan case 2: 72e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1); 73e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1); 74e5136372SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1); 75e5136372SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1); 76e5136372SVijay Mahadevan break; 77e5136372SVijay Mahadevan case 3: 78e5136372SVijay Mahadevan default: 79e5136372SVijay Mahadevan connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 80e5136372SVijay Mahadevan connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k); 81e5136372SVijay Mahadevan connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 82e5136372SVijay Mahadevan connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k); 83e5136372SVijay Mahadevan connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1)); 84e5136372SVijay Mahadevan connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 85e5136372SVijay Mahadevan connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1)); 86e5136372SVijay Mahadevan connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 87e5136372SVijay Mahadevan break; 88e5136372SVijay Mahadevan } 89e5136372SVijay Mahadevan } 9051d15aeeSVijay Mahadevan 9151d15aeeSVijay Mahadevan #undef __FUNCT__ 9251d15aeeSVijay Mahadevan #define __FUNCT__ "DMMoabCreateBoxMesh" 93*aa859218SVijay Mahadevan /*@ 94*aa859218SVijay Mahadevan DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with user specified bounds. 95*aa859218SVijay Mahadevan 96*aa859218SVijay Mahadevan Collective on MPI_Comm 97*aa859218SVijay Mahadevan 98*aa859218SVijay Mahadevan Input Parameters: 99*aa859218SVijay Mahadevan + comm - The communicator for the DM object 100*aa859218SVijay Mahadevan . dim - The spatial dimension 101*aa859218SVijay Mahadevan - bounds - The bounds of the box specified with [x-left, x-right, y-bottom, y-top, z-bottom, z-top] depending on the spatial dimension 102*aa859218SVijay Mahadevan - nele - The number of discrete elements in each direction 103*aa859218SVijay Mahadevan - user_nghost - The number of ghosted layers needed in the partitioned mesh 104*aa859218SVijay Mahadevan 105*aa859218SVijay Mahadevan Output Parameter: 106*aa859218SVijay Mahadevan . dm - The DM object 107*aa859218SVijay Mahadevan 108*aa859218SVijay Mahadevan Level: beginner 109*aa859218SVijay Mahadevan 110*aa859218SVijay Mahadevan .keywords: DM, create 111*aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile() 112*aa859218SVijay Mahadevan @*/ 11366f88a78SVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, const PetscReal* bounds, PetscInt nele, PetscInt user_nghost, DM *dm) 11451d15aeeSVijay Mahadevan { 11551d15aeeSVijay Mahadevan PetscErrorCode ierr; 11651d15aeeSVijay Mahadevan moab::ErrorCode merr; 1173a4aeca1SVijay Mahadevan PetscInt i,j,k,n,nprocs; 11851d15aeeSVijay Mahadevan DM_Moab *dmmoab; 11951d15aeeSVijay Mahadevan moab::Interface *mbiface; 12051d15aeeSVijay Mahadevan moab::ParallelComm *pcomm; 121a4717116SVijay Mahadevan moab::ReadUtilIface* readMeshIface; 122a4717116SVijay Mahadevan 12351d15aeeSVijay Mahadevan moab::Tag id_tag=PETSC_NULL; 124a4717116SVijay Mahadevan moab::Range ownedvtx,ownedelms; 1253a4aeca1SVijay Mahadevan moab::EntityHandle vfirst,efirst,regionset,faceset,edgeset,vtxset; 126a4717116SVijay Mahadevan std::vector<double*> vcoords; 127a4717116SVijay Mahadevan moab::EntityHandle *connectivity = 0; 12851d15aeeSVijay Mahadevan moab::EntityType etype; 129c8d5365dSVijay Mahadevan PetscInt ise[6]; 13066f88a78SVijay Mahadevan PetscReal xse[6],defbounds[6]; 1313a4aeca1SVijay Mahadevan /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */ 1323a4aeca1SVijay Mahadevan const PetscInt nghost=0; 1333a4aeca1SVijay Mahadevan 1343a4aeca1SVijay Mahadevan moab::Tag geom_tag; 1353a4aeca1SVijay Mahadevan 1363a4aeca1SVijay Mahadevan moab::Range adj,dim3,dim2; 1373a4aeca1SVijay Mahadevan bool build_adjacencies=false; 13851d15aeeSVijay Mahadevan 139a4717116SVijay Mahadevan const PetscInt npts=nele+1; /* Number of points in every dimension */ 140e5136372SVijay Mahadevan PetscInt vpere,locnele,locnpts,ghnele,ghnpts; /* Number of verts/element, vertices, elements owned by this process */ 14151d15aeeSVijay Mahadevan 14251d15aeeSVijay Mahadevan PetscFunctionBegin; 143e5136372SVijay Mahadevan if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n"); 144e5136372SVijay Mahadevan 145c8d5365dSVijay Mahadevan ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr); 146e5136372SVijay Mahadevan /* total number of vertices in all dimensions */ 147e5136372SVijay Mahadevan n=pow(npts,dim); 148e5136372SVijay Mahadevan 149e5136372SVijay Mahadevan /* do some error checking */ 150e5136372SVijay Mahadevan if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n"); 151e5136372SVijay 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"); 152e5136372SVijay Mahadevan if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n"); 153e5136372SVijay Mahadevan 154a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 15551d15aeeSVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 15651d15aeeSVijay Mahadevan 157a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 15851d15aeeSVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 15951d15aeeSVijay Mahadevan mbiface = dmmoab->mbiface; 16051d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm; 16151d15aeeSVijay Mahadevan id_tag = dmmoab->ltog_tag; 16251d15aeeSVijay Mahadevan nprocs = pcomm->size(); 163c8d5365dSVijay Mahadevan dmmoab->dim = dim; 16451d15aeeSVijay Mahadevan 165b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */ 166b5410836SVijay Mahadevan merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 167b5410836SVijay Mahadevan 168a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */ 169a4717116SVijay Mahadevan merr = mbiface->query_interface(readMeshIface);MBERRNM(merr); 17051d15aeeSVijay Mahadevan 17151d15aeeSVijay Mahadevan ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr); 172a4717116SVijay Mahadevan 173a4717116SVijay Mahadevan /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */ 17451d15aeeSVijay Mahadevan ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr); 17551d15aeeSVijay Mahadevan 176a4717116SVijay Mahadevan /* set some variables based on dimension */ 17751d15aeeSVijay Mahadevan switch(dim) { 17851d15aeeSVijay Mahadevan case 1: 17951d15aeeSVijay Mahadevan vpere = 2; 180a4717116SVijay Mahadevan locnele = (ise[1]-ise[0]); 181a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1); 182c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 ); 183c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0); 18451d15aeeSVijay Mahadevan etype = moab::MBEDGE; 18551d15aeeSVijay Mahadevan break; 18651d15aeeSVijay Mahadevan case 2: 18751d15aeeSVijay Mahadevan vpere = 4; 188a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2]); 189a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1); 190c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0); 191c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0); 19251d15aeeSVijay Mahadevan etype = moab::MBQUAD; 19351d15aeeSVijay Mahadevan break; 19451d15aeeSVijay Mahadevan case 3: 19551d15aeeSVijay Mahadevan vpere = 8; 196a4717116SVijay Mahadevan locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]); 197a4717116SVijay Mahadevan locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1); 198c8d5365dSVijay Mahadevan ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0); 199c8d5365dSVijay Mahadevan ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0); 20051d15aeeSVijay Mahadevan etype = moab::MBHEX; 20151d15aeeSVijay Mahadevan break; 20251d15aeeSVijay Mahadevan } 20351d15aeeSVijay Mahadevan 20451d15aeeSVijay Mahadevan /* we have a domain of size [1,1,1] - now compute local co-ordinate box */ 20551d15aeeSVijay Mahadevan ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr); 206c8d5365dSVijay Mahadevan for(i=0; i<6; ++i) { 20751d15aeeSVijay Mahadevan xse[i]=(PetscReal)ise[i]/nele; 20851d15aeeSVijay Mahadevan } 20951d15aeeSVijay Mahadevan 210a4717116SVijay Mahadevan /* Create vertexes and set the coodinate of each vertex */ 211c8d5365dSVijay Mahadevan merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr); 21251d15aeeSVijay Mahadevan 213a4717116SVijay Mahadevan /* Compute the co-ordinates of vertices and global IDs */ 214c8d5365dSVijay Mahadevan std::vector<int> vgid(locnpts+ghnpts); 21551d15aeeSVijay Mahadevan int vcount=0; 21666f88a78SVijay Mahadevan 21766f88a78SVijay Mahadevan if (!bounds) { /* default box mesh is defined on a unit-cube */ 21866f88a78SVijay Mahadevan defbounds[0]=0.0; defbounds[1]=1.0; 21966f88a78SVijay Mahadevan defbounds[2]=0.0; defbounds[3]=1.0; 22066f88a78SVijay Mahadevan defbounds[4]=0.0; defbounds[5]=1.0; 22166f88a78SVijay Mahadevan bounds=defbounds; 22266f88a78SVijay Mahadevan } 22366f88a78SVijay Mahadevan else { 22466f88a78SVijay Mahadevan /* validate the bounds data */ 22566f88a78SVijay Mahadevan if(bounds[0] >= bounds[1]) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"X-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[0],bounds[1]); 22666f88a78SVijay Mahadevan if(dim > 1 && (bounds[2] >= bounds[3])) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Y-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[2],bounds[3]); 22766f88a78SVijay Mahadevan if(dim > 2 && (bounds[4] >= bounds[5])) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Z-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[4],bounds[5]); 22866f88a78SVijay Mahadevan } 22966f88a78SVijay Mahadevan 23066f88a78SVijay Mahadevan const double hx=(bounds[1]-bounds[0])/nele; 23166f88a78SVijay Mahadevan const double hy=(dim > 1 ? (bounds[3]-bounds[2])/nele : 0.0); 23266f88a78SVijay Mahadevan const double hz=(dim > 2 ? (bounds[5]-bounds[4])/nele : 0.0); 233e5136372SVijay Mahadevan 234c8d5365dSVijay Mahadevan /* create all the owned vertices */ 235c8d5365dSVijay Mahadevan for (k = ise[4]; k <= ise[5]; k++) { 236c8d5365dSVijay Mahadevan for (j = ise[2]; j <= ise[3]; j++) { 237c8d5365dSVijay Mahadevan for (i = ise[0]; i <= ise[1]; i++, vcount++) { 238*aa859218SVijay Mahadevan DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 239c8d5365dSVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 240c8d5365dSVijay Mahadevan } 241c8d5365dSVijay Mahadevan } 242c8d5365dSVijay Mahadevan } 243c8d5365dSVijay Mahadevan 244e5136372SVijay Mahadevan /* create ghosted vertices requested by user - below the current plane */ 245e5136372SVijay Mahadevan if (ise[2*dim-2] > 0) { 246e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) { 247e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) { 248e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) { 249*aa859218SVijay Mahadevan DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 250e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 251e5136372SVijay Mahadevan } 252e5136372SVijay Mahadevan } 253e5136372SVijay Mahadevan } 254e5136372SVijay Mahadevan } 255e5136372SVijay Mahadevan 256e5136372SVijay Mahadevan /* create ghosted vertices requested by user - above the current plane */ 257e5136372SVijay Mahadevan if (ise[2*dim-1] < nele) { 258e5136372SVijay Mahadevan for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) { 259e5136372SVijay Mahadevan for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) { 260e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) { 261*aa859218SVijay Mahadevan DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 262e5136372SVijay Mahadevan vgid[vcount] = (k*npts+j)*npts+i+1; 263e5136372SVijay Mahadevan } 26451d15aeeSVijay Mahadevan } 26551d15aeeSVijay Mahadevan } 26651d15aeeSVijay Mahadevan } 26751d15aeeSVijay Mahadevan 268c8d5365dSVijay Mahadevan if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount); 269c8d5365dSVijay Mahadevan 27051d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 27151d15aeeSVijay Mahadevan 272a4717116SVijay Mahadevan /* The global ID tag is applied to each owned 273a4717116SVijay Mahadevan vertex. It acts as an global identifier which MOAB uses to 274a4717116SVijay Mahadevan assemble the individual pieces of the mesh: 275a4717116SVijay Mahadevan Set the global ID indices */ 27651d15aeeSVijay Mahadevan merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr); 27751d15aeeSVijay Mahadevan 278a4717116SVijay Mahadevan /* Create elements between mesh points using the ReadUtilInterface 279a4717116SVijay Mahadevan get the reference to element connectivities for all local elements from the ReadUtilInterface */ 280e5136372SVijay Mahadevan merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr); 28151d15aeeSVijay Mahadevan 282a4717116SVijay Mahadevan /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */ 283c8d5365dSVijay Mahadevan // PetscPrintf(PETSC_COMM_SELF, "[%D] first local handle %D\n", rank, vfirst); 284e5136372SVijay Mahadevan vfirst-=vgid[0]-1; 28551d15aeeSVijay Mahadevan 286a4717116SVijay Mahadevan /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k): 287a4717116SVijay Mahadevan and then set the connectivity for each element appropriately */ 28851d15aeeSVijay Mahadevan int ecount=0; 28951d15aeeSVijay Mahadevan 290e5136372SVijay Mahadevan /* create ghosted elements requested by user - below the current plane */ 291e5136372SVijay Mahadevan if (ise[2*dim-2] >= nghost) { 292e5136372SVijay Mahadevan for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) { 293e5136372SVijay Mahadevan for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) { 294e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++, ecount++) { 295*aa859218SVijay Mahadevan DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 29651d15aeeSVijay Mahadevan } 29751d15aeeSVijay Mahadevan } 29851d15aeeSVijay Mahadevan } 29951d15aeeSVijay Mahadevan } 300e5136372SVijay Mahadevan 301e5136372SVijay Mahadevan /* create owned elements requested by user */ 302e5136372SVijay Mahadevan for (k = ise[4]; k < std::max(ise[5],1); k++) { 303e5136372SVijay Mahadevan for (j = ise[2]; j < std::max(ise[3],1); j++) { 304e5136372SVijay Mahadevan for (i = ise[0]; i < std::max(ise[1],1); i++,ecount++) { 305*aa859218SVijay Mahadevan DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 306e5136372SVijay Mahadevan } 307e5136372SVijay Mahadevan } 308e5136372SVijay Mahadevan } 309e5136372SVijay Mahadevan 310e5136372SVijay Mahadevan /* create ghosted elements requested by user - above the current plane */ 311e5136372SVijay Mahadevan if (ise[2*dim-1] <= nele-nghost) { 312e5136372SVijay Mahadevan for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) { 313e5136372SVijay Mahadevan for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) { 314e5136372SVijay Mahadevan for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++, ecount++) { 315*aa859218SVijay Mahadevan DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 316e5136372SVijay Mahadevan } 317e5136372SVijay Mahadevan } 318e5136372SVijay Mahadevan } 319e5136372SVijay Mahadevan } 320e5136372SVijay Mahadevan 321e5136372SVijay Mahadevan merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr); 32251d15aeeSVijay Mahadevan 323a4717116SVijay 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. 324a4717116SVijay Mahadevan first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */ 32551d15aeeSVijay Mahadevan merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr); 326a4717116SVijay Mahadevan 327e5136372SVijay Mahadevan if (locnele+ghnele != (int) ownedelms.size()) 328e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size()); 329e5136372SVijay Mahadevan else if(locnpts+ghnpts != (int) ownedvtx.size()) 330e5136372SVijay Mahadevan SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size()); 33151d15aeeSVijay Mahadevan else 332e427d9c9SVijay Mahadevan PetscInfo2(NULL, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size()); 33351d15aeeSVijay Mahadevan 3343a4aeca1SVijay Mahadevan /* lets create some sets */ 3353a4aeca1SVijay 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); 3363a4aeca1SVijay Mahadevan 3373a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr); 3383a4aeca1SVijay Mahadevan merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr); 3393a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, ®ionset, 1, &dmmoab->dim);MBERRNM(merr); 340b5410836SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr); 3413a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr); 3423a4aeca1SVijay Mahadevan 3433a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr); 3443a4aeca1SVijay Mahadevan merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr); 3453a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr); 3463a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr); 3473a4aeca1SVijay Mahadevan 3483a4aeca1SVijay Mahadevan if (build_adjacencies) { 3493a4aeca1SVijay Mahadevan // generate all lower dimensional adjacencies 3503a4aeca1SVijay Mahadevan merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr); 3513a4aeca1SVijay Mahadevan merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr); 3523a4aeca1SVijay Mahadevan adj.merge(dim2); 3533a4aeca1SVijay Mahadevan 3543a4aeca1SVijay Mahadevan /* create face sets */ 3553a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr); 3563a4aeca1SVijay Mahadevan merr = mbiface->add_entities(faceset, adj);MBERRNM(merr); 3573a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr); 3583a4aeca1SVijay Mahadevan i=dim-1; 3593a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr); 3603a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr); 3613a4aeca1SVijay Mahadevan PetscInfo2(dm, "Found %d %d-Dim quantities.\n", adj.size(), dim-1); 3623a4aeca1SVijay Mahadevan 3633a4aeca1SVijay Mahadevan if (dim > 2) { 3643a4aeca1SVijay Mahadevan dim2.clear(); 3653a4aeca1SVijay Mahadevan /* create edge sets, if appropriate i.e., if dim=3 */ 3663a4aeca1SVijay Mahadevan merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr); 3673a4aeca1SVijay Mahadevan merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr); 3683a4aeca1SVijay Mahadevan merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr); 3693a4aeca1SVijay Mahadevan merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr); 3703a4aeca1SVijay Mahadevan i=dim-2; 3713a4aeca1SVijay Mahadevan merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr); 3723a4aeca1SVijay Mahadevan merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr); 3733a4aeca1SVijay Mahadevan PetscInfo2(dm, "Found %d %d-Dim quantities.\n", adj.size(), dim-2); 3743a4aeca1SVijay Mahadevan } 3753a4aeca1SVijay Mahadevan } 3763a4aeca1SVijay Mahadevan 377a4717116SVijay Mahadevan /* check the handles */ 378a4717116SVijay Mahadevan merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr); 37951d15aeeSVijay Mahadevan 380c8d5365dSVijay Mahadevan #if 0 381e5136372SVijay Mahadevan std::stringstream sstr; 382e5136372SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 383e5136372SVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 384e5136372SVijay Mahadevan #endif 385e5136372SVijay Mahadevan 386a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */ 387a4717116SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr); 388ce1fd009SVijay Mahadevan merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,NULL,&id_tag);MBERRV(mbiface,merr); 38951d15aeeSVijay Mahadevan 3903a4aeca1SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr); 391c8d5365dSVijay Mahadevan 392e427d9c9SVijay Mahadevan /* Reassign global IDs on all entities. */ 393e427d9c9SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true,false);MBERRNM(merr); 394e427d9c9SVijay Mahadevan 395a4717116SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 396a4717116SVijay Mahadevan on all of the ghost vertexes */ 397c8d5365dSVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 398c8d5365dSVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr); 399c8d5365dSVijay Mahadevan 400a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr); 401a4717116SVijay Mahadevan merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr); 40251d15aeeSVijay Mahadevan PetscFunctionReturn(0); 40351d15aeeSVijay Mahadevan } 40451d15aeeSVijay Mahadevan 40551d15aeeSVijay Mahadevan 406a4717116SVijay Mahadevan #undef __FUNCT__ 407a4717116SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private" 4087023aa44SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* extra_opts, const char** read_opts) 40951d15aeeSVijay Mahadevan { 410a4717116SVijay Mahadevan std::ostringstream str; 41151d15aeeSVijay Mahadevan 412a4717116SVijay Mahadevan PetscFunctionBegin; 413e23c60ebSVijay Mahadevan /* do parallel read unless using only one processor */ 414a4717116SVijay Mahadevan if (numproc > 1) { 4157023aa44SVijay Mahadevan str << "PARALLEL=" << mode << ";PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"; 4167023aa44SVijay Mahadevan str << "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;"; 417a4717116SVijay Mahadevan if (by_rank) 418a4717116SVijay Mahadevan str << "PARTITION_BY_RANK;"; 419a4717116SVijay Mahadevan } 42051d15aeeSVijay Mahadevan 421c8d5365dSVijay Mahadevan if (dbglevel) { 422c8d5365dSVijay Mahadevan str << "CPUTIME;DEBUG_IO=" << dbglevel << ";"; 423c8d5365dSVijay Mahadevan if (numproc>1) 424c8d5365dSVijay Mahadevan str << "DEBUG_PIO=" << dbglevel << ";"; 425c8d5365dSVijay Mahadevan } 42651d15aeeSVijay Mahadevan 427c8d5365dSVijay Mahadevan if (extra_opts) 428c8d5365dSVijay Mahadevan str << extra_opts; 42951d15aeeSVijay Mahadevan 4307023aa44SVijay Mahadevan *read_opts = str.str().c_str(); 43151d15aeeSVijay Mahadevan PetscFunctionReturn(0); 43251d15aeeSVijay Mahadevan } 43351d15aeeSVijay Mahadevan 43451d15aeeSVijay Mahadevan 43551d15aeeSVijay Mahadevan #undef __FUNCT__ 436a4717116SVijay Mahadevan #define __FUNCT__ "DMMoabLoadFromFile" 437*aa859218SVijay Mahadevan /*@ 438*aa859218SVijay Mahadevan DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file. 439*aa859218SVijay Mahadevan 440*aa859218SVijay Mahadevan Collective on MPI_Comm 441*aa859218SVijay Mahadevan 442*aa859218SVijay Mahadevan Input Parameters: 443*aa859218SVijay Mahadevan + comm - The communicator for the DM object 444*aa859218SVijay Mahadevan . dim - The spatial dimension 445*aa859218SVijay Mahadevan - filename - The name of the mesh file to be loaded 446*aa859218SVijay Mahadevan - usrreadopts - The options string to read a MOAB mesh. 447*aa859218SVijay Mahadevan Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo) 448*aa859218SVijay Mahadevan 449*aa859218SVijay Mahadevan Output Parameter: 450*aa859218SVijay Mahadevan . dm - The DM object 451*aa859218SVijay Mahadevan 452*aa859218SVijay Mahadevan Level: beginner 453*aa859218SVijay Mahadevan 454*aa859218SVijay Mahadevan .keywords: DM, create 455*aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh() 456*aa859218SVijay Mahadevan @*/ 457a4717116SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm) 45851d15aeeSVijay Mahadevan { 459a4717116SVijay Mahadevan moab::ErrorCode merr; 4607023aa44SVijay Mahadevan PetscInt nprocs,dbglevel=0; 4617023aa44SVijay Mahadevan PetscBool part_by_rank=PETSC_FALSE; 462a4717116SVijay Mahadevan DM_Moab *dmmoab; 463a4717116SVijay Mahadevan moab::Interface *mbiface; 464a4717116SVijay Mahadevan moab::ParallelComm *pcomm; 465a4717116SVijay Mahadevan moab::Range verts,elems; 4667023aa44SVijay Mahadevan const char *readopts; 467a4717116SVijay Mahadevan PetscErrorCode ierr; 46851d15aeeSVijay Mahadevan 46951d15aeeSVijay Mahadevan PetscFunctionBegin; 470a4717116SVijay Mahadevan PetscValidPointer(dm,5); 47151d15aeeSVijay Mahadevan 472a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 473a4717116SVijay Mahadevan ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 47451d15aeeSVijay Mahadevan 475a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 476a4717116SVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 477a4717116SVijay Mahadevan mbiface = dmmoab->mbiface; 478a4717116SVijay Mahadevan pcomm = dmmoab->pcomm; 479a4717116SVijay Mahadevan nprocs = pcomm->size(); 480*aa859218SVijay Mahadevan /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */ 481c8d5365dSVijay Mahadevan dmmoab->dim = dim; 482a4717116SVijay Mahadevan 483b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */ 484b5410836SVijay Mahadevan merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 485b5410836SVijay Mahadevan 4867023aa44SVijay Mahadevan /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */ 4877023aa44SVijay Mahadevan ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab"); 4887023aa44SVijay Mahadevan ierr = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr); 4897023aa44SVijay 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); 4907023aa44SVijay Mahadevan ierr = PetscOptionsEnd(); 4917023aa44SVijay Mahadevan 492a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */ 4937023aa44SVijay Mahadevan ierr = DMMoab_GetReadOptions_Private(part_by_rank, nprocs, dim, MOAB_PARROPTS_READ_PART, dbglevel, usrreadopts, &readopts);CHKERRQ(ierr); 4947023aa44SVijay Mahadevan 495e5136372SVijay Mahadevan PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts); 496a4717116SVijay Mahadevan 497a4717116SVijay Mahadevan /* Load the mesh from a file. */ 4987023aa44SVijay Mahadevan merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr); 499a4717116SVijay Mahadevan 5006e40195eSVijay Mahadevan /* Reassign global IDs on all entities. */ 501e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 502e5136372SVijay Mahadevan 503e5136372SVijay Mahadevan /* load the local vertices */ 504e5136372SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 505e5136372SVijay Mahadevan /* load the local elements */ 506e5136372SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 507e5136372SVijay Mahadevan 508e5136372SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 509e5136372SVijay Mahadevan on all of the ghost vertexes */ 510e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr); 511e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr); 512e5136372SVijay Mahadevan 513e5136372SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr); 5146e40195eSVijay Mahadevan 515a4717116SVijay Mahadevan merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 516a4717116SVijay Mahadevan 517e5136372SVijay Mahadevan PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size()); 518e5136372SVijay Mahadevan 519e427d9c9SVijay Mahadevan #if 0 520e5136372SVijay Mahadevan std::stringstream sstr; 521e5136372SVijay Mahadevan sstr << "test_" << pcomm->rank() << ".vtk"; 522e5136372SVijay Mahadevan mbiface->write_mesh(sstr.str().c_str()); 523e5136372SVijay Mahadevan #endif 52451d15aeeSVijay Mahadevan PetscFunctionReturn(0); 52551d15aeeSVijay Mahadevan } 52651d15aeeSVijay Mahadevan 527