1 #include <petsc/private/dmmbimpl.h> /*I "petscdmmoab.h" I*/ 2 #include <petsc/private/vecimpl.h> 3 4 #include <petscdmmoab.h> 5 #include <MBTagConventions.hpp> 6 #include <moab/ReadUtilIface.hpp> 7 #include <moab/MergeMesh.hpp> 8 #include <moab/CN.hpp> 9 10 typedef struct { 11 // options 12 PetscInt A, B, C, M, N, K, dim; 13 PetscInt blockSizeVertexXYZ[3]; // Number of element blocks per partition 14 PetscInt blockSizeElementXYZ[3]; 15 PetscReal xyzbounds[6]; // the physical size of the domain 16 bool newMergeMethod, keep_skins, simplex, adjEnts; 17 18 // compute params 19 PetscReal dx, dy, dz; 20 PetscInt NX, NY, NZ, nex, ney, nez; 21 PetscInt q, xstride, ystride, zstride; 22 PetscBool usrxyzgrid, usrprocgrid, usrrefgrid; 23 PetscInt fraction, remainder, cumfraction; 24 PetscLogEvent generateMesh, generateElements, generateVertices, parResolve; 25 } DMMoabMeshGeneratorCtx; 26 27 PetscInt DMMoab_SetTensorElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt offset, PetscInt corner, std::vector<PetscInt>& subent_conn, moab::EntityHandle *connectivity) 28 { 29 switch (genCtx.dim) { 30 case 1: 31 subent_conn.resize(2); 32 moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data()); 33 connectivity[offset + subent_conn[0]] = corner; 34 connectivity[offset + subent_conn[1]] = corner + 1; 35 break; 36 case 2: 37 subent_conn.resize(4); 38 moab::CN::SubEntityVertexIndices(moab::MBQUAD, 2, 0, subent_conn.data()); 39 connectivity[offset + subent_conn[0]] = corner; 40 connectivity[offset + subent_conn[1]] = corner + 1; 41 connectivity[offset + subent_conn[2]] = corner + 1 + genCtx.ystride; 42 connectivity[offset + subent_conn[3]] = corner + genCtx.ystride; 43 break; 44 case 3: 45 default: 46 subent_conn.resize(8); 47 moab::CN::SubEntityVertexIndices(moab::MBHEX, 3, 0, subent_conn.data()); 48 connectivity[offset + subent_conn[0]] = corner; 49 connectivity[offset + subent_conn[1]] = corner + 1; 50 connectivity[offset + subent_conn[2]] = corner + 1 + genCtx.ystride; 51 connectivity[offset + subent_conn[3]] = corner + genCtx.ystride; 52 connectivity[offset + subent_conn[4]] = corner + genCtx.zstride; 53 connectivity[offset + subent_conn[5]] = corner + 1 + genCtx.zstride; 54 connectivity[offset + subent_conn[6]] = corner + 1 + genCtx.ystride + genCtx.zstride; 55 connectivity[offset + subent_conn[7]] = corner + genCtx.ystride + genCtx.zstride; 56 break; 57 } 58 return subent_conn.size(); 59 } 60 61 PetscInt DMMoab_SetSimplexElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt subelem, PetscInt offset, PetscInt corner, std::vector<PetscInt>& subent_conn, moab::EntityHandle *connectivity) 62 { 63 PetscInt A, B, C, D, E, F, G, H, M; 64 const PetscInt trigen_opts = 1; /* 1 - Aligned diagonally to right, 2 - Aligned diagonally to left, 3 - 4 elements per quad */ 65 A = corner; 66 B = corner + 1; 67 switch (genCtx.dim) { 68 case 1: 69 subent_conn.resize(2); /* only linear EDGE supported now */ 70 moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data()); 71 connectivity[offset + subent_conn[0]] = A; 72 connectivity[offset + subent_conn[1]] = B; 73 break; 74 case 2: 75 C = corner + 1 + genCtx.ystride; 76 D = corner + genCtx.ystride; 77 M = corner + 0.5; /* technically -- need to modify vertex generation */ 78 subent_conn.resize(3); /* only linear TRI supported */ 79 moab::CN::SubEntityVertexIndices(moab::MBTRI, 2, 0, subent_conn.data()); 80 if (trigen_opts == 1) { 81 if (subelem) { /* 0 1 2 of a QUAD */ 82 connectivity[offset + subent_conn[0]] = B; 83 connectivity[offset + subent_conn[1]] = C; 84 connectivity[offset + subent_conn[2]] = A; 85 } 86 else { /* 2 3 0 of a QUAD */ 87 connectivity[offset + subent_conn[0]] = D; 88 connectivity[offset + subent_conn[1]] = A; 89 connectivity[offset + subent_conn[2]] = C; 90 } 91 } 92 else if (trigen_opts == 2) { 93 if (subelem) { /* 0 1 2 of a QUAD */ 94 connectivity[offset + subent_conn[0]] = A; 95 connectivity[offset + subent_conn[1]] = B; 96 connectivity[offset + subent_conn[2]] = D; 97 } 98 else { /* 2 3 0 of a QUAD */ 99 connectivity[offset + subent_conn[0]] = C; 100 connectivity[offset + subent_conn[1]] = D; 101 connectivity[offset + subent_conn[2]] = B; 102 } 103 } 104 else { 105 switch (subelem) { /* 0 1 2 of a QUAD */ 106 case 0: 107 connectivity[offset + subent_conn[0]] = A; 108 connectivity[offset + subent_conn[1]] = B; 109 connectivity[offset + subent_conn[2]] = M; 110 break; 111 case 1: 112 connectivity[offset + subent_conn[0]] = B; 113 connectivity[offset + subent_conn[1]] = C; 114 connectivity[offset + subent_conn[2]] = M; 115 break; 116 case 2: 117 connectivity[offset + subent_conn[0]] = C; 118 connectivity[offset + subent_conn[1]] = D; 119 connectivity[offset + subent_conn[2]] = M; 120 break; 121 case 3: 122 connectivity[offset + subent_conn[0]] = D; 123 connectivity[offset + subent_conn[1]] = A; 124 connectivity[offset + subent_conn[2]] = M; 125 break; 126 } 127 } 128 break; 129 case 3: 130 default: 131 C = corner + 1 + genCtx.ystride; 132 D = corner + genCtx.ystride; 133 E = corner + genCtx.zstride; 134 F = corner + 1 + genCtx.zstride; 135 G = corner + 1 + genCtx.ystride + genCtx.zstride; 136 H = corner + genCtx.ystride + genCtx.zstride; 137 subent_conn.resize(4); /* only linear TET supported */ 138 moab::CN::SubEntityVertexIndices(moab::MBTET, 3, 0, subent_conn.data()); 139 switch (subelem) { 140 case 0: /* 4 3 7 6 of a HEX */ 141 connectivity[offset + subent_conn[0]] = E; 142 connectivity[offset + subent_conn[1]] = D; 143 connectivity[offset + subent_conn[2]] = H; 144 connectivity[offset + subent_conn[3]] = G; 145 break; 146 case 1: /* 0 1 2 5 of a HEX */ 147 connectivity[offset + subent_conn[0]] = A; 148 connectivity[offset + subent_conn[1]] = B; 149 connectivity[offset + subent_conn[2]] = C; 150 connectivity[offset + subent_conn[3]] = F; 151 break; 152 case 2: /* 0 3 4 5 of a HEX */ 153 connectivity[offset + subent_conn[0]] = A; 154 connectivity[offset + subent_conn[1]] = D; 155 connectivity[offset + subent_conn[2]] = E; 156 connectivity[offset + subent_conn[3]] = F; 157 break; 158 case 3: /* 2 6 3 5 of a HEX */ 159 connectivity[offset + subent_conn[0]] = C; 160 connectivity[offset + subent_conn[1]] = G; 161 connectivity[offset + subent_conn[2]] = D; 162 connectivity[offset + subent_conn[3]] = F; 163 break; 164 case 4: /* 0 2 3 5 of a HEX */ 165 connectivity[offset + subent_conn[0]] = A; 166 connectivity[offset + subent_conn[1]] = C; 167 connectivity[offset + subent_conn[2]] = D; 168 connectivity[offset + subent_conn[3]] = F; 169 break; 170 case 5: /* 3 6 4 5 of a HEX */ 171 connectivity[offset + subent_conn[0]] = D; 172 connectivity[offset + subent_conn[1]] = G; 173 connectivity[offset + subent_conn[2]] = E; 174 connectivity[offset + subent_conn[3]] = F; 175 break; 176 } 177 break; 178 } 179 return subent_conn.size(); 180 } 181 182 std::pair<PetscInt, PetscInt> DMMoab_SetElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt offset, PetscInt corner, moab::EntityHandle *connectivity) 183 { 184 PetscInt vcount = 0; 185 PetscInt simplices_per_tensor[4] = {0, 1, 2, 6}; 186 std::vector<PetscInt> subent_conn; /* only linear edge, tri, tet supported now */ 187 subent_conn.reserve(27); 188 PetscInt m, subelem; 189 if (genCtx.simplex) { 190 subelem = simplices_per_tensor[genCtx.dim]; 191 for (m = 0; m < subelem; m++) { 192 vcount = DMMoab_SetSimplexElementConnectivity_Private(genCtx, m, offset, corner, subent_conn, connectivity); 193 offset += vcount; 194 } 195 } 196 else { 197 subelem = 1; 198 vcount = DMMoab_SetTensorElementConnectivity_Private(genCtx, offset, corner, subent_conn, connectivity); 199 } 200 return std::pair<PetscInt, PetscInt>(vcount * subelem, subelem); 201 } 202 203 PetscErrorCode DMMoab_GenerateVertices_Private(moab::Interface *mbImpl, moab::ReadUtilIface *iface, DMMoabMeshGeneratorCtx& genCtx, PetscInt m, PetscInt n, PetscInt k, 204 PetscInt a, PetscInt b, PetscInt c, moab::Tag& global_id_tag, moab::EntityHandle& startv, moab::Range& uverts) 205 { 206 PetscInt x, y, z, ix, nnodes; 207 PetscInt ii, jj, kk; 208 std::vector<PetscReal*> arrays; 209 PetscInt* gids; 210 moab::ErrorCode merr; 211 212 PetscFunctionBegin; 213 /* we will generate (q*block+1)^3 vertices, and block^3 hexas; q is 1 for linear, 2 for quadratic 214 * the global id of the vertices will come from m, n, k, a, b, c 215 * x will vary from m*A*q*block + a*q*block to m*A*q*block+(a+1)*q*block etc. 216 */ 217 nnodes = genCtx.blockSizeVertexXYZ[0] * (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] * (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1) : 1); 218 PetscCall(PetscMalloc1(nnodes, &gids)); 219 220 merr = iface->get_node_coords(3, nnodes, 0, startv, arrays);MBERR("Can't get node coords.", merr); 221 222 /* will start with the lower corner: */ 223 /* x = ( m * genCtx.A + a) * genCtx.q * genCtx.blockSizeElementXYZ[0]; */ 224 /* y = ( n * genCtx.B + b) * genCtx.q * genCtx.blockSizeElementXYZ[1]; */ 225 /* z = ( k * genCtx.C + c) * genCtx.q * genCtx.blockSizeElementXYZ[2]; */ 226 227 x = ( m * genCtx.A + a) * genCtx.q; 228 y = ( n * genCtx.B + b) * genCtx.q; 229 z = ( k * genCtx.C + c) * genCtx.q; 230 PetscCall(PetscInfo(NULL, "Starting offset for coordinates := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", x, y, z)); 231 ix = 0; 232 moab::Range verts(startv, startv + nnodes - 1); 233 for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1); kk++) { 234 for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] : 1); jj++) { 235 for (ii = 0; ii < genCtx.blockSizeVertexXYZ[0]; ii++, ix++) { 236 /* set coordinates for the vertices */ 237 arrays[0][ix] = (x + ii) * genCtx.dx + genCtx.xyzbounds[0]; 238 arrays[1][ix] = (y + jj) * genCtx.dy + genCtx.xyzbounds[2]; 239 arrays[2][ix] = (z + kk) * genCtx.dz + genCtx.xyzbounds[4]; 240 PetscCall(PetscInfo(NULL, "Creating vertex with coordinates := %f, %f, %f\n", arrays[0][ix], arrays[1][ix], arrays[2][ix])); 241 242 /* If we want to set some tags on the vertices -> use the following entity handle definition: 243 moab::EntityHandle v = startv + ix; 244 */ 245 /* compute the global ID for vertex */ 246 gids[ix] = 1 + (x + ii) + (y + jj) * genCtx.NX + (z + kk) * (genCtx.NX * genCtx.NY); 247 } 248 } 249 } 250 /* set global ID data on vertices */ 251 mbImpl->tag_set_data(global_id_tag, verts, &gids[0]); 252 verts.swap(uverts); 253 PetscCall(PetscFree(gids)); 254 PetscFunctionReturn(0); 255 } 256 257 PetscErrorCode DMMoab_GenerateElements_Private(moab::Interface* mbImpl, moab::ReadUtilIface* iface, DMMoabMeshGeneratorCtx& genCtx, PetscInt m, PetscInt n, PetscInt k, 258 PetscInt a, PetscInt b, PetscInt c, moab::Tag& global_id_tag, moab::EntityHandle startv, moab::Range& cells) 259 { 260 moab::ErrorCode merr; 261 PetscInt ix, ie, xe, ye, ze; 262 PetscInt ii, jj, kk, nvperelem; 263 PetscInt simplices_per_tensor[4] = {0, 1, 2, 6}; 264 PetscInt ntensorelems = genCtx.blockSizeElementXYZ[0] * (genCtx.dim > 1 ? genCtx.blockSizeElementXYZ[1] * (genCtx.dim > 2 ? genCtx.blockSizeElementXYZ[2] : 1) : 1); /*pow(genCtx.blockSizeElement,genCtx.dim);*/ 265 PetscInt nelems = ntensorelems; 266 moab::EntityHandle starte; /* connectivity */ 267 moab::EntityHandle* conn; 268 269 PetscFunctionBegin; 270 switch (genCtx.dim) { 271 case 1: 272 nvperelem = 2; 273 merr = iface->get_element_connect(nelems, 2, moab::MBEDGE, 0, starte, conn);MBERR("Can't get EDGE2 element connectivity.", merr); 274 break; 275 case 2: 276 if (genCtx.simplex) { 277 nvperelem = 3; 278 nelems = ntensorelems * simplices_per_tensor[genCtx.dim]; 279 merr = iface->get_element_connect(nelems, 3, moab::MBTRI, 0, starte, conn);MBERR("Can't get TRI3 element connectivity.", merr); 280 } 281 else { 282 nvperelem = 4; 283 merr = iface->get_element_connect(nelems, 4, moab::MBQUAD, 0, starte, conn);MBERR("Can't get QUAD4 element connectivity.", merr); 284 } 285 break; 286 case 3: 287 default: 288 if (genCtx.simplex) { 289 nvperelem = 4; 290 nelems = ntensorelems * simplices_per_tensor[genCtx.dim]; 291 merr = iface->get_element_connect(nelems, 4, moab::MBTET, 0, starte, conn);MBERR("Can't get TET4 element connectivity.", merr); 292 } 293 else { 294 nvperelem = 8; 295 merr = iface->get_element_connect(nelems, 8, moab::MBHEX, 0, starte, conn);MBERR("Can't get HEX8 element connectivity.", merr); 296 } 297 break; 298 } 299 300 ix = ie = 0; /* index now in the elements, for global ids */ 301 302 /* create a temporary range to store local element handles */ 303 moab::Range tmp(starte, starte + nelems - 1); 304 std::vector<PetscInt> gids(nelems); 305 306 /* identify the elements at the lower corner, for their global ids */ 307 xe = m * genCtx.A * genCtx.blockSizeElementXYZ[0] + a * genCtx.blockSizeElementXYZ[0]; 308 ye = (genCtx.dim > 1 ? n * genCtx.B * genCtx.blockSizeElementXYZ[1] + b * genCtx.blockSizeElementXYZ[1] : 0); 309 ze = (genCtx.dim > 2 ? k * genCtx.C * genCtx.blockSizeElementXYZ[2] + c * genCtx.blockSizeElementXYZ[2] : 0); 310 311 /* create owned elements requested by genCtx */ 312 for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeElementXYZ[2] : 1); kk++) { 313 for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeElementXYZ[1] : 1); jj++) { 314 for (ii = 0; ii < genCtx.blockSizeElementXYZ[0]; ii++) { 315 316 moab::EntityHandle corner = startv + genCtx.q * ii + genCtx.q * jj * genCtx.ystride + genCtx.q * kk * genCtx.zstride; 317 318 std::pair<PetscInt, PetscInt> entoffset = DMMoab_SetElementConnectivity_Private(genCtx, ix, corner, conn); 319 320 for (PetscInt j = 0; j < entoffset.second; j++) { 321 /* The entity handle for the particular element -> if we want to set some tags is 322 moab::EntityHandle eh = starte + ie + j; 323 */ 324 gids[ie + j] = 1 + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)); 325 /* gids[ie+j] = ie + j + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)); */ 326 /* gids[ie+j] = 1 + ie; */ 327 /* ie++; */ 328 } 329 330 ix += entoffset.first; 331 ie += entoffset.second; 332 } 333 } 334 } 335 if (genCtx.adjEnts) { /* we need to update adjacencies now, because some elements are new */ 336 merr = iface->update_adjacencies(starte, nelems, nvperelem, conn);MBERR("Can't update adjacencies", merr); 337 } 338 tmp.swap(cells); 339 merr = mbImpl->tag_set_data(global_id_tag, cells, &gids[0]);MBERR("Can't set global ids to elements.", merr); 340 PetscFunctionReturn(0); 341 } 342 343 PetscErrorCode DMMBUtil_InitializeOptions(DMMoabMeshGeneratorCtx& genCtx, PetscInt dim, PetscBool simplex, PetscInt rank, PetscInt nprocs, const PetscReal* bounds, PetscInt nelems) 344 { 345 PetscFunctionBegin; 346 /* Initialize all genCtx data */ 347 genCtx.dim = dim; 348 genCtx.simplex = simplex; 349 genCtx.newMergeMethod = genCtx.keep_skins = genCtx.adjEnts = true; 350 /* determine other global quantities for the mesh used for nodes increments */ 351 genCtx.q = 1; 352 genCtx.fraction = genCtx.remainder = genCtx.cumfraction = 0; 353 354 if (!genCtx.usrxyzgrid) { /* not overridden by genCtx - assume nele equally and that genCtx wants a uniform cube mesh */ 355 356 genCtx.fraction = nelems / nprocs; /* partition only by the largest dimension */ 357 genCtx.remainder = nelems % nprocs; /* remainder after partition which gets evenly distributed by round-robin */ 358 genCtx.cumfraction = (rank > 0 ? (genCtx.fraction) * (rank) + (rank - 1 < genCtx.remainder ? rank : genCtx.remainder) : 0); 359 if (rank < genCtx.remainder) /* This process gets "fraction+1" elements */ 360 genCtx.fraction++; 361 362 PetscCall(PetscInfo(NULL, "Fraction = %" PetscInt_FMT ", Remainder = %" PetscInt_FMT ", Cumulative fraction = %" PetscInt_FMT "\n", genCtx.fraction, genCtx.remainder, genCtx.cumfraction)); 363 switch (genCtx.dim) { 364 case 1: 365 genCtx.blockSizeElementXYZ[0] = genCtx.fraction; 366 genCtx.blockSizeElementXYZ[1] = 1; 367 genCtx.blockSizeElementXYZ[2] = 1; 368 break; 369 case 2: 370 genCtx.blockSizeElementXYZ[0] = nelems; 371 genCtx.blockSizeElementXYZ[1] = genCtx.fraction; 372 genCtx.blockSizeElementXYZ[2] = 1; 373 break; 374 case 3: 375 default: 376 genCtx.blockSizeElementXYZ[0] = nelems; 377 genCtx.blockSizeElementXYZ[1] = nelems; 378 genCtx.blockSizeElementXYZ[2] = genCtx.fraction; 379 break; 380 } 381 } 382 383 /* partition only by the largest dimension */ 384 /* Total number of local elements := genCtx.blockSizeElementXYZ[0]*(genCtx.dim>1? genCtx.blockSizeElementXYZ[1]*(genCtx.dim>2 ? genCtx.blockSizeElementXYZ[2]:1) :1); */ 385 if (bounds) { 386 for (PetscInt i = 0; i < 6; i++) genCtx.xyzbounds[i] = bounds[i]; 387 } 388 else { 389 genCtx.xyzbounds[0] = genCtx.xyzbounds[2] = genCtx.xyzbounds[4] = 0.0; 390 genCtx.xyzbounds[1] = genCtx.xyzbounds[3] = genCtx.xyzbounds[5] = 1.0; 391 } 392 393 if (!genCtx.usrprocgrid) { 394 switch (genCtx.dim) { 395 case 1: 396 genCtx.M = nprocs; 397 genCtx.N = genCtx.K = 1; 398 break; 399 case 2: 400 genCtx.N = nprocs; 401 genCtx.M = genCtx.K = 1; 402 break; 403 default: 404 genCtx.K = nprocs; 405 genCtx.M = genCtx.N = 1; 406 break; 407 } 408 } 409 410 if (!genCtx.usrrefgrid) { 411 genCtx.A = genCtx.B = genCtx.C = 1; 412 } 413 414 /* more default values */ 415 genCtx.nex = genCtx.ney = genCtx.nez = 0; 416 genCtx.xstride = genCtx.ystride = genCtx.zstride = 0; 417 genCtx.NX = genCtx.NY = genCtx.NZ = 0; 418 genCtx.nex = genCtx.ney = genCtx.nez = 0; 419 genCtx.blockSizeVertexXYZ[0] = genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 1; 420 421 switch (genCtx.dim) { 422 case 3: 423 genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 424 genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1; 425 genCtx.blockSizeVertexXYZ[2] = genCtx.q * genCtx.blockSizeElementXYZ[2] + 1; 426 427 genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 428 genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */ 429 genCtx.NX = (genCtx.q * genCtx.nex + 1); 430 genCtx.xstride = 1; 431 genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1]; /* number of elements in y direction .... */ 432 genCtx.dy = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */ 433 genCtx.NY = (genCtx.q * genCtx.ney + 1); 434 genCtx.ystride = genCtx.blockSizeVertexXYZ[0]; 435 genCtx.nez = genCtx.K * genCtx.C * genCtx.blockSizeElementXYZ[2]; /* number of elements in z direction .... */ 436 genCtx.dz = (genCtx.xyzbounds[5] - genCtx.xyzbounds[4]) / (nelems * genCtx.q); /* distance between 2 nodes in z direction */ 437 genCtx.NZ = (genCtx.q * genCtx.nez + 1); 438 genCtx.zstride = genCtx.blockSizeVertexXYZ[0] * genCtx.blockSizeVertexXYZ[1]; 439 break; 440 case 2: 441 genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 442 genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1; 443 genCtx.blockSizeVertexXYZ[2] = 0; 444 445 genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 446 genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (genCtx.nex * genCtx.q); /* distance between 2 nodes in x direction */ 447 genCtx.NX = (genCtx.q * genCtx.nex + 1); 448 genCtx.xstride = 1; 449 genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1]; /* number of elements in y direction .... */ 450 genCtx.dy = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */ 451 genCtx.NY = (genCtx.q * genCtx.ney + 1); 452 genCtx.ystride = genCtx.blockSizeVertexXYZ[0]; 453 break; 454 case 1: 455 genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 0; 456 genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 457 458 genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 459 genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */ 460 genCtx.NX = (genCtx.q * genCtx.nex + 1); 461 genCtx.xstride = 1; 462 break; 463 } 464 465 /* Lets check for some valid input */ 466 PetscCheck(genCtx.dim >= 1 && genCtx.dim <= 3,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid topological dimension specified: %" PetscInt_FMT ".", genCtx.dim); 467 PetscCheck(genCtx.M * genCtx.N * genCtx.K == nprocs,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid [m, n, k] data: %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ". Product must be equal to global size = %" PetscInt_FMT ".", genCtx.M, genCtx.N, genCtx.K, nprocs); 468 /* validate the bounds data */ 469 PetscCheck(genCtx.xyzbounds[0] < genCtx.xyzbounds[1],PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "X-dim: Left boundary cannot be greater than right. [%G >= %G]", genCtx.xyzbounds[0], genCtx.xyzbounds[1]); 470 PetscCheck(genCtx.dim <= 1 || genCtx.xyzbounds[2] < genCtx.xyzbounds[3],PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Y-dim: Left boundary cannot be greater than right. [%G >= %G]", genCtx.xyzbounds[2], genCtx.xyzbounds[3]); 471 PetscCheck(genCtx.dim <= 2 || genCtx.xyzbounds[4] < genCtx.xyzbounds[5],PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Z-dim: Left boundary cannot be greater than right. [%G >= %G]", genCtx.xyzbounds[4], genCtx.xyzbounds[5]); 472 473 PetscCall(PetscInfo(NULL, "Local elements:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.blockSizeElementXYZ[0], genCtx.blockSizeElementXYZ[1], genCtx.blockSizeElementXYZ[2])); 474 PetscCall(PetscInfo(NULL, "Local vertices:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.blockSizeVertexXYZ[0], genCtx.blockSizeVertexXYZ[1], genCtx.blockSizeVertexXYZ[2])); 475 PetscCall(PetscInfo(NULL, "Local blocks/processors := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.A, genCtx.B, genCtx.C)); 476 PetscCall(PetscInfo(NULL, "Local processors := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.M, genCtx.N, genCtx.K)); 477 PetscCall(PetscInfo(NULL, "Local nexyz:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.nex, genCtx.ney, genCtx.nez)); 478 PetscCall(PetscInfo(NULL, "Local delxyz:= %g, %g, %g\n", genCtx.dx, genCtx.dy, genCtx.dz)); 479 PetscCall(PetscInfo(NULL, "Local strides:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.xstride, genCtx.ystride, genCtx.zstride)); 480 PetscFunctionReturn(0); 481 } 482 483 /*@C 484 DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with genCtx specified bounds. 485 486 Collective 487 488 Input Parameters: 489 + comm - The communicator for the DM object 490 . dim - The spatial dimension 491 . 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 492 . nele - The number of discrete elements in each direction 493 - user_nghost - The number of ghosted layers needed in the partitioned mesh 494 495 Output Parameter: 496 . dm - The DM object 497 498 Level: beginner 499 500 .seealso: `DMSetType()`, `DMCreate()`, `DMMoabLoadFromFile()` 501 @*/ 502 PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal* bounds, PetscInt nele, PetscInt nghost, DM *dm) 503 { 504 moab::ErrorCode merr; 505 PetscInt a, b, c, n, global_size, global_rank; 506 DM_Moab *dmmoab; 507 moab::Interface *mbImpl; 508 #ifdef MOAB_HAVE_MPI 509 moab::ParallelComm *pcomm; 510 #endif 511 moab::ReadUtilIface *readMeshIface; 512 moab::Range verts, cells, edges, faces, adj, dim3, dim2; 513 DMMoabMeshGeneratorCtx genCtx; 514 const PetscInt npts = nele + 1; /* Number of points in every dimension */ 515 516 moab::Tag global_id_tag, part_tag, geom_tag, mat_tag, dir_tag, neu_tag; 517 moab::Range ownedvtx, ownedelms, localvtxs, localelms; 518 moab::EntityHandle regionset; 519 PetscInt ml = 0, nl = 0, kl = 0; 520 521 PetscFunctionBegin; 522 PetscCheck(dim >= 1 && dim <= 3,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension argument for mesh: dim=[1,3]."); 523 524 PetscCall(PetscLogEventRegister("GenerateMesh", DM_CLASSID, &genCtx.generateMesh)); 525 PetscCall(PetscLogEventRegister("AddVertices", DM_CLASSID, &genCtx.generateVertices)); 526 PetscCall(PetscLogEventRegister("AddElements", DM_CLASSID, &genCtx.generateElements)); 527 PetscCall(PetscLogEventRegister("ParResolve", DM_CLASSID, &genCtx.parResolve)); 528 PetscCall(PetscLogEventBegin(genCtx.generateMesh, 0, 0, 0, 0)); 529 PetscCallMPI(MPI_Comm_size(comm, &global_size)); 530 /* total number of vertices in all dimensions */ 531 n = pow(npts, dim); 532 533 /* do some error checking */ 534 PetscCheck(n >= 2,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be >= 2."); 535 PetscCheck(global_size <= n,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of processors must be less than or equal to number of elements."); 536 PetscCheck(nghost >= 0,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of ghost layers cannot be negative."); 537 538 /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 539 PetscCall(DMMoabCreateMoab(comm, NULL, NULL, NULL, dm)); 540 541 /* get all the necessary handles from the private DM object */ 542 dmmoab = (DM_Moab*)(*dm)->data; 543 mbImpl = dmmoab->mbiface; 544 #ifdef MOAB_HAVE_MPI 545 pcomm = dmmoab->pcomm; 546 global_rank = pcomm->rank(); 547 #else 548 global_rank = 0; 549 global_size = 1; 550 #endif 551 global_id_tag = dmmoab->ltog_tag; 552 dmmoab->dim = dim; 553 dmmoab->nghostrings = nghost; 554 dmmoab->refct = 1; 555 556 /* create a file set to associate all entities in current mesh */ 557 merr = mbImpl->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 558 559 /* No errors yet; proceed with building the mesh */ 560 merr = mbImpl->query_interface(readMeshIface);MBERRNM(merr); 561 562 genCtx.M = genCtx.N = genCtx.K = 1; 563 genCtx.A = genCtx.B = genCtx.C = 1; 564 genCtx.blockSizeElementXYZ[0] = 0; 565 genCtx.blockSizeElementXYZ[1] = 0; 566 genCtx.blockSizeElementXYZ[2] = 0; 567 568 PetscOptionsBegin(comm, "", "DMMoab Creation Options", "DMMOAB"); 569 /* Handle DMMoab spatial resolution */ 570 PetscCall(PetscOptionsInt("-dmb_grid_x", "Number of grid points in x direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[0], &genCtx.blockSizeElementXYZ[0], &genCtx.usrxyzgrid)); 571 if (dim > 1) PetscCall(PetscOptionsInt("-dmb_grid_y", "Number of grid points in y direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[1], &genCtx.blockSizeElementXYZ[1], &genCtx.usrxyzgrid)); 572 if (dim > 2) PetscCall(PetscOptionsInt("-dmb_grid_z", "Number of grid points in z direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[2], &genCtx.blockSizeElementXYZ[2], &genCtx.usrxyzgrid)); 573 574 /* Handle DMMoab parallel distribution */ 575 PetscCall(PetscOptionsInt("-dmb_processors_x", "Number of processors in x direction", "DMMoabSetNumProcs", genCtx.M, &genCtx.M, &genCtx.usrprocgrid)); 576 if (dim > 1) PetscCall(PetscOptionsInt("-dmb_processors_y", "Number of processors in y direction", "DMMoabSetNumProcs", genCtx.N, &genCtx.N, &genCtx.usrprocgrid)); 577 if (dim > 2) PetscCall(PetscOptionsInt("-dmb_processors_z", "Number of processors in z direction", "DMMoabSetNumProcs", genCtx.K, &genCtx.K, &genCtx.usrprocgrid)); 578 579 /* Handle DMMoab block refinement */ 580 PetscCall(PetscOptionsInt("-dmb_refine_x", "Number of refinement blocks in x direction", "DMMoabSetRefinement", genCtx.A, &genCtx.A, &genCtx.usrrefgrid)); 581 if (dim > 1) PetscCall(PetscOptionsInt("-dmb_refine_y", "Number of refinement blocks in y direction", "DMMoabSetRefinement", genCtx.B, &genCtx.B, &genCtx.usrrefgrid)); 582 if (dim > 2) PetscCall(PetscOptionsInt("-dmb_refine_z", "Number of refinement blocks in z direction", "DMMoabSetRefinement", genCtx.C, &genCtx.C, &genCtx.usrrefgrid)); 583 PetscOptionsEnd(); 584 585 PetscCall(DMMBUtil_InitializeOptions(genCtx, dim, useSimplex, global_rank, global_size, bounds, nele)); 586 587 //PetscCheck(nele>=nprocs,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The dimensional discretization size should be greater or equal to number of processors: %" PetscInt_FMT " < %" PetscInt_FMT,nele,nprocs); 588 589 if (genCtx.adjEnts) genCtx.keep_skins = true; /* do not delete anything - consumes more memory */ 590 591 /* determine m, n, k for processor rank */ 592 ml = nl = kl = 0; 593 switch (genCtx.dim) { 594 case 1: 595 ml = (genCtx.cumfraction); 596 break; 597 case 2: 598 nl = (genCtx.cumfraction); 599 break; 600 default: 601 kl = (genCtx.cumfraction) / genCtx.q / genCtx.blockSizeElementXYZ[2] / genCtx.C; //genCtx.K 602 break; 603 } 604 605 /* 606 * so there are a total of M * A * blockSizeElement elements in x direction (so M * A * blockSizeElement + 1 verts in x direction) 607 * so there are a total of N * B * blockSizeElement elements in y direction (so N * B * blockSizeElement + 1 verts in y direction) 608 * so there are a total of K * C * blockSizeElement elements in z direction (so K * C * blockSizeElement + 1 verts in z direction) 609 610 * there are ( M * A blockSizeElement) * ( N * B * blockSizeElement) * (K * C * blockSizeElement) hexas 611 * there are ( M * A * blockSizeElement + 1) * ( N * B * blockSizeElement + 1) * (K * C * blockSizeElement + 1) vertices 612 * x is the first dimension that varies 613 */ 614 615 /* generate the block at (a, b, c); it will represent a partition , it will get a partition tag */ 616 PetscInt dum_id = -1; 617 merr = mbImpl->tag_get_handle("GLOBAL_ID", 1, moab::MB_TYPE_INTEGER, global_id_tag);MBERR("Getting Global_ID Tag handle failed", merr); 618 619 merr = mbImpl->tag_get_handle(MATERIAL_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, mat_tag);MBERR("Getting Material set Tag handle failed", merr); 620 merr = mbImpl->tag_get_handle(DIRICHLET_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, dir_tag);MBERR("Getting Dirichlet set Tag handle failed", merr); 621 merr = mbImpl->tag_get_handle(NEUMANN_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, neu_tag);MBERR("Getting Neumann set Tag handle failed", merr); 622 623 merr = mbImpl->tag_get_handle("PARALLEL_PARTITION", 1, moab::MB_TYPE_INTEGER, part_tag, moab::MB_TAG_CREAT | moab::MB_TAG_SPARSE, &dum_id);MBERR("Getting Partition Tag handle failed", merr); 624 625 /* lets create some sets */ 626 merr = mbImpl->tag_get_handle(GEOM_DIMENSION_TAG_NAME, 1, moab::MB_TYPE_INTEGER, geom_tag, moab::MB_TAG_CREAT | moab::MB_TAG_SPARSE, &dum_id);MBERRNM(merr); 627 merr = mbImpl->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr); 628 PetscCall(PetscLogEventEnd(genCtx.generateMesh, 0, 0, 0, 0)); 629 630 for (a = 0; a < (genCtx.dim > 0 ? genCtx.A : genCtx.A); a++) { 631 for (b = 0; b < (genCtx.dim > 1 ? genCtx.B : 1); b++) { 632 for (c = 0; c < (genCtx.dim > 2 ? genCtx.C : 1); c++) { 633 634 moab::EntityHandle startv; 635 636 PetscCall(PetscLogEventBegin(genCtx.generateVertices, 0, 0, 0, 0)); 637 PetscCall(DMMoab_GenerateVertices_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, verts)); 638 PetscCall(PetscLogEventEnd(genCtx.generateVertices, 0, 0, 0, 0)); 639 640 PetscCall(PetscLogEventBegin(genCtx.generateElements, 0, 0, 0, 0)); 641 PetscCall(DMMoab_GenerateElements_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, cells)); 642 PetscCall(PetscLogEventEnd(genCtx.generateElements, 0, 0, 0, 0)); 643 644 PetscInt part_num = 0; 645 switch (genCtx.dim) { 646 case 3: 647 part_num += (c + kl * genCtx.C) * (genCtx.M * genCtx.A * genCtx.N * genCtx.B); 648 case 2: 649 part_num += (b + nl * genCtx.B) * (genCtx.M * genCtx.A); 650 case 1: 651 part_num += (a + ml * genCtx.A); 652 break; 653 } 654 655 moab::EntityHandle part_set; 656 merr = mbImpl->create_meshset(moab::MESHSET_SET, part_set);MBERR("Can't create mesh set.", merr); 657 658 merr = mbImpl->add_entities(part_set, verts);MBERR("Can't add vertices to set.", merr); 659 merr = mbImpl->add_entities(part_set, cells);MBERR("Can't add entities to set.", merr); 660 merr = mbImpl->add_entities(regionset, cells);MBERR("Can't add entities to set.", merr); 661 662 /* if needed, add all edges and faces */ 663 if (genCtx.adjEnts) 664 { 665 if (genCtx.dim > 1) { 666 merr = mbImpl->get_adjacencies(cells, 1, true, edges, moab::Interface::UNION);MBERR("Can't get edges", merr); 667 merr = mbImpl->add_entities(part_set, edges);MBERR("Can't add edges to partition set.", merr); 668 } 669 if (genCtx.dim > 2) { 670 merr = mbImpl->get_adjacencies(cells, 2, true, faces, moab::Interface::UNION);MBERR("Can't get faces", merr); 671 merr = mbImpl->add_entities(part_set, faces);MBERR("Can't add faces to partition set.", merr); 672 } 673 edges.clear(); 674 faces.clear(); 675 } 676 verts.clear(); cells.clear(); 677 678 merr = mbImpl->tag_set_data(part_tag, &part_set, 1, &part_num);MBERR("Can't set part tag on set", merr); 679 if (dmmoab->fileset) { 680 merr = mbImpl->add_parent_child(dmmoab->fileset, part_set);MBERR("Can't add part set to file set.", merr); 681 merr = mbImpl->unite_meshset(dmmoab->fileset, part_set);MBERRNM(merr); 682 } 683 merr = mbImpl->add_entities(dmmoab->fileset, &part_set, 1);MBERRNM(merr); 684 } 685 } 686 } 687 688 merr = mbImpl->add_parent_child(dmmoab->fileset, regionset);MBERRNM(merr); 689 690 /* Only in parallel: resolve shared entities between processors and exchange ghost layers */ 691 if (global_size > 1) { 692 693 PetscCall(PetscLogEventBegin(genCtx.parResolve, 0, 0, 0, 0)); 694 695 merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, genCtx.dim, cells);MBERR("Can't get all d-dimensional elements.", merr); 696 merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 0, verts);MBERR("Can't get all vertices.", merr); 697 698 if (genCtx.A * genCtx.B * genCtx.C != 1) { // merge needed 699 moab::MergeMesh mm(mbImpl); 700 if (genCtx.newMergeMethod) { 701 merr = mm.merge_using_integer_tag(verts, global_id_tag);MBERR("Can't merge with GLOBAL_ID tag", merr); 702 } 703 else { 704 merr = mm.merge_entities(cells, 0.0001);MBERR("Can't merge with coordinates", merr); 705 } 706 } 707 708 #ifdef MOAB_HAVE_MPI 709 /* check the handles */ 710 merr = pcomm->check_all_shared_handles();MBERRV(mbImpl, merr); 711 712 /* resolve the shared entities by exchanging information to adjacent processors */ 713 merr = pcomm->resolve_shared_ents(dmmoab->fileset, cells, dim, dim - 1, NULL, &global_id_tag);MBERRV(mbImpl, merr); 714 if (dmmoab->fileset) { 715 merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false, &dmmoab->fileset);MBERRV(mbImpl, merr); 716 } 717 else { 718 merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false);MBERRV(mbImpl, merr); 719 } 720 721 /* Reassign global IDs on all entities. */ 722 merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, false, true, false);MBERRNM(merr); 723 #endif 724 725 PetscCall(PetscLogEventEnd(genCtx.parResolve, 0, 0, 0, 0)); 726 } 727 728 if (!genCtx.keep_skins) { // default is to delete the 1- and 2-dimensional entities 729 // delete all quads and edges 730 moab::Range toDelete; 731 if (genCtx.dim > 1) { 732 merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 1, toDelete);MBERR("Can't get edges", merr); 733 } 734 735 if (genCtx.dim > 2) { 736 merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 2, toDelete);MBERR("Can't get faces", merr); 737 } 738 739 #ifdef MOAB_HAVE_MPI 740 merr = dmmoab->pcomm->delete_entities(toDelete) ;MBERR("Can't delete entities", merr); 741 #endif 742 } 743 744 /* set geometric dimension tag for regions */ 745 merr = mbImpl->tag_set_data(geom_tag, ®ionset, 1, &dmmoab->dim);MBERRNM(merr); 746 /* set default material ID for regions */ 747 int default_material = 1; 748 merr = mbImpl->tag_set_data(mat_tag, ®ionset, 1, &default_material);MBERRNM(merr); 749 /* 750 int default_dbc = 0; 751 merr = mbImpl->tag_set_data(dir_tag, &vertexset, 1, &default_dbc);MBERRNM(merr); 752 */ 753 PetscFunctionReturn(0); 754 } 755 756 PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, PetscInt nghost, MoabReadMode mode, PetscInt dbglevel, const char* dm_opts, const char* extra_opts, const char** read_opts) 757 { 758 char *ropts; 759 char ropts_par[PETSC_MAX_PATH_LEN], ropts_pargh[PETSC_MAX_PATH_LEN]; 760 char ropts_dbg[PETSC_MAX_PATH_LEN]; 761 762 PetscFunctionBegin; 763 PetscCall(PetscMalloc1(PETSC_MAX_PATH_LEN, &ropts)); 764 PetscCall(PetscMemzero(&ropts_par, PETSC_MAX_PATH_LEN)); 765 PetscCall(PetscMemzero(&ropts_pargh, PETSC_MAX_PATH_LEN)); 766 PetscCall(PetscMemzero(&ropts_dbg, PETSC_MAX_PATH_LEN)); 767 768 /* do parallel read unless using only one processor */ 769 if (numproc > 1) { 770 // PetscCall(PetscSNPrintf(ropts_par, PETSC_MAX_PATH_LEN, "PARALLEL=%s;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=%d.0.1%s;",MoabReadModes[mode],dim,(by_rank ? ";PARTITION_BY_RANK":""))); 771 PetscCall(PetscSNPrintf(ropts_par, PETSC_MAX_PATH_LEN, "PARALLEL=%s;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS;%s", MoabReadModes[mode], (by_rank ? "PARTITION_BY_RANK;" : ""))); 772 if (nghost) { 773 PetscCall(PetscSNPrintf(ropts_pargh, PETSC_MAX_PATH_LEN, "PARALLEL_GHOSTS=%" PetscInt_FMT ".0.%" PetscInt_FMT ";", dim, nghost)); 774 } 775 } 776 777 if (dbglevel) { 778 if (numproc > 1) { 779 PetscCall(PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%" PetscInt_FMT ";DEBUG_PIO=%" PetscInt_FMT ";", dbglevel, dbglevel)); 780 } 781 else PetscCall(PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%" PetscInt_FMT ";", dbglevel)); 782 } 783 784 PetscCall(PetscSNPrintf(ropts, PETSC_MAX_PATH_LEN, "%s%s%s%s%s", ropts_par, (nghost ? ropts_pargh : ""), ropts_dbg, (extra_opts ? extra_opts : ""), (dm_opts ? dm_opts : ""))); 785 *read_opts = ropts; 786 PetscFunctionReturn(0); 787 } 788 789 /*@C 790 DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file. 791 792 Collective 793 794 Input Parameters: 795 + comm - The communicator for the DM object 796 . dim - The spatial dimension 797 . filename - The name of the mesh file to be loaded 798 - usrreadopts - The options string to read a MOAB mesh. 799 800 Reference (Parallel Mesh Initialization: https://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo) 801 802 Output Parameter: 803 . dm - The DM object 804 805 Level: beginner 806 807 .seealso: `DMSetType()`, `DMCreate()`, `DMMoabCreateBoxMesh()` 808 @*/ 809 PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm, PetscInt dim, PetscInt nghost, const char* filename, const char* usrreadopts, DM *dm) 810 { 811 moab::ErrorCode merr; 812 PetscInt nprocs; 813 DM_Moab *dmmoab; 814 moab::Interface *mbiface; 815 #ifdef MOAB_HAVE_MPI 816 moab::ParallelComm *pcomm; 817 #endif 818 moab::Range verts, elems; 819 const char *readopts; 820 821 PetscFunctionBegin; 822 PetscValidPointer(dm, 6); 823 824 /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 825 PetscCall(DMMoabCreateMoab(comm, NULL, NULL, NULL, dm)); 826 827 /* get all the necessary handles from the private DM object */ 828 dmmoab = (DM_Moab*)(*dm)->data; 829 mbiface = dmmoab->mbiface; 830 #ifdef MOAB_HAVE_MPI 831 pcomm = dmmoab->pcomm; 832 nprocs = pcomm->size(); 833 #else 834 nprocs = 1; 835 #endif 836 /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */ 837 dmmoab->dim = dim; 838 dmmoab->nghostrings = nghost; 839 dmmoab->refct = 1; 840 841 /* create a file set to associate all entities in current mesh */ 842 merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 843 844 /* add mesh loading options specific to the DM */ 845 PetscCall(DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, nghost, dmmoab->read_mode, 846 dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts)); 847 848 PetscCall(PetscInfo(*dm, "Reading file %s with options: %s\n", filename, readopts)); 849 850 /* Load the mesh from a file. */ 851 if (dmmoab->fileset) { 852 merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface, "Reading MOAB file failed.", merr); 853 } 854 else { 855 merr = mbiface->load_file(filename, 0, readopts);MBERRVM(mbiface, "Reading MOAB file failed.", merr); 856 } 857 858 #ifdef MOAB_HAVE_MPI 859 /* Reassign global IDs on all entities. */ 860 /* merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, true, true, true);MBERRNM(merr); */ 861 #endif 862 863 /* load the local vertices */ 864 merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 865 /* load the local elements */ 866 merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 867 868 #ifdef MOAB_HAVE_MPI 869 /* Everything is set up, now just do a tag exchange to update tags 870 on all of the ghost vertexes */ 871 merr = pcomm->exchange_tags(dmmoab->ltog_tag, verts);MBERRV(mbiface, merr); 872 merr = pcomm->exchange_tags(dmmoab->ltog_tag, elems);MBERRV(mbiface, merr); 873 merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 874 #endif 875 876 PetscCall(PetscInfo(*dm, "MOAB file '%s' was successfully loaded. Found %zu vertices and %zu elements.\n", filename, verts.size(), elems.size())); 877 PetscCall(PetscFree(readopts)); 878 PetscFunctionReturn(0); 879 } 880 881 /*@C 882 DMMoabRenumberMeshEntities - Order and number all entities (vertices->elements) to be contiguously ordered 883 in parallel 884 885 Collective 886 887 Input Parameters: 888 . dm - The DM object 889 890 Level: advanced 891 892 .seealso: `DMSetUp()`, `DMCreate()` 893 @*/ 894 PetscErrorCode DMMoabRenumberMeshEntities(DM dm) 895 { 896 moab::Range verts; 897 898 PetscFunctionBegin; 899 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 900 901 #ifdef MOAB_HAVE_MPI 902 /* Insert new points */ 903 moab::ErrorCode merr; 904 merr = ((DM_Moab*) dm->data)->pcomm->assign_global_ids(((DM_Moab*) dm->data)->fileset, 3, 0, false, true, false);MBERRNM(merr); 905 #endif 906 PetscFunctionReturn(0); 907 } 908