1 #include <petsc-private/dmmbimpl.h> /*I "petscdm.h" I*/ 2 #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/ 3 4 #include <petscdmmoab.h> 5 #include <MBTagConventions.hpp> 6 #include <moab/ReadUtilIface.hpp> 7 #include <moab/ScdInterface.hpp> 8 #include <moab/CN.hpp> 9 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "DMMoabComputeDomainBounds_Private" 13 static PetscErrorCode DMMoabComputeDomainBounds_Private(moab::ParallelComm* pcomm, PetscInt dim, PetscInt neleglob, PetscInt *ise) 14 { 15 PetscInt size,rank; 16 PetscInt fraction,remainder; 17 PetscInt starts[3],sizes[3]; 18 19 PetscFunctionBegin; 20 if(dim<1 && dim>3) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The problem dimension is invalid: %D",dim); 21 22 size=pcomm->size(); 23 rank=pcomm->rank(); 24 fraction=neleglob/size; /* partition only by the largest dimension */ 25 remainder=neleglob%size; /* remainder after partition which gets evenly distributed by round-robin */ 26 27 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); 28 29 starts[0]=starts[1]=starts[2]=0; /* default dimensional element = 1 */ 30 sizes[0]=sizes[1]=sizes[2]=neleglob; /* default dimensional element = 1 */ 31 32 if(rank < remainder) { 33 /* This process gets "fraction+1" elements */ 34 sizes[dim-1] = (fraction + 1); 35 starts[dim-1] = rank * (fraction+1); 36 } else { 37 /* This process gets "fraction" elements */ 38 sizes[dim-1] = fraction; 39 starts[dim-1] = (remainder*(fraction+1) + fraction*(rank-remainder)); 40 } 41 42 for(int i=dim-1; i>=0; --i) { 43 ise[2*i]=starts[i];ise[2*i+1]=starts[i]+sizes[i]; 44 } 45 46 PetscFunctionReturn(0); 47 } 48 49 #undef __FUNCT__ 50 #define __FUNCT__ "DMMoab_SetStructuredCoords_Private" 51 static void DMMoab_SetStructuredCoords_Private(PetscInt i, PetscInt j, PetscInt k, PetscReal hx, PetscReal hy, PetscReal hz, PetscInt vcount, std::vector<double*>& vcoords) 52 { 53 vcoords[0][vcount] = i*hx; 54 vcoords[1][vcount] = j*hy; 55 vcoords[2][vcount] = k*hz; 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "DMMoab_SetElementConnectivity_Private" 60 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) 61 { 62 std::vector<int> subent_conn(pow(2,dim)); /* only linear edge, quad, hex supported now */ 63 64 moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data()); 65 66 switch(dim) { 67 case 1: 68 connectivity[offset+subent_conn[0]] = vfirst+i; 69 connectivity[offset+subent_conn[1]] = vfirst+(i+1); 70 break; 71 case 2: 72 connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1); 73 connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1); 74 connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1); 75 connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1); 76 break; 77 case 3: 78 default: 79 connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k); 80 connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k); 81 connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k); 82 connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k); 83 connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1)); 84 connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1)); 85 connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1)); 86 connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1)); 87 break; 88 } 89 } 90 91 #undef __FUNCT__ 92 #define __FUNCT__ "DMMoabCreateBoxMesh" 93 /*@ 94 DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with user specified bounds. 95 96 Collective on MPI_Comm 97 98 Input Parameters: 99 + comm - The communicator for the DM object 100 . dim - The spatial dimension 101 - 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 - nele - The number of discrete elements in each direction 103 - user_nghost - The number of ghosted layers needed in the partitioned mesh 104 105 Output Parameter: 106 . dm - The DM object 107 108 Level: beginner 109 110 .keywords: DM, create 111 .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile() 112 @*/ 113 PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, const PetscReal* bounds, PetscInt nele, PetscInt user_nghost, DM *dm) 114 { 115 PetscErrorCode ierr; 116 moab::ErrorCode merr; 117 PetscInt i,j,k,n,nprocs; 118 DM_Moab *dmmoab; 119 moab::Interface *mbiface; 120 moab::ParallelComm *pcomm; 121 moab::ReadUtilIface* readMeshIface; 122 123 moab::Tag id_tag=PETSC_NULL; 124 moab::Range ownedvtx,ownedelms; 125 moab::EntityHandle vfirst,efirst,regionset,faceset,edgeset,vtxset; 126 std::vector<double*> vcoords; 127 moab::EntityHandle *connectivity = 0; 128 moab::EntityType etype; 129 PetscInt ise[6]; 130 PetscReal xse[6],defbounds[6]; 131 /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */ 132 const PetscInt nghost=0; 133 134 moab::Tag geom_tag; 135 136 moab::Range adj,dim3,dim2; 137 bool build_adjacencies=false; 138 139 const PetscInt npts=nele+1; /* Number of points in every dimension */ 140 PetscInt vpere,locnele,locnpts,ghnele,ghnpts; /* Number of verts/element, vertices, elements owned by this process */ 141 142 PetscFunctionBegin; 143 if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n"); 144 145 ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr); 146 /* total number of vertices in all dimensions */ 147 n=pow(npts,dim); 148 149 /* do some error checking */ 150 if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n"); 151 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"); 152 if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n"); 153 154 /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 155 ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 156 157 /* get all the necessary handles from the private DM object */ 158 dmmoab = (DM_Moab*)(*dm)->data; 159 mbiface = dmmoab->mbiface; 160 pcomm = dmmoab->pcomm; 161 id_tag = dmmoab->ltog_tag; 162 nprocs = pcomm->size(); 163 dmmoab->dim = dim; 164 165 /* create a file set to associate all entities in current mesh */ 166 merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 167 168 /* No errors yet; proceed with building the mesh */ 169 merr = mbiface->query_interface(readMeshIface);MBERRNM(merr); 170 171 ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr); 172 173 /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */ 174 ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr); 175 176 /* set some variables based on dimension */ 177 switch(dim) { 178 case 1: 179 vpere = 2; 180 locnele = (ise[1]-ise[0]); 181 locnpts = (ise[1]-ise[0]+1); 182 ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 ); 183 ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0); 184 etype = moab::MBEDGE; 185 break; 186 case 2: 187 vpere = 4; 188 locnele = (ise[1]-ise[0])*(ise[3]-ise[2]); 189 locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1); 190 ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0); 191 ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0); 192 etype = moab::MBQUAD; 193 break; 194 case 3: 195 vpere = 8; 196 locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]); 197 locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1); 198 ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0); 199 ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0); 200 etype = moab::MBHEX; 201 break; 202 } 203 204 /* we have a domain of size [1,1,1] - now compute local co-ordinate box */ 205 ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr); 206 for(i=0; i<6; ++i) { 207 xse[i]=(PetscReal)ise[i]/nele; 208 } 209 210 /* Create vertexes and set the coodinate of each vertex */ 211 merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr); 212 213 /* Compute the co-ordinates of vertices and global IDs */ 214 std::vector<int> vgid(locnpts+ghnpts); 215 int vcount=0; 216 217 if (!bounds) { /* default box mesh is defined on a unit-cube */ 218 defbounds[0]=0.0; defbounds[1]=1.0; 219 defbounds[2]=0.0; defbounds[3]=1.0; 220 defbounds[4]=0.0; defbounds[5]=1.0; 221 bounds=defbounds; 222 } 223 else { 224 /* validate the bounds data */ 225 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]); 226 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]); 227 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]); 228 } 229 230 const double hx=(bounds[1]-bounds[0])/nele; 231 const double hy=(dim > 1 ? (bounds[3]-bounds[2])/nele : 0.0); 232 const double hz=(dim > 2 ? (bounds[5]-bounds[4])/nele : 0.0); 233 234 /* create all the owned vertices */ 235 for (k = ise[4]; k <= ise[5]; k++) { 236 for (j = ise[2]; j <= ise[3]; j++) { 237 for (i = ise[0]; i <= ise[1]; i++, vcount++) { 238 DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 239 vgid[vcount] = (k*npts+j)*npts+i+1; 240 } 241 } 242 } 243 244 /* create ghosted vertices requested by user - below the current plane */ 245 if (ise[2*dim-2] > 0) { 246 for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) { 247 for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) { 248 for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) { 249 DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 250 vgid[vcount] = (k*npts+j)*npts+i+1; 251 } 252 } 253 } 254 } 255 256 /* create ghosted vertices requested by user - above the current plane */ 257 if (ise[2*dim-1] < nele) { 258 for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) { 259 for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) { 260 for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) { 261 DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords); 262 vgid[vcount] = (k*npts+j)*npts+i+1; 263 } 264 } 265 } 266 } 267 268 if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount); 269 270 merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 271 272 /* The global ID tag is applied to each owned 273 vertex. It acts as an global identifier which MOAB uses to 274 assemble the individual pieces of the mesh: 275 Set the global ID indices */ 276 merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr); 277 278 /* Create elements between mesh points using the ReadUtilInterface 279 get the reference to element connectivities for all local elements from the ReadUtilInterface */ 280 merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr); 281 282 /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */ 283 // PetscPrintf(PETSC_COMM_SELF, "[%D] first local handle %D\n", rank, vfirst); 284 vfirst-=vgid[0]-1; 285 286 /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k): 287 and then set the connectivity for each element appropriately */ 288 int ecount=0; 289 290 /* create ghosted elements requested by user - below the current plane */ 291 if (ise[2*dim-2] >= nghost) { 292 for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) { 293 for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) { 294 for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++, ecount++) { 295 DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 296 } 297 } 298 } 299 } 300 301 /* create owned elements requested by user */ 302 for (k = ise[4]; k < std::max(ise[5],1); k++) { 303 for (j = ise[2]; j < std::max(ise[3],1); j++) { 304 for (i = ise[0]; i < std::max(ise[1],1); i++,ecount++) { 305 DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 306 } 307 } 308 } 309 310 /* create ghosted elements requested by user - above the current plane */ 311 if (ise[2*dim-1] <= nele-nghost) { 312 for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) { 313 for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) { 314 for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++, ecount++) { 315 DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity); 316 } 317 } 318 } 319 } 320 321 merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr); 322 323 /* 2. Get the vertices and hexes from moab and check their numbers against I*J*K and (I-1)*(J-1)*(K-1), resp. 324 first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */ 325 merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr); 326 327 if (locnele+ghnele != (int) ownedelms.size()) 328 SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size()); 329 else if(locnpts+ghnpts != (int) ownedvtx.size()) 330 SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size()); 331 else 332 PetscInfo2(NULL, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size()); 333 334 /* lets create some sets */ 335 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); 336 337 merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr); 338 merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr); 339 merr = mbiface->tag_set_data(geom_tag, ®ionset, 1, &dmmoab->dim);MBERRNM(merr); 340 merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr); 341 merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr); 342 343 merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr); 344 merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr); 345 merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr); 346 merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr); 347 348 if (build_adjacencies) { 349 // generate all lower dimensional adjacencies 350 merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr); 351 merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr); 352 adj.merge(dim2); 353 354 /* create face sets */ 355 merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr); 356 merr = mbiface->add_entities(faceset, adj);MBERRNM(merr); 357 merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr); 358 i=dim-1; 359 merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr); 360 merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr); 361 PetscInfo2(dm, "Found %d %d-Dim quantities.\n", adj.size(), dim-1); 362 363 if (dim > 2) { 364 dim2.clear(); 365 /* create edge sets, if appropriate i.e., if dim=3 */ 366 merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr); 367 merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr); 368 merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr); 369 merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr); 370 i=dim-2; 371 merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr); 372 merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr); 373 PetscInfo2(dm, "Found %d %d-Dim quantities.\n", adj.size(), dim-2); 374 } 375 } 376 377 /* check the handles */ 378 merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr); 379 380 #if 0 381 std::stringstream sstr; 382 sstr << "test_" << pcomm->rank() << ".vtk"; 383 mbiface->write_mesh(sstr.str().c_str()); 384 #endif 385 386 /* resolve the shared entities by exchanging information to adjacent processors */ 387 merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr); 388 merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,NULL,&id_tag);MBERRV(mbiface,merr); 389 390 merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr); 391 392 /* Reassign global IDs on all entities. */ 393 merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true,false);MBERRNM(merr); 394 395 /* Everything is set up, now just do a tag exchange to update tags 396 on all of the ghost vertexes */ 397 merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr); 398 merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr); 399 400 merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr); 401 merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr); 402 PetscFunctionReturn(0); 403 } 404 405 406 #undef __FUNCT__ 407 #define __FUNCT__ "DMMoab_GetReadOptions_Private" 408 PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, MoabReadMode mode, PetscInt dbglevel, const char* extra_opts, const char** read_opts) 409 { 410 std::ostringstream str; 411 412 PetscFunctionBegin; 413 /* do parallel read unless using only one processor */ 414 if (numproc > 1) { 415 str << "PARALLEL=" << mode << ";PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;"; 416 str << "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;"; 417 if (by_rank) 418 str << "PARTITION_BY_RANK;"; 419 } 420 421 if (dbglevel) { 422 str << "CPUTIME;DEBUG_IO=" << dbglevel << ";"; 423 if (numproc>1) 424 str << "DEBUG_PIO=" << dbglevel << ";"; 425 } 426 427 if (extra_opts) 428 str << extra_opts; 429 430 *read_opts = str.str().c_str(); 431 PetscFunctionReturn(0); 432 } 433 434 435 #undef __FUNCT__ 436 #define __FUNCT__ "DMMoabLoadFromFile" 437 /*@ 438 DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file. 439 440 Collective on MPI_Comm 441 442 Input Parameters: 443 + comm - The communicator for the DM object 444 . dim - The spatial dimension 445 - filename - The name of the mesh file to be loaded 446 - usrreadopts - The options string to read a MOAB mesh. 447 Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo) 448 449 Output Parameter: 450 . dm - The DM object 451 452 Level: beginner 453 454 .keywords: DM, create 455 .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh() 456 @*/ 457 PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm) 458 { 459 moab::ErrorCode merr; 460 PetscInt nprocs,dbglevel=0; 461 PetscBool part_by_rank=PETSC_FALSE; 462 DM_Moab *dmmoab; 463 moab::Interface *mbiface; 464 moab::ParallelComm *pcomm; 465 moab::Range verts,elems; 466 const char *readopts; 467 PetscErrorCode ierr; 468 469 PetscFunctionBegin; 470 PetscValidPointer(dm,5); 471 472 /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 473 ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr); 474 475 /* get all the necessary handles from the private DM object */ 476 dmmoab = (DM_Moab*)(*dm)->data; 477 mbiface = dmmoab->mbiface; 478 pcomm = dmmoab->pcomm; 479 nprocs = pcomm->size(); 480 /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */ 481 dmmoab->dim = dim; 482 483 /* create a file set to associate all entities in current mesh */ 484 merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 485 486 /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */ 487 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab"); 488 ierr = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr); 489 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); 490 ierr = PetscOptionsEnd(); 491 492 /* add mesh loading options specific to the DM */ 493 ierr = DMMoab_GetReadOptions_Private(part_by_rank, nprocs, dim, MOAB_PARROPTS_READ_PART, dbglevel, usrreadopts, &readopts);CHKERRQ(ierr); 494 495 PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts); 496 497 /* Load the mesh from a file. */ 498 merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr); 499 500 /* Reassign global IDs on all entities. */ 501 merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr); 502 503 /* load the local vertices */ 504 merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 505 /* load the local elements */ 506 merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 507 508 /* Everything is set up, now just do a tag exchange to update tags 509 on all of the ghost vertexes */ 510 merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr); 511 merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr); 512 513 merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr); 514 515 merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 516 517 PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size()); 518 519 #if 0 520 std::stringstream sstr; 521 sstr << "test_" << pcomm->rank() << ".vtk"; 522 mbiface->write_mesh(sstr.str().c_str()); 523 #endif 524 PetscFunctionReturn(0); 525 } 526 527