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