1af0996ceSBarry Smith #include <petsc/private/dmmbimpl.h> /*I "petscdmmoab.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/vecimpl.h> 351d15aeeSVijay Mahadevan 451d15aeeSVijay Mahadevan #include <petscdmmoab.h> 551d15aeeSVijay Mahadevan #include <MBTagConventions.hpp> 651d15aeeSVijay Mahadevan #include <moab/ReadUtilIface.hpp> 7ce27a4eeSVijay Mahadevan #include <moab/MergeMesh.hpp> 851d15aeeSVijay Mahadevan #include <moab/CN.hpp> 951d15aeeSVijay Mahadevan 10ce27a4eeSVijay Mahadevan typedef struct { 11ce27a4eeSVijay Mahadevan // options 1249d66b22SVijay Mahadevan PetscInt A, B, C, M, N, K, dim; 13ce27a4eeSVijay Mahadevan PetscInt blockSizeVertexXYZ[3]; // Number of element blocks per partition 14ce27a4eeSVijay Mahadevan PetscInt blockSizeElementXYZ[3]; 15ce27a4eeSVijay Mahadevan PetscReal xyzbounds[6]; // the physical size of the domain 16ce27a4eeSVijay Mahadevan bool newMergeMethod, keep_skins, simplex, adjEnts; 1751d15aeeSVijay Mahadevan 18ce27a4eeSVijay Mahadevan // compute params 19ce27a4eeSVijay Mahadevan PetscReal dx, dy, dz; 20ce27a4eeSVijay Mahadevan PetscInt NX, NY, NZ, nex, ney, nez; 21ce27a4eeSVijay Mahadevan PetscInt q, xstride, ystride, zstride; 22ce27a4eeSVijay Mahadevan PetscBool usrxyzgrid, usrprocgrid, usrrefgrid; 23755f3dfbSVijay Mahadevan PetscInt fraction, remainder, cumfraction; 24755f3dfbSVijay Mahadevan PetscLogEvent generateMesh, generateElements, generateVertices, parResolve; 25ce27a4eeSVijay Mahadevan } DMMoabMeshGeneratorCtx; 26ce27a4eeSVijay Mahadevan 2749d66b22SVijay Mahadevan PetscInt DMMoab_SetTensorElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt offset, PetscInt corner, std::vector<PetscInt>& subent_conn, moab::EntityHandle *connectivity) 2851d15aeeSVijay Mahadevan { 29ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 30e5136372SVijay Mahadevan case 1: 31ce27a4eeSVijay Mahadevan subent_conn.resize(2); 32ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data()); 33ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[0]] = corner; 34ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[1]] = corner + 1; 35e5136372SVijay Mahadevan break; 36e5136372SVijay Mahadevan case 2: 37ce27a4eeSVijay Mahadevan subent_conn.resize(4); 38ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBQUAD, 2, 0, subent_conn.data()); 39ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[0]] = corner; 40ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[1]] = corner + 1; 41ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[2]] = corner + 1 + genCtx.ystride; 42ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[3]] = corner + genCtx.ystride; 43e5136372SVijay Mahadevan break; 44e5136372SVijay Mahadevan case 3: 45e5136372SVijay Mahadevan default: 46ce27a4eeSVijay Mahadevan subent_conn.resize(8); 47ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBHEX, 3, 0, subent_conn.data()); 48ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[0]] = corner; 49ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[1]] = corner + 1; 50ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[2]] = corner + 1 + genCtx.ystride; 51ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[3]] = corner + genCtx.ystride; 52ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[4]] = corner + genCtx.zstride; 53ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[5]] = corner + 1 + genCtx.zstride; 54ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[6]] = corner + 1 + genCtx.ystride + genCtx.zstride; 55ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[7]] = corner + genCtx.ystride + genCtx.zstride; 56e5136372SVijay Mahadevan break; 57e5136372SVijay Mahadevan } 58ce27a4eeSVijay Mahadevan return subent_conn.size(); 59e5136372SVijay Mahadevan } 6051d15aeeSVijay Mahadevan 6149d66b22SVijay Mahadevan PetscInt DMMoab_SetSimplexElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt subelem, PetscInt offset, PetscInt corner, std::vector<PetscInt>& subent_conn, moab::EntityHandle *connectivity) 62c3dd00c7SVijay Mahadevan { 63755f3dfbSVijay Mahadevan PetscInt A, B, C, D, E, F, G, H, M; 64755f3dfbSVijay Mahadevan const PetscInt trigen_opts = 1; /* 1 - Aligned diagonally to right, 2 - Aligned diagonally to left, 3 - 4 elements per quad */ 6549d66b22SVijay Mahadevan A = corner; 6649d66b22SVijay Mahadevan B = corner + 1; 67ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 68c3dd00c7SVijay Mahadevan case 1: 69ce27a4eeSVijay Mahadevan subent_conn.resize(2); /* only linear EDGE supported now */ 70ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data()); 7149d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 7249d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = B; 73c3dd00c7SVijay Mahadevan break; 74c3dd00c7SVijay Mahadevan case 2: 7549d66b22SVijay Mahadevan C = corner + 1 + genCtx.ystride; 7649d66b22SVijay Mahadevan D = corner + genCtx.ystride; 7763d025dbSVijay Mahadevan M = corner + 0.5; /* technically -- need to modify vertex generation */ 78ce27a4eeSVijay Mahadevan subent_conn.resize(3); /* only linear TRI supported */ 79ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBTRI, 2, 0, subent_conn.data()); 80755f3dfbSVijay Mahadevan if (trigen_opts == 1) { 81ce27a4eeSVijay Mahadevan if (subelem) { /* 0 1 2 of a QUAD */ 8277a8c067SVijay Mahadevan connectivity[offset + subent_conn[0]] = B; 83755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[1]] = C; 8477a8c067SVijay Mahadevan connectivity[offset + subent_conn[2]] = A; 85*1baa6e33SBarry Smith } else { /* 2 3 0 of a QUAD */ 8677a8c067SVijay Mahadevan connectivity[offset + subent_conn[0]] = D; 8777a8c067SVijay Mahadevan connectivity[offset + subent_conn[1]] = A; 88755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = C; 89755f3dfbSVijay Mahadevan } 90*1baa6e33SBarry Smith } else if (trigen_opts == 2) { 91755f3dfbSVijay Mahadevan if (subelem) { /* 0 1 2 of a QUAD */ 92755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 9377a8c067SVijay Mahadevan connectivity[offset + subent_conn[1]] = B; 9477a8c067SVijay Mahadevan connectivity[offset + subent_conn[2]] = D; 95*1baa6e33SBarry Smith } else { /* 2 3 0 of a QUAD */ 9677a8c067SVijay Mahadevan connectivity[offset + subent_conn[0]] = C; 9777a8c067SVijay Mahadevan connectivity[offset + subent_conn[1]] = D; 98755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = B; 99755f3dfbSVijay Mahadevan } 100*1baa6e33SBarry Smith } else { 101755f3dfbSVijay Mahadevan switch (subelem) { /* 0 1 2 of a QUAD */ 102755f3dfbSVijay Mahadevan case 0: 103755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 104755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[1]] = B; 105755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = M; 106755f3dfbSVijay Mahadevan break; 107755f3dfbSVijay Mahadevan case 1: 108755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[0]] = B; 109755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[1]] = C; 110755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = M; 111755f3dfbSVijay Mahadevan break; 112755f3dfbSVijay Mahadevan case 2: 11349d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = C; 11449d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = D; 115755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = M; 116755f3dfbSVijay Mahadevan break; 117755f3dfbSVijay Mahadevan case 3: 118755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[0]] = D; 119755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[1]] = A; 120755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = M; 121755f3dfbSVijay Mahadevan break; 122755f3dfbSVijay Mahadevan } 123c3dd00c7SVijay Mahadevan } 124c3dd00c7SVijay Mahadevan break; 125c3dd00c7SVijay Mahadevan case 3: 126c3dd00c7SVijay Mahadevan default: 12749d66b22SVijay Mahadevan C = corner + 1 + genCtx.ystride; 12849d66b22SVijay Mahadevan D = corner + genCtx.ystride; 12949d66b22SVijay Mahadevan E = corner + genCtx.zstride; 13049d66b22SVijay Mahadevan F = corner + 1 + genCtx.zstride; 13149d66b22SVijay Mahadevan G = corner + 1 + genCtx.ystride + genCtx.zstride; 13249d66b22SVijay Mahadevan H = corner + genCtx.ystride + genCtx.zstride; 133ce27a4eeSVijay Mahadevan subent_conn.resize(4); /* only linear TET supported */ 134ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBTET, 3, 0, subent_conn.data()); 135c3dd00c7SVijay Mahadevan switch (subelem) { 13649d66b22SVijay Mahadevan case 0: /* 4 3 7 6 of a HEX */ 13749d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = E; 13849d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = D; 13949d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = H; 14049d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = G; 141c3dd00c7SVijay Mahadevan break; 14249d66b22SVijay Mahadevan case 1: /* 0 1 2 5 of a HEX */ 14349d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 14449d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = B; 14549d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = C; 14649d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 147c3dd00c7SVijay Mahadevan break; 14849d66b22SVijay Mahadevan case 2: /* 0 3 4 5 of a HEX */ 14949d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 15049d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = D; 15149d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = E; 15249d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 153c3dd00c7SVijay Mahadevan break; 15449d66b22SVijay Mahadevan case 3: /* 2 6 3 5 of a HEX */ 15549d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = C; 15649d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = G; 15749d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = D; 15849d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 159c3dd00c7SVijay Mahadevan break; 16049d66b22SVijay Mahadevan case 4: /* 0 2 3 5 of a HEX */ 16149d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 16249d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = C; 16349d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = D; 16449d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 16549d66b22SVijay Mahadevan break; 16649d66b22SVijay Mahadevan case 5: /* 3 6 4 5 of a HEX */ 16749d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = D; 16849d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = G; 16949d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = E; 17049d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 171c3dd00c7SVijay Mahadevan break; 172c3dd00c7SVijay Mahadevan } 173c3dd00c7SVijay Mahadevan break; 174c3dd00c7SVijay Mahadevan } 175ce27a4eeSVijay Mahadevan return subent_conn.size(); 176c3dd00c7SVijay Mahadevan } 177c3dd00c7SVijay Mahadevan 17849d66b22SVijay Mahadevan std::pair<PetscInt, PetscInt> DMMoab_SetElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt offset, PetscInt corner, moab::EntityHandle *connectivity) 179c3dd00c7SVijay Mahadevan { 18049d66b22SVijay Mahadevan PetscInt vcount = 0; 18149d66b22SVijay Mahadevan PetscInt simplices_per_tensor[4] = {0, 1, 2, 6}; 18249d66b22SVijay Mahadevan std::vector<PetscInt> subent_conn; /* only linear edge, tri, tet supported now */ 183ce27a4eeSVijay Mahadevan subent_conn.reserve(27); 18449d66b22SVijay Mahadevan PetscInt m, subelem; 185ce27a4eeSVijay Mahadevan if (genCtx.simplex) { 18649d66b22SVijay Mahadevan subelem = simplices_per_tensor[genCtx.dim]; 187ce27a4eeSVijay Mahadevan for (m = 0; m < subelem; m++) { 188ce27a4eeSVijay Mahadevan vcount = DMMoab_SetSimplexElementConnectivity_Private(genCtx, m, offset, corner, subent_conn, connectivity); 189ce27a4eeSVijay Mahadevan offset += vcount; 190ce27a4eeSVijay Mahadevan } 191*1baa6e33SBarry Smith } else { 192c3dd00c7SVijay Mahadevan subelem = 1; 193ce27a4eeSVijay Mahadevan vcount = DMMoab_SetTensorElementConnectivity_Private(genCtx, offset, corner, subent_conn, connectivity); 194c3dd00c7SVijay Mahadevan } 19549d66b22SVijay Mahadevan return std::pair<PetscInt, PetscInt>(vcount * subelem, subelem); 196c3dd00c7SVijay Mahadevan } 197c3dd00c7SVijay Mahadevan 19849d66b22SVijay Mahadevan PetscErrorCode DMMoab_GenerateVertices_Private(moab::Interface *mbImpl, moab::ReadUtilIface *iface, DMMoabMeshGeneratorCtx& genCtx, PetscInt m, PetscInt n, PetscInt k, 19949d66b22SVijay Mahadevan PetscInt a, PetscInt b, PetscInt c, moab::Tag& global_id_tag, moab::EntityHandle& startv, moab::Range& uverts) 200ce27a4eeSVijay Mahadevan { 20149d66b22SVijay Mahadevan PetscInt x, y, z, ix, nnodes; 20249d66b22SVijay Mahadevan PetscInt ii, jj, kk; 20349d66b22SVijay Mahadevan std::vector<PetscReal*> arrays; 204755f3dfbSVijay Mahadevan PetscInt* gids; 205ce27a4eeSVijay Mahadevan moab::ErrorCode merr; 206ce27a4eeSVijay Mahadevan 207ce27a4eeSVijay Mahadevan PetscFunctionBegin; 208ce27a4eeSVijay Mahadevan /* we will generate (q*block+1)^3 vertices, and block^3 hexas; q is 1 for linear, 2 for quadratic 209ce27a4eeSVijay Mahadevan * the global id of the vertices will come from m, n, k, a, b, c 210ce27a4eeSVijay Mahadevan * x will vary from m*A*q*block + a*q*block to m*A*q*block+(a+1)*q*block etc. 211ce27a4eeSVijay Mahadevan */ 212ce27a4eeSVijay Mahadevan nnodes = genCtx.blockSizeVertexXYZ[0] * (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] * (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1) : 1); 2139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnodes, &gids)); 214ce27a4eeSVijay Mahadevan 215ce27a4eeSVijay Mahadevan merr = iface->get_node_coords(3, nnodes, 0, startv, arrays);MBERR("Can't get node coords.", merr); 216ce27a4eeSVijay Mahadevan 217ce27a4eeSVijay Mahadevan /* will start with the lower corner: */ 21863d025dbSVijay Mahadevan /* x = ( m * genCtx.A + a) * genCtx.q * genCtx.blockSizeElementXYZ[0]; */ 21963d025dbSVijay Mahadevan /* y = ( n * genCtx.B + b) * genCtx.q * genCtx.blockSizeElementXYZ[1]; */ 22063d025dbSVijay Mahadevan /* z = ( k * genCtx.C + c) * genCtx.q * genCtx.blockSizeElementXYZ[2]; */ 221755f3dfbSVijay Mahadevan 222755f3dfbSVijay Mahadevan x = ( m * genCtx.A + a) * genCtx.q; 223755f3dfbSVijay Mahadevan y = ( n * genCtx.B + b) * genCtx.q; 224755f3dfbSVijay Mahadevan z = ( k * genCtx.C + c) * genCtx.q; 2259566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Starting offset for coordinates := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", x, y, z)); 226ce27a4eeSVijay Mahadevan ix = 0; 227ce27a4eeSVijay Mahadevan moab::Range verts(startv, startv + nnodes - 1); 228ce27a4eeSVijay Mahadevan for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1); kk++) { 229ce27a4eeSVijay Mahadevan for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] : 1); jj++) { 230ce27a4eeSVijay Mahadevan for (ii = 0; ii < genCtx.blockSizeVertexXYZ[0]; ii++, ix++) { 231ce27a4eeSVijay Mahadevan /* set coordinates for the vertices */ 232ce27a4eeSVijay Mahadevan arrays[0][ix] = (x + ii) * genCtx.dx + genCtx.xyzbounds[0]; 233ce27a4eeSVijay Mahadevan arrays[1][ix] = (y + jj) * genCtx.dy + genCtx.xyzbounds[2]; 234ce27a4eeSVijay Mahadevan arrays[2][ix] = (z + kk) * genCtx.dz + genCtx.xyzbounds[4]; 2359566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Creating vertex with coordinates := %f, %f, %f\n", arrays[0][ix], arrays[1][ix], arrays[2][ix])); 236ce27a4eeSVijay Mahadevan 237ce27a4eeSVijay Mahadevan /* If we want to set some tags on the vertices -> use the following entity handle definition: 238ce27a4eeSVijay Mahadevan moab::EntityHandle v = startv + ix; 239ce27a4eeSVijay Mahadevan */ 240ce27a4eeSVijay Mahadevan /* compute the global ID for vertex */ 241ce27a4eeSVijay Mahadevan gids[ix] = 1 + (x + ii) + (y + jj) * genCtx.NX + (z + kk) * (genCtx.NX * genCtx.NY); 242ce27a4eeSVijay Mahadevan } 243ce27a4eeSVijay Mahadevan } 244ce27a4eeSVijay Mahadevan } 245ce27a4eeSVijay Mahadevan /* set global ID data on vertices */ 246ce27a4eeSVijay Mahadevan mbImpl->tag_set_data(global_id_tag, verts, &gids[0]); 247ce27a4eeSVijay Mahadevan verts.swap(uverts); 2489566063dSJacob Faibussowitsch PetscCall(PetscFree(gids)); 249ce27a4eeSVijay Mahadevan PetscFunctionReturn(0); 250ce27a4eeSVijay Mahadevan } 251ce27a4eeSVijay Mahadevan 25249d66b22SVijay Mahadevan PetscErrorCode DMMoab_GenerateElements_Private(moab::Interface* mbImpl, moab::ReadUtilIface* iface, DMMoabMeshGeneratorCtx& genCtx, PetscInt m, PetscInt n, PetscInt k, 25349d66b22SVijay Mahadevan PetscInt a, PetscInt b, PetscInt c, moab::Tag& global_id_tag, moab::EntityHandle startv, moab::Range& cells) 254ce27a4eeSVijay Mahadevan { 255ce27a4eeSVijay Mahadevan moab::ErrorCode merr; 256ce27a4eeSVijay Mahadevan PetscInt ix, ie, xe, ye, ze; 257ce27a4eeSVijay Mahadevan PetscInt ii, jj, kk, nvperelem; 25849d66b22SVijay Mahadevan PetscInt simplices_per_tensor[4] = {0, 1, 2, 6}; 259ce27a4eeSVijay Mahadevan PetscInt ntensorelems = genCtx.blockSizeElementXYZ[0] * (genCtx.dim > 1 ? genCtx.blockSizeElementXYZ[1] * (genCtx.dim > 2 ? genCtx.blockSizeElementXYZ[2] : 1) : 1); /*pow(genCtx.blockSizeElement,genCtx.dim);*/ 260ce27a4eeSVijay Mahadevan PetscInt nelems = ntensorelems; 26163d025dbSVijay Mahadevan moab::EntityHandle starte; /* connectivity */ 262ce27a4eeSVijay Mahadevan moab::EntityHandle* conn; 263ce27a4eeSVijay Mahadevan 264ce27a4eeSVijay Mahadevan PetscFunctionBegin; 265ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 266ce27a4eeSVijay Mahadevan case 1: 267ce27a4eeSVijay Mahadevan nvperelem = 2; 268ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 2, moab::MBEDGE, 0, starte, conn);MBERR("Can't get EDGE2 element connectivity.", merr); 269ce27a4eeSVijay Mahadevan break; 270ce27a4eeSVijay Mahadevan case 2: 271ce27a4eeSVijay Mahadevan if (genCtx.simplex) { 272ce27a4eeSVijay Mahadevan nvperelem = 3; 27349d66b22SVijay Mahadevan nelems = ntensorelems * simplices_per_tensor[genCtx.dim]; 274ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 3, moab::MBTRI, 0, starte, conn);MBERR("Can't get TRI3 element connectivity.", merr); 275*1baa6e33SBarry Smith } else { 276ce27a4eeSVijay Mahadevan nvperelem = 4; 277ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 4, moab::MBQUAD, 0, starte, conn);MBERR("Can't get QUAD4 element connectivity.", merr); 278ce27a4eeSVijay Mahadevan } 279ce27a4eeSVijay Mahadevan break; 280ce27a4eeSVijay Mahadevan case 3: 281ce27a4eeSVijay Mahadevan default: 282ce27a4eeSVijay Mahadevan if (genCtx.simplex) { 283ce27a4eeSVijay Mahadevan nvperelem = 4; 28449d66b22SVijay Mahadevan nelems = ntensorelems * simplices_per_tensor[genCtx.dim]; 285ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 4, moab::MBTET, 0, starte, conn);MBERR("Can't get TET4 element connectivity.", merr); 286*1baa6e33SBarry Smith } else { 287ce27a4eeSVijay Mahadevan nvperelem = 8; 288ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 8, moab::MBHEX, 0, starte, conn);MBERR("Can't get HEX8 element connectivity.", merr); 289ce27a4eeSVijay Mahadevan } 290ce27a4eeSVijay Mahadevan break; 291ce27a4eeSVijay Mahadevan } 292ce27a4eeSVijay Mahadevan 29363d025dbSVijay Mahadevan ix = ie = 0; /* index now in the elements, for global ids */ 294ce27a4eeSVijay Mahadevan 295ce27a4eeSVijay Mahadevan /* create a temporary range to store local element handles */ 296ce27a4eeSVijay Mahadevan moab::Range tmp(starte, starte + nelems - 1); 29749d66b22SVijay Mahadevan std::vector<PetscInt> gids(nelems); 298ce27a4eeSVijay Mahadevan 299ce27a4eeSVijay Mahadevan /* identify the elements at the lower corner, for their global ids */ 300ce27a4eeSVijay Mahadevan xe = m * genCtx.A * genCtx.blockSizeElementXYZ[0] + a * genCtx.blockSizeElementXYZ[0]; 301ce27a4eeSVijay Mahadevan ye = (genCtx.dim > 1 ? n * genCtx.B * genCtx.blockSizeElementXYZ[1] + b * genCtx.blockSizeElementXYZ[1] : 0); 302ce27a4eeSVijay Mahadevan ze = (genCtx.dim > 2 ? k * genCtx.C * genCtx.blockSizeElementXYZ[2] + c * genCtx.blockSizeElementXYZ[2] : 0); 303ce27a4eeSVijay Mahadevan 304ce27a4eeSVijay Mahadevan /* create owned elements requested by genCtx */ 305ce27a4eeSVijay Mahadevan for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeElementXYZ[2] : 1); kk++) { 306ce27a4eeSVijay Mahadevan for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeElementXYZ[1] : 1); jj++) { 307ce27a4eeSVijay Mahadevan for (ii = 0; ii < genCtx.blockSizeElementXYZ[0]; ii++) { 308ce27a4eeSVijay Mahadevan 309ce27a4eeSVijay Mahadevan moab::EntityHandle corner = startv + genCtx.q * ii + genCtx.q * jj * genCtx.ystride + genCtx.q * kk * genCtx.zstride; 310ce27a4eeSVijay Mahadevan 31149d66b22SVijay Mahadevan std::pair<PetscInt, PetscInt> entoffset = DMMoab_SetElementConnectivity_Private(genCtx, ix, corner, conn); 312ce27a4eeSVijay Mahadevan 31349d66b22SVijay Mahadevan for (PetscInt j = 0; j < entoffset.second; j++) { 314ce27a4eeSVijay Mahadevan /* The entity handle for the particular element -> if we want to set some tags is 315ce27a4eeSVijay Mahadevan moab::EntityHandle eh = starte + ie + j; 316ce27a4eeSVijay Mahadevan */ 317ce27a4eeSVijay Mahadevan gids[ie + j] = 1 + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)); 31863d025dbSVijay Mahadevan /* gids[ie+j] = ie + j + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)); */ 31963d025dbSVijay Mahadevan /* gids[ie+j] = 1 + ie; */ 32063d025dbSVijay Mahadevan /* ie++; */ 321ce27a4eeSVijay Mahadevan } 322ce27a4eeSVijay Mahadevan 323ce27a4eeSVijay Mahadevan ix += entoffset.first; 32449d66b22SVijay Mahadevan ie += entoffset.second; 325ce27a4eeSVijay Mahadevan } 326ce27a4eeSVijay Mahadevan } 327ce27a4eeSVijay Mahadevan } 328ce27a4eeSVijay Mahadevan if (genCtx.adjEnts) { /* we need to update adjacencies now, because some elements are new */ 329ce27a4eeSVijay Mahadevan merr = iface->update_adjacencies(starte, nelems, nvperelem, conn);MBERR("Can't update adjacencies", merr); 330ce27a4eeSVijay Mahadevan } 331ce27a4eeSVijay Mahadevan tmp.swap(cells); 332ce27a4eeSVijay Mahadevan merr = mbImpl->tag_set_data(global_id_tag, cells, &gids[0]);MBERR("Can't set global ids to elements.", merr); 333ce27a4eeSVijay Mahadevan PetscFunctionReturn(0); 334ce27a4eeSVijay Mahadevan } 335ce27a4eeSVijay Mahadevan 336755f3dfbSVijay Mahadevan PetscErrorCode DMMBUtil_InitializeOptions(DMMoabMeshGeneratorCtx& genCtx, PetscInt dim, PetscBool simplex, PetscInt rank, PetscInt nprocs, const PetscReal* bounds, PetscInt nelems) 337ce27a4eeSVijay Mahadevan { 338ce27a4eeSVijay Mahadevan PetscFunctionBegin; 339ce27a4eeSVijay Mahadevan /* Initialize all genCtx data */ 340ce27a4eeSVijay Mahadevan genCtx.dim = dim; 341ce27a4eeSVijay Mahadevan genCtx.simplex = simplex; 34249d66b22SVijay Mahadevan genCtx.newMergeMethod = genCtx.keep_skins = genCtx.adjEnts = true; 343ce27a4eeSVijay Mahadevan /* determine other global quantities for the mesh used for nodes increments */ 344ce27a4eeSVijay Mahadevan genCtx.q = 1; 345755f3dfbSVijay Mahadevan genCtx.fraction = genCtx.remainder = genCtx.cumfraction = 0; 346ce27a4eeSVijay Mahadevan 347ce27a4eeSVijay Mahadevan if (!genCtx.usrxyzgrid) { /* not overridden by genCtx - assume nele equally and that genCtx wants a uniform cube mesh */ 348755f3dfbSVijay Mahadevan 349755f3dfbSVijay Mahadevan genCtx.fraction = nelems / nprocs; /* partition only by the largest dimension */ 350755f3dfbSVijay Mahadevan genCtx.remainder = nelems % nprocs; /* remainder after partition which gets evenly distributed by round-robin */ 351755f3dfbSVijay Mahadevan genCtx.cumfraction = (rank > 0 ? (genCtx.fraction) * (rank) + (rank - 1 < genCtx.remainder ? rank : genCtx.remainder) : 0); 352755f3dfbSVijay Mahadevan if (rank < genCtx.remainder) /* This process gets "fraction+1" elements */ 353755f3dfbSVijay Mahadevan genCtx.fraction++; 354755f3dfbSVijay Mahadevan 3559566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Fraction = %" PetscInt_FMT ", Remainder = %" PetscInt_FMT ", Cumulative fraction = %" PetscInt_FMT "\n", genCtx.fraction, genCtx.remainder, genCtx.cumfraction)); 356755f3dfbSVijay Mahadevan switch (genCtx.dim) { 357755f3dfbSVijay Mahadevan case 1: 358755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = genCtx.fraction; 359755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = 1; 360755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = 1; 361755f3dfbSVijay Mahadevan break; 362755f3dfbSVijay Mahadevan case 2: 363755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = nelems; 364755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = genCtx.fraction; 365755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = 1; 366755f3dfbSVijay Mahadevan break; 367755f3dfbSVijay Mahadevan case 3: 368755f3dfbSVijay Mahadevan default: 369755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = nelems; 370755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = nelems; 371755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = genCtx.fraction; 372755f3dfbSVijay Mahadevan break; 37349d66b22SVijay Mahadevan } 374ce27a4eeSVijay Mahadevan } 375ce27a4eeSVijay Mahadevan 376ce27a4eeSVijay Mahadevan /* partition only by the largest dimension */ 377ce27a4eeSVijay Mahadevan /* Total number of local elements := genCtx.blockSizeElementXYZ[0]*(genCtx.dim>1? genCtx.blockSizeElementXYZ[1]*(genCtx.dim>2 ? genCtx.blockSizeElementXYZ[2]:1) :1); */ 378ce27a4eeSVijay Mahadevan if (bounds) { 3795f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < 6; i++) genCtx.xyzbounds[i] = bounds[i]; 380*1baa6e33SBarry Smith } else { 381ce27a4eeSVijay Mahadevan genCtx.xyzbounds[0] = genCtx.xyzbounds[2] = genCtx.xyzbounds[4] = 0.0; 382ce27a4eeSVijay Mahadevan genCtx.xyzbounds[1] = genCtx.xyzbounds[3] = genCtx.xyzbounds[5] = 1.0; 383ce27a4eeSVijay Mahadevan } 384ce27a4eeSVijay Mahadevan 385ce27a4eeSVijay Mahadevan if (!genCtx.usrprocgrid) { 386ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 387ce27a4eeSVijay Mahadevan case 1: 388ce27a4eeSVijay Mahadevan genCtx.M = nprocs; 389ce27a4eeSVijay Mahadevan genCtx.N = genCtx.K = 1; 390ce27a4eeSVijay Mahadevan break; 391ce27a4eeSVijay Mahadevan case 2: 392755f3dfbSVijay Mahadevan genCtx.N = nprocs; 39363d025dbSVijay Mahadevan genCtx.M = genCtx.K = 1; 394ce27a4eeSVijay Mahadevan break; 395ce27a4eeSVijay Mahadevan default: 396755f3dfbSVijay Mahadevan genCtx.K = nprocs; 39763d025dbSVijay Mahadevan genCtx.M = genCtx.N = 1; 398ce27a4eeSVijay Mahadevan break; 399ce27a4eeSVijay Mahadevan } 400ce27a4eeSVijay Mahadevan } 401ce27a4eeSVijay Mahadevan 402ce27a4eeSVijay Mahadevan if (!genCtx.usrrefgrid) { 403ce27a4eeSVijay Mahadevan genCtx.A = genCtx.B = genCtx.C = 1; 404ce27a4eeSVijay Mahadevan } 405ce27a4eeSVijay Mahadevan 406ce27a4eeSVijay Mahadevan /* more default values */ 407ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.ney = genCtx.nez = 0; 408ce27a4eeSVijay Mahadevan genCtx.xstride = genCtx.ystride = genCtx.zstride = 0; 409ce27a4eeSVijay Mahadevan genCtx.NX = genCtx.NY = genCtx.NZ = 0; 410ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.ney = genCtx.nez = 0; 411ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 1; 412ce27a4eeSVijay Mahadevan 413ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 414ce27a4eeSVijay Mahadevan case 3: 415ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 416ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1; 417ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[2] = genCtx.q * genCtx.blockSizeElementXYZ[2] + 1; 418ce27a4eeSVijay Mahadevan 419ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 420755f3dfbSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */ 421ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1); 422ce27a4eeSVijay Mahadevan genCtx.xstride = 1; 423ce27a4eeSVijay Mahadevan genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1]; /* number of elements in y direction .... */ 424755f3dfbSVijay Mahadevan genCtx.dy = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */ 425ce27a4eeSVijay Mahadevan genCtx.NY = (genCtx.q * genCtx.ney + 1); 426ce27a4eeSVijay Mahadevan genCtx.ystride = genCtx.blockSizeVertexXYZ[0]; 427ce27a4eeSVijay Mahadevan genCtx.nez = genCtx.K * genCtx.C * genCtx.blockSizeElementXYZ[2]; /* number of elements in z direction .... */ 428755f3dfbSVijay Mahadevan genCtx.dz = (genCtx.xyzbounds[5] - genCtx.xyzbounds[4]) / (nelems * genCtx.q); /* distance between 2 nodes in z direction */ 429ce27a4eeSVijay Mahadevan genCtx.NZ = (genCtx.q * genCtx.nez + 1); 430ce27a4eeSVijay Mahadevan genCtx.zstride = genCtx.blockSizeVertexXYZ[0] * genCtx.blockSizeVertexXYZ[1]; 431ce27a4eeSVijay Mahadevan break; 432ce27a4eeSVijay Mahadevan case 2: 433ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 434ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1; 4359daf19fdSVijay Mahadevan genCtx.blockSizeVertexXYZ[2] = 0; 436ce27a4eeSVijay Mahadevan 437ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 438ce27a4eeSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (genCtx.nex * genCtx.q); /* distance between 2 nodes in x direction */ 439ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1); 440ce27a4eeSVijay Mahadevan genCtx.xstride = 1; 441ce27a4eeSVijay Mahadevan genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1]; /* number of elements in y direction .... */ 442755f3dfbSVijay Mahadevan genCtx.dy = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */ 443ce27a4eeSVijay Mahadevan genCtx.NY = (genCtx.q * genCtx.ney + 1); 444ce27a4eeSVijay Mahadevan genCtx.ystride = genCtx.blockSizeVertexXYZ[0]; 445ce27a4eeSVijay Mahadevan break; 446ce27a4eeSVijay Mahadevan case 1: 4479daf19fdSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 0; 448ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 449ce27a4eeSVijay Mahadevan 450ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 451755f3dfbSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */ 452ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1); 453ce27a4eeSVijay Mahadevan genCtx.xstride = 1; 454ce27a4eeSVijay Mahadevan break; 455ce27a4eeSVijay Mahadevan } 456ce27a4eeSVijay Mahadevan 457755f3dfbSVijay Mahadevan /* Lets check for some valid input */ 4585f80ce2aSJacob Faibussowitsch PetscCheck(genCtx.dim >= 1 && genCtx.dim <= 3,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid topological dimension specified: %" PetscInt_FMT ".", genCtx.dim); 4595f80ce2aSJacob Faibussowitsch 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); 460755f3dfbSVijay Mahadevan /* validate the bounds data */ 4615f80ce2aSJacob Faibussowitsch 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]); 4621dca8a05SBarry Smith 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]); 4631dca8a05SBarry Smith 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]); 464755f3dfbSVijay Mahadevan 4659566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local elements:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.blockSizeElementXYZ[0], genCtx.blockSizeElementXYZ[1], genCtx.blockSizeElementXYZ[2])); 4669566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local vertices:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.blockSizeVertexXYZ[0], genCtx.blockSizeVertexXYZ[1], genCtx.blockSizeVertexXYZ[2])); 4679566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local blocks/processors := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.A, genCtx.B, genCtx.C)); 4689566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local processors := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.M, genCtx.N, genCtx.K)); 4699566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local nexyz:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.nex, genCtx.ney, genCtx.nez)); 4709566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local delxyz:= %g, %g, %g\n", genCtx.dx, genCtx.dy, genCtx.dz)); 4719566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local strides:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.xstride, genCtx.ystride, genCtx.zstride)); 472ce27a4eeSVijay Mahadevan PetscFunctionReturn(0); 473ce27a4eeSVijay Mahadevan } 474ce27a4eeSVijay Mahadevan 475cab5ea25SPierre Jolivet /*@C 476ce27a4eeSVijay Mahadevan DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with genCtx specified bounds. 477aa859218SVijay Mahadevan 478d083f849SBarry Smith Collective 479aa859218SVijay Mahadevan 480aa859218SVijay Mahadevan Input Parameters: 481aa859218SVijay Mahadevan + comm - The communicator for the DM object 482aa859218SVijay Mahadevan . dim - The spatial dimension 483b8ecf6d3SVijay Mahadevan . bounds - The bounds of the box specified with [x-left, x-right, y-bottom, y-top, z-bottom, z-top] depending on the spatial dimension 484b8ecf6d3SVijay Mahadevan . nele - The number of discrete elements in each direction 485a2b725a8SWilliam Gropp - user_nghost - The number of ghosted layers needed in the partitioned mesh 486aa859218SVijay Mahadevan 487aa859218SVijay Mahadevan Output Parameter: 488aa859218SVijay Mahadevan . dm - The DM object 489aa859218SVijay Mahadevan 490aa859218SVijay Mahadevan Level: beginner 491aa859218SVijay Mahadevan 492db781477SPatrick Sanan .seealso: `DMSetType()`, `DMCreate()`, `DMMoabLoadFromFile()` 493aa859218SVijay Mahadevan @*/ 494ce27a4eeSVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal* bounds, PetscInt nele, PetscInt nghost, DM *dm) 49551d15aeeSVijay Mahadevan { 49651d15aeeSVijay Mahadevan moab::ErrorCode merr; 4979daf19fdSVijay Mahadevan PetscInt a, b, c, n, global_size, global_rank; 49851d15aeeSVijay Mahadevan DM_Moab *dmmoab; 499ce27a4eeSVijay Mahadevan moab::Interface *mbImpl; 5009daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 50151d15aeeSVijay Mahadevan moab::ParallelComm *pcomm; 5029daf19fdSVijay Mahadevan #endif 503a4717116SVijay Mahadevan moab::ReadUtilIface *readMeshIface; 504ce27a4eeSVijay Mahadevan moab::Range verts, cells, edges, faces, adj, dim3, dim2; 505ce27a4eeSVijay Mahadevan DMMoabMeshGeneratorCtx genCtx; 506a4717116SVijay Mahadevan const PetscInt npts = nele + 1; /* Number of points in every dimension */ 507ce27a4eeSVijay Mahadevan 5086ea892caSVijay Mahadevan moab::Tag global_id_tag, part_tag, geom_tag, mat_tag, dir_tag, neu_tag; 509ce27a4eeSVijay Mahadevan moab::Range ownedvtx, ownedelms, localvtxs, localelms; 510ce27a4eeSVijay Mahadevan moab::EntityHandle regionset; 511755f3dfbSVijay Mahadevan PetscInt ml = 0, nl = 0, kl = 0; 51251d15aeeSVijay Mahadevan 51351d15aeeSVijay Mahadevan PetscFunctionBegin; 5141dca8a05SBarry Smith PetscCheck(dim >= 1 && dim <= 3,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension argument for mesh: dim=[1,3]."); 515e5136372SVijay Mahadevan 5169566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("GenerateMesh", DM_CLASSID, &genCtx.generateMesh)); 5179566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("AddVertices", DM_CLASSID, &genCtx.generateVertices)); 5189566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("AddElements", DM_CLASSID, &genCtx.generateElements)); 5199566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("ParResolve", DM_CLASSID, &genCtx.parResolve)); 5209566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(genCtx.generateMesh, 0, 0, 0, 0)); 5219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &global_size)); 522e5136372SVijay Mahadevan /* total number of vertices in all dimensions */ 523e5136372SVijay Mahadevan n = pow(npts, dim); 524e5136372SVijay Mahadevan 525e5136372SVijay Mahadevan /* do some error checking */ 52608401ef6SPierre Jolivet PetscCheck(n >= 2,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be >= 2."); 52708401ef6SPierre Jolivet PetscCheck(global_size <= n,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of processors must be less than or equal to number of elements."); 52808401ef6SPierre Jolivet PetscCheck(nghost >= 0,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of ghost layers cannot be negative."); 529e5136372SVijay Mahadevan 530a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 5319566063dSJacob Faibussowitsch PetscCall(DMMoabCreateMoab(comm, NULL, NULL, NULL, dm)); 53251d15aeeSVijay Mahadevan 533a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 53451d15aeeSVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 535ce27a4eeSVijay Mahadevan mbImpl = dmmoab->mbiface; 5369daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 53751d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm; 538ce27a4eeSVijay Mahadevan global_rank = pcomm->rank(); 5399daf19fdSVijay Mahadevan #else 5409daf19fdSVijay Mahadevan global_rank = 0; 5419daf19fdSVijay Mahadevan global_size = 1; 5429daf19fdSVijay Mahadevan #endif 5439daf19fdSVijay Mahadevan global_id_tag = dmmoab->ltog_tag; 544c8d5365dSVijay Mahadevan dmmoab->dim = dim; 54549d66b22SVijay Mahadevan dmmoab->nghostrings = nghost; 546304006b3SVijay Mahadevan dmmoab->refct = 1; 54751d15aeeSVijay Mahadevan 548b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */ 549ce27a4eeSVijay Mahadevan merr = mbImpl->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 550b5410836SVijay Mahadevan 551a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */ 552ce27a4eeSVijay Mahadevan merr = mbImpl->query_interface(readMeshIface);MBERRNM(merr); 55351d15aeeSVijay Mahadevan 554755f3dfbSVijay Mahadevan genCtx.M = genCtx.N = genCtx.K = 1; 555755f3dfbSVijay Mahadevan genCtx.A = genCtx.B = genCtx.C = 1; 556da8b1a3eSBarry Smith genCtx.blockSizeElementXYZ[0] = 0; 557da8b1a3eSBarry Smith genCtx.blockSizeElementXYZ[1] = 0; 558da8b1a3eSBarry Smith genCtx.blockSizeElementXYZ[2] = 0; 559755f3dfbSVijay Mahadevan 560d0609cedSBarry Smith PetscOptionsBegin(comm, "", "DMMoab Creation Options", "DMMOAB"); 561ce27a4eeSVijay Mahadevan /* Handle DMMoab spatial resolution */ 5629566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dmb_grid_x", "Number of grid points in x direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[0], &genCtx.blockSizeElementXYZ[0], &genCtx.usrxyzgrid)); 5639566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsInt("-dmb_grid_y", "Number of grid points in y direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[1], &genCtx.blockSizeElementXYZ[1], &genCtx.usrxyzgrid)); 5649566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsInt("-dmb_grid_z", "Number of grid points in z direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[2], &genCtx.blockSizeElementXYZ[2], &genCtx.usrxyzgrid)); 565a4717116SVijay Mahadevan 5666aad120cSJose E. Roman /* Handle DMMoab parallel distribution */ 5679566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dmb_processors_x", "Number of processors in x direction", "DMMoabSetNumProcs", genCtx.M, &genCtx.M, &genCtx.usrprocgrid)); 5689566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsInt("-dmb_processors_y", "Number of processors in y direction", "DMMoabSetNumProcs", genCtx.N, &genCtx.N, &genCtx.usrprocgrid)); 5699566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsInt("-dmb_processors_z", "Number of processors in z direction", "DMMoabSetNumProcs", genCtx.K, &genCtx.K, &genCtx.usrprocgrid)); 57051d15aeeSVijay Mahadevan 571ce27a4eeSVijay Mahadevan /* Handle DMMoab block refinement */ 5729566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dmb_refine_x", "Number of refinement blocks in x direction", "DMMoabSetRefinement", genCtx.A, &genCtx.A, &genCtx.usrrefgrid)); 5739566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsInt("-dmb_refine_y", "Number of refinement blocks in y direction", "DMMoabSetRefinement", genCtx.B, &genCtx.B, &genCtx.usrrefgrid)); 5749566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsInt("-dmb_refine_z", "Number of refinement blocks in z direction", "DMMoabSetRefinement", genCtx.C, &genCtx.C, &genCtx.usrrefgrid)); 575d0609cedSBarry Smith PetscOptionsEnd(); 57651d15aeeSVijay Mahadevan 5779566063dSJacob Faibussowitsch PetscCall(DMMBUtil_InitializeOptions(genCtx, dim, useSimplex, global_rank, global_size, bounds, nele)); 57851d15aeeSVijay Mahadevan 57963a3b9bcSJacob Faibussowitsch //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); 58064e1c140SVijay Mahadevan 581ce27a4eeSVijay Mahadevan if (genCtx.adjEnts) genCtx.keep_skins = true; /* do not delete anything - consumes more memory */ 58266f88a78SVijay Mahadevan 583ce27a4eeSVijay Mahadevan /* determine m, n, k for processor rank */ 584755f3dfbSVijay Mahadevan ml = nl = kl = 0; 585755f3dfbSVijay Mahadevan switch (genCtx.dim) { 586755f3dfbSVijay Mahadevan case 1: 587755f3dfbSVijay Mahadevan ml = (genCtx.cumfraction); 588755f3dfbSVijay Mahadevan break; 589755f3dfbSVijay Mahadevan case 2: 590755f3dfbSVijay Mahadevan nl = (genCtx.cumfraction); 591755f3dfbSVijay Mahadevan break; 592755f3dfbSVijay Mahadevan default: 593755f3dfbSVijay Mahadevan kl = (genCtx.cumfraction) / genCtx.q / genCtx.blockSizeElementXYZ[2] / genCtx.C; //genCtx.K 594755f3dfbSVijay Mahadevan break; 595755f3dfbSVijay Mahadevan } 596e5136372SVijay Mahadevan 597ce27a4eeSVijay Mahadevan /* 598ce27a4eeSVijay Mahadevan * so there are a total of M * A * blockSizeElement elements in x direction (so M * A * blockSizeElement + 1 verts in x direction) 599ce27a4eeSVijay Mahadevan * so there are a total of N * B * blockSizeElement elements in y direction (so N * B * blockSizeElement + 1 verts in y direction) 600ce27a4eeSVijay Mahadevan * so there are a total of K * C * blockSizeElement elements in z direction (so K * C * blockSizeElement + 1 verts in z direction) 601c8d5365dSVijay Mahadevan 602ce27a4eeSVijay Mahadevan * there are ( M * A blockSizeElement) * ( N * B * blockSizeElement) * (K * C * blockSizeElement) hexas 603ce27a4eeSVijay Mahadevan * there are ( M * A * blockSizeElement + 1) * ( N * B * blockSizeElement + 1) * (K * C * blockSizeElement + 1) vertices 604ce27a4eeSVijay Mahadevan * x is the first dimension that varies 605ce27a4eeSVijay Mahadevan */ 606e5136372SVijay Mahadevan 607ce27a4eeSVijay Mahadevan /* generate the block at (a, b, c); it will represent a partition , it will get a partition tag */ 60849d66b22SVijay Mahadevan PetscInt dum_id = -1; 6096ea892caSVijay Mahadevan merr = mbImpl->tag_get_handle("GLOBAL_ID", 1, moab::MB_TYPE_INTEGER, global_id_tag);MBERR("Getting Global_ID Tag handle failed", merr); 61051d15aeeSVijay Mahadevan 6116ea892caSVijay Mahadevan merr = mbImpl->tag_get_handle(MATERIAL_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, mat_tag);MBERR("Getting Material set Tag handle failed", merr); 6126ea892caSVijay Mahadevan merr = mbImpl->tag_get_handle(DIRICHLET_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, dir_tag);MBERR("Getting Dirichlet set Tag handle failed", merr); 6136ea892caSVijay Mahadevan merr = mbImpl->tag_get_handle(NEUMANN_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, neu_tag);MBERR("Getting Neumann set Tag handle failed", merr); 6146ea892caSVijay Mahadevan 6156ea892caSVijay Mahadevan 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); 61651d15aeeSVijay Mahadevan 6173a4aeca1SVijay Mahadevan /* lets create some sets */ 618ce27a4eeSVijay Mahadevan 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); 61949d66b22SVijay Mahadevan merr = mbImpl->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr); 6209566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(genCtx.generateMesh, 0, 0, 0, 0)); 6213a4aeca1SVijay Mahadevan 622ce27a4eeSVijay Mahadevan for (a = 0; a < (genCtx.dim > 0 ? genCtx.A : genCtx.A); a++) { 623ce27a4eeSVijay Mahadevan for (b = 0; b < (genCtx.dim > 1 ? genCtx.B : 1); b++) { 624ce27a4eeSVijay Mahadevan for (c = 0; c < (genCtx.dim > 2 ? genCtx.C : 1); c++) { 62549d66b22SVijay Mahadevan 62649d66b22SVijay Mahadevan moab::EntityHandle startv; 62749d66b22SVijay Mahadevan 6289566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(genCtx.generateVertices, 0, 0, 0, 0)); 6299566063dSJacob Faibussowitsch PetscCall(DMMoab_GenerateVertices_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, verts)); 6309566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(genCtx.generateVertices, 0, 0, 0, 0)); 6313a4aeca1SVijay Mahadevan 6329566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(genCtx.generateElements, 0, 0, 0, 0)); 6339566063dSJacob Faibussowitsch PetscCall(DMMoab_GenerateElements_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, cells)); 6349566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(genCtx.generateElements, 0, 0, 0, 0)); 6353a4aeca1SVijay Mahadevan 63649d66b22SVijay Mahadevan PetscInt part_num = 0; 63749d66b22SVijay Mahadevan switch (genCtx.dim) { 63849d66b22SVijay Mahadevan case 3: 63949d66b22SVijay Mahadevan part_num += (c + kl * genCtx.C) * (genCtx.M * genCtx.A * genCtx.N * genCtx.B); 64049d66b22SVijay Mahadevan case 2: 64149d66b22SVijay Mahadevan part_num += (b + nl * genCtx.B) * (genCtx.M * genCtx.A); 64249d66b22SVijay Mahadevan case 1: 64349d66b22SVijay Mahadevan part_num += (a + ml * genCtx.A); 64449d66b22SVijay Mahadevan break; 64549d66b22SVijay Mahadevan } 64649d66b22SVijay Mahadevan 64749d66b22SVijay Mahadevan moab::EntityHandle part_set; 648ce27a4eeSVijay Mahadevan merr = mbImpl->create_meshset(moab::MESHSET_SET, part_set);MBERR("Can't create mesh set.", merr); 6493a4aeca1SVijay Mahadevan 650ce27a4eeSVijay Mahadevan merr = mbImpl->add_entities(part_set, verts);MBERR("Can't add vertices to set.", merr); 65149d66b22SVijay Mahadevan merr = mbImpl->add_entities(part_set, cells);MBERR("Can't add entities to set.", merr); 652ce27a4eeSVijay Mahadevan merr = mbImpl->add_entities(regionset, cells);MBERR("Can't add entities to set.", merr); 653ce27a4eeSVijay Mahadevan 654ce27a4eeSVijay Mahadevan /* if needed, add all edges and faces */ 655ce27a4eeSVijay Mahadevan if (genCtx.adjEnts) 656ce27a4eeSVijay Mahadevan { 657ce27a4eeSVijay Mahadevan if (genCtx.dim > 1) { 658ce27a4eeSVijay Mahadevan merr = mbImpl->get_adjacencies(cells, 1, true, edges, moab::Interface::UNION);MBERR("Can't get edges", merr); 659ce27a4eeSVijay Mahadevan merr = mbImpl->add_entities(part_set, edges);MBERR("Can't add edges to partition set.", merr); 660ce27a4eeSVijay Mahadevan } 661ce27a4eeSVijay Mahadevan if (genCtx.dim > 2) { 662ce27a4eeSVijay Mahadevan merr = mbImpl->get_adjacencies(cells, 2, true, faces, moab::Interface::UNION);MBERR("Can't get faces", merr); 663ce27a4eeSVijay Mahadevan merr = mbImpl->add_entities(part_set, faces);MBERR("Can't add faces to partition set.", merr); 664ce27a4eeSVijay Mahadevan } 665ce27a4eeSVijay Mahadevan edges.clear(); 666ce27a4eeSVijay Mahadevan faces.clear(); 667ce27a4eeSVijay Mahadevan } 66849d66b22SVijay Mahadevan verts.clear(); cells.clear(); 669ce27a4eeSVijay Mahadevan 670ce27a4eeSVijay Mahadevan merr = mbImpl->tag_set_data(part_tag, &part_set, 1, &part_num);MBERR("Can't set part tag on set", merr); 67149d66b22SVijay Mahadevan if (dmmoab->fileset) { 672ce27a4eeSVijay Mahadevan merr = mbImpl->add_parent_child(dmmoab->fileset, part_set);MBERR("Can't add part set to file set.", merr); 673ce27a4eeSVijay Mahadevan merr = mbImpl->unite_meshset(dmmoab->fileset, part_set);MBERRNM(merr); 67449d66b22SVijay Mahadevan } 67549d66b22SVijay Mahadevan merr = mbImpl->add_entities(dmmoab->fileset, &part_set, 1);MBERRNM(merr); 676ce27a4eeSVijay Mahadevan } 677ce27a4eeSVijay Mahadevan } 678ce27a4eeSVijay Mahadevan } 679ce27a4eeSVijay Mahadevan 68049d66b22SVijay Mahadevan merr = mbImpl->add_parent_child(dmmoab->fileset, regionset);MBERRNM(merr); 681ce27a4eeSVijay Mahadevan 68249d66b22SVijay Mahadevan /* Only in parallel: resolve shared entities between processors and exchange ghost layers */ 683ce27a4eeSVijay Mahadevan if (global_size > 1) { 684ce27a4eeSVijay Mahadevan 6859566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(genCtx.parResolve, 0, 0, 0, 0)); 68677a8c067SVijay Mahadevan 687ce27a4eeSVijay Mahadevan merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, genCtx.dim, cells);MBERR("Can't get all d-dimensional elements.", merr); 688ce27a4eeSVijay Mahadevan merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 0, verts);MBERR("Can't get all vertices.", merr); 689ce27a4eeSVijay Mahadevan 690ce27a4eeSVijay Mahadevan if (genCtx.A * genCtx.B * genCtx.C != 1) { // merge needed 69149d66b22SVijay Mahadevan moab::MergeMesh mm(mbImpl); 692ce27a4eeSVijay Mahadevan if (genCtx.newMergeMethod) { 693ce27a4eeSVijay Mahadevan merr = mm.merge_using_integer_tag(verts, global_id_tag);MBERR("Can't merge with GLOBAL_ID tag", merr); 694ce27a4eeSVijay Mahadevan } 695ce27a4eeSVijay Mahadevan else { 696ce27a4eeSVijay Mahadevan merr = mm.merge_entities(cells, 0.0001);MBERR("Can't merge with coordinates", merr); 6973a4aeca1SVijay Mahadevan } 6983a4aeca1SVijay Mahadevan } 6993a4aeca1SVijay Mahadevan 7009daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 701a4717116SVijay Mahadevan /* check the handles */ 702ce27a4eeSVijay Mahadevan merr = pcomm->check_all_shared_handles();MBERRV(mbImpl, merr); 70351d15aeeSVijay Mahadevan 704a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */ 705ce27a4eeSVijay Mahadevan merr = pcomm->resolve_shared_ents(dmmoab->fileset, cells, dim, dim - 1, NULL, &global_id_tag);MBERRV(mbImpl, merr); 70649d66b22SVijay Mahadevan if (dmmoab->fileset) { 707ce27a4eeSVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false, &dmmoab->fileset);MBERRV(mbImpl, merr); 70849d66b22SVijay Mahadevan } 70949d66b22SVijay Mahadevan else { 71049d66b22SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false);MBERRV(mbImpl, merr); 71149d66b22SVijay Mahadevan } 712c8d5365dSVijay Mahadevan 713e427d9c9SVijay Mahadevan /* Reassign global IDs on all entities. */ 714e427d9c9SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, false, true, false);MBERRNM(merr); 7159daf19fdSVijay Mahadevan #endif 71677a8c067SVijay Mahadevan 7179566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(genCtx.parResolve, 0, 0, 0, 0)); 71849d66b22SVijay Mahadevan } 7193f1c6e43SVijay Mahadevan 720ce27a4eeSVijay Mahadevan if (!genCtx.keep_skins) { // default is to delete the 1- and 2-dimensional entities 721ce27a4eeSVijay Mahadevan // delete all quads and edges 722ce27a4eeSVijay Mahadevan moab::Range toDelete; 723ce27a4eeSVijay Mahadevan if (genCtx.dim > 1) { 724ce27a4eeSVijay Mahadevan merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 1, toDelete);MBERR("Can't get edges", merr); 725ce27a4eeSVijay Mahadevan } 7263f1c6e43SVijay Mahadevan 727ce27a4eeSVijay Mahadevan if (genCtx.dim > 2) { 728ce27a4eeSVijay Mahadevan merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 2, toDelete);MBERR("Can't get faces", merr); 729ce27a4eeSVijay Mahadevan } 730c8d5365dSVijay Mahadevan 7319daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 732ce27a4eeSVijay Mahadevan merr = dmmoab->pcomm->delete_entities(toDelete) ;MBERR("Can't delete entities", merr); 7339daf19fdSVijay Mahadevan #endif 734ce27a4eeSVijay Mahadevan } 7356ea892caSVijay Mahadevan 7366ea892caSVijay Mahadevan /* set geometric dimension tag for regions */ 7376ea892caSVijay Mahadevan merr = mbImpl->tag_set_data(geom_tag, ®ionset, 1, &dmmoab->dim);MBERRNM(merr); 7386ea892caSVijay Mahadevan /* set default material ID for regions */ 7396ea892caSVijay Mahadevan int default_material = 1; 7406ea892caSVijay Mahadevan merr = mbImpl->tag_set_data(mat_tag, ®ionset, 1, &default_material);MBERRNM(merr); 7416ea892caSVijay Mahadevan /* 7426ea892caSVijay Mahadevan int default_dbc = 0; 7436ea892caSVijay Mahadevan merr = mbImpl->tag_set_data(dir_tag, &vertexset, 1, &default_dbc);MBERRNM(merr); 7446ea892caSVijay Mahadevan */ 74551d15aeeSVijay Mahadevan PetscFunctionReturn(0); 74651d15aeeSVijay Mahadevan } 74751d15aeeSVijay Mahadevan 74849d66b22SVijay Mahadevan 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) 74951d15aeeSVijay Mahadevan { 7506acfe860SVijay Mahadevan char *ropts; 75149d66b22SVijay Mahadevan char ropts_par[PETSC_MAX_PATH_LEN], ropts_pargh[PETSC_MAX_PATH_LEN]; 75261eb6e27SVijay Mahadevan char ropts_dbg[PETSC_MAX_PATH_LEN]; 75351d15aeeSVijay Mahadevan 754a4717116SVijay Mahadevan PetscFunctionBegin; 7559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PETSC_MAX_PATH_LEN, &ropts)); 7569566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ropts_par, PETSC_MAX_PATH_LEN)); 7579566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ropts_pargh, PETSC_MAX_PATH_LEN)); 7589566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ropts_dbg, PETSC_MAX_PATH_LEN)); 75961eb6e27SVijay Mahadevan 760e23c60ebSVijay Mahadevan /* do parallel read unless using only one processor */ 761a4717116SVijay Mahadevan if (numproc > 1) { 7629566063dSJacob Faibussowitsch // 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":""))); 7639566063dSJacob Faibussowitsch 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;" : ""))); 76449d66b22SVijay Mahadevan if (nghost) { 7659566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(ropts_pargh, PETSC_MAX_PATH_LEN, "PARALLEL_GHOSTS=%" PetscInt_FMT ".0.%" PetscInt_FMT ";", dim, nghost)); 76649d66b22SVijay Mahadevan } 7672e4e7c01SVijay Mahadevan } 7682e4e7c01SVijay Mahadevan 769c8d5365dSVijay Mahadevan if (dbglevel) { 770f90c3b0eSVijay Mahadevan if (numproc > 1) { 7719566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%" PetscInt_FMT ";DEBUG_PIO=%" PetscInt_FMT ";", dbglevel, dbglevel)); 77261eb6e27SVijay Mahadevan } 7739566063dSJacob Faibussowitsch else PetscCall(PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%" PetscInt_FMT ";", dbglevel)); 774c8d5365dSVijay Mahadevan } 77551d15aeeSVijay Mahadevan 7769566063dSJacob Faibussowitsch 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 : ""))); 777f90c3b0eSVijay Mahadevan *read_opts = ropts; 77851d15aeeSVijay Mahadevan PetscFunctionReturn(0); 77951d15aeeSVijay Mahadevan } 78051d15aeeSVijay Mahadevan 781cab5ea25SPierre Jolivet /*@C 782aa859218SVijay Mahadevan DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file. 783aa859218SVijay Mahadevan 784d083f849SBarry Smith Collective 785aa859218SVijay Mahadevan 786aa859218SVijay Mahadevan Input Parameters: 787aa859218SVijay Mahadevan + comm - The communicator for the DM object 788aa859218SVijay Mahadevan . dim - The spatial dimension 789b8ecf6d3SVijay Mahadevan . filename - The name of the mesh file to be loaded 790a2b725a8SWilliam Gropp - usrreadopts - The options string to read a MOAB mesh. 791a2b725a8SWilliam Gropp 792a8d69d7bSBarry Smith Reference (Parallel Mesh Initialization: https://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo) 793aa859218SVijay Mahadevan 794aa859218SVijay Mahadevan Output Parameter: 795aa859218SVijay Mahadevan . dm - The DM object 796aa859218SVijay Mahadevan 797aa859218SVijay Mahadevan Level: beginner 798aa859218SVijay Mahadevan 799db781477SPatrick Sanan .seealso: `DMSetType()`, `DMCreate()`, `DMMoabCreateBoxMesh()` 800aa859218SVijay Mahadevan @*/ 80149d66b22SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm, PetscInt dim, PetscInt nghost, const char* filename, const char* usrreadopts, DM *dm) 80251d15aeeSVijay Mahadevan { 803a4717116SVijay Mahadevan moab::ErrorCode merr; 8042e4e7c01SVijay Mahadevan PetscInt nprocs; 805a4717116SVijay Mahadevan DM_Moab *dmmoab; 806a4717116SVijay Mahadevan moab::Interface *mbiface; 8079daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 808a4717116SVijay Mahadevan moab::ParallelComm *pcomm; 8099daf19fdSVijay Mahadevan #endif 810a4717116SVijay Mahadevan moab::Range verts, elems; 8117023aa44SVijay Mahadevan const char *readopts; 81251d15aeeSVijay Mahadevan 81351d15aeeSVijay Mahadevan PetscFunctionBegin; 81449d66b22SVijay Mahadevan PetscValidPointer(dm, 6); 81551d15aeeSVijay Mahadevan 816a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 8179566063dSJacob Faibussowitsch PetscCall(DMMoabCreateMoab(comm, NULL, NULL, NULL, dm)); 81851d15aeeSVijay Mahadevan 819a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 820a4717116SVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 821a4717116SVijay Mahadevan mbiface = dmmoab->mbiface; 8229daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 823a4717116SVijay Mahadevan pcomm = dmmoab->pcomm; 824a4717116SVijay Mahadevan nprocs = pcomm->size(); 8259daf19fdSVijay Mahadevan #else 8269daf19fdSVijay Mahadevan nprocs = 1; 8279daf19fdSVijay Mahadevan #endif 828aa859218SVijay Mahadevan /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */ 829c8d5365dSVijay Mahadevan dmmoab->dim = dim; 83049d66b22SVijay Mahadevan dmmoab->nghostrings = nghost; 831304006b3SVijay Mahadevan dmmoab->refct = 1; 832a4717116SVijay Mahadevan 833b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */ 834b5410836SVijay Mahadevan merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr); 835b5410836SVijay Mahadevan 836a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */ 8379566063dSJacob Faibussowitsch PetscCall(DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, nghost, dmmoab->read_mode, 8385f80ce2aSJacob Faibussowitsch dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts)); 8397023aa44SVijay Mahadevan 8409566063dSJacob Faibussowitsch PetscCall(PetscInfo(*dm, "Reading file %s with options: %s\n", filename, readopts)); 841a4717116SVijay Mahadevan 842a4717116SVijay Mahadevan /* Load the mesh from a file. */ 84349d66b22SVijay Mahadevan if (dmmoab->fileset) { 8447023aa44SVijay Mahadevan merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface, "Reading MOAB file failed.", merr); 84549d66b22SVijay Mahadevan } 84649d66b22SVijay Mahadevan else { 84749d66b22SVijay Mahadevan merr = mbiface->load_file(filename, 0, readopts);MBERRVM(mbiface, "Reading MOAB file failed.", merr); 84849d66b22SVijay Mahadevan } 849a4717116SVijay Mahadevan 8509daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 8516e40195eSVijay Mahadevan /* Reassign global IDs on all entities. */ 8526ea892caSVijay Mahadevan /* merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, true, true, true);MBERRNM(merr); */ 8539daf19fdSVijay Mahadevan #endif 854e5136372SVijay Mahadevan 855e5136372SVijay Mahadevan /* load the local vertices */ 856e5136372SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr); 857e5136372SVijay Mahadevan /* load the local elements */ 858e5136372SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr); 859e5136372SVijay Mahadevan 8609daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 861e5136372SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 862e5136372SVijay Mahadevan on all of the ghost vertexes */ 863e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag, verts);MBERRV(mbiface, merr); 864e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag, elems);MBERRV(mbiface, merr); 865a4717116SVijay Mahadevan merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr); 8669daf19fdSVijay Mahadevan #endif 867a4717116SVijay Mahadevan 86863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(*dm, "MOAB file '%s' was successfully loaded. Found %zu vertices and %zu elements.\n", filename, verts.size(), elems.size())); 8699566063dSJacob Faibussowitsch PetscCall(PetscFree(readopts)); 87051d15aeeSVijay Mahadevan PetscFunctionReturn(0); 87151d15aeeSVijay Mahadevan } 87251d15aeeSVijay Mahadevan 873cab5ea25SPierre Jolivet /*@C 874304006b3SVijay Mahadevan DMMoabRenumberMeshEntities - Order and number all entities (vertices->elements) to be contiguously ordered 875304006b3SVijay Mahadevan in parallel 876304006b3SVijay Mahadevan 877d083f849SBarry Smith Collective 878304006b3SVijay Mahadevan 879304006b3SVijay Mahadevan Input Parameters: 880a2b725a8SWilliam Gropp . dm - The DM object 881304006b3SVijay Mahadevan 882304006b3SVijay Mahadevan Level: advanced 883304006b3SVijay Mahadevan 884db781477SPatrick Sanan .seealso: `DMSetUp()`, `DMCreate()` 885304006b3SVijay Mahadevan @*/ 886304006b3SVijay Mahadevan PetscErrorCode DMMoabRenumberMeshEntities(DM dm) 887304006b3SVijay Mahadevan { 888304006b3SVijay Mahadevan moab::Range verts; 889304006b3SVijay Mahadevan 890304006b3SVijay Mahadevan PetscFunctionBegin; 891304006b3SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 892304006b3SVijay Mahadevan 8939daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 894304006b3SVijay Mahadevan /* Insert new points */ 8959daf19fdSVijay Mahadevan moab::ErrorCode merr; 8969daf19fdSVijay Mahadevan merr = ((DM_Moab*) dm->data)->pcomm->assign_global_ids(((DM_Moab*) dm->data)->fileset, 3, 0, false, true, false);MBERRNM(merr); 8979daf19fdSVijay Mahadevan #endif 898304006b3SVijay Mahadevan PetscFunctionReturn(0); 899304006b3SVijay Mahadevan } 900