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; 25755f3dfbSVijay Mahadevan 26ce27a4eeSVijay Mahadevan } DMMoabMeshGeneratorCtx; 27ce27a4eeSVijay Mahadevan 28ce27a4eeSVijay Mahadevan 29ce27a4eeSVijay Mahadevan #undef __FUNCT__ 30ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_SetTensorElementConnectivity_Private" 3149d66b22SVijay Mahadevan PetscInt DMMoab_SetTensorElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt offset, PetscInt corner, std::vector<PetscInt>& subent_conn, moab::EntityHandle *connectivity) 3251d15aeeSVijay Mahadevan { 33ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 34e5136372SVijay Mahadevan case 1: 35ce27a4eeSVijay Mahadevan subent_conn.resize(2); 36ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data()); 37ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[0]] = corner; 38ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[1]] = corner + 1; 39e5136372SVijay Mahadevan break; 40e5136372SVijay Mahadevan case 2: 41ce27a4eeSVijay Mahadevan subent_conn.resize(4); 42ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBQUAD, 2, 0, subent_conn.data()); 43ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[0]] = corner; 44ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[1]] = corner + 1; 45ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[2]] = corner + 1 + genCtx.ystride; 46ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[3]] = corner + genCtx.ystride; 47e5136372SVijay Mahadevan break; 48e5136372SVijay Mahadevan case 3: 49e5136372SVijay Mahadevan default: 50ce27a4eeSVijay Mahadevan subent_conn.resize(8); 51ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBHEX, 3, 0, subent_conn.data()); 52ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[0]] = corner; 53ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[1]] = corner + 1; 54ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[2]] = corner + 1 + genCtx.ystride; 55ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[3]] = corner + genCtx.ystride; 56ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[4]] = corner + genCtx.zstride; 57ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[5]] = corner + 1 + genCtx.zstride; 58ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[6]] = corner + 1 + genCtx.ystride + genCtx.zstride; 59ce27a4eeSVijay Mahadevan connectivity[offset + subent_conn[7]] = corner + genCtx.ystride + genCtx.zstride; 60e5136372SVijay Mahadevan break; 61e5136372SVijay Mahadevan } 62ce27a4eeSVijay Mahadevan return subent_conn.size(); 63e5136372SVijay Mahadevan } 6451d15aeeSVijay Mahadevan 65ce27a4eeSVijay Mahadevan 66ce27a4eeSVijay Mahadevan #undef __FUNCT__ 67ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_SetSimplexElementConnectivity_Private" 6849d66b22SVijay Mahadevan PetscInt DMMoab_SetSimplexElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt subelem, PetscInt offset, PetscInt corner, std::vector<PetscInt>& subent_conn, moab::EntityHandle *connectivity) 69c3dd00c7SVijay Mahadevan { 70755f3dfbSVijay Mahadevan PetscInt A, B, C, D, E, F, G, H, M; 71755f3dfbSVijay Mahadevan const PetscInt trigen_opts = 1; /* 1 - Aligned diagonally to right, 2 - Aligned diagonally to left, 3 - 4 elements per quad */ 7249d66b22SVijay Mahadevan A = corner; 7349d66b22SVijay Mahadevan B = corner + 1; 74ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 75c3dd00c7SVijay Mahadevan case 1: 76ce27a4eeSVijay Mahadevan subent_conn.resize(2); /* only linear EDGE supported now */ 77ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data()); 7849d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 7949d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = B; 80c3dd00c7SVijay Mahadevan break; 81c3dd00c7SVijay Mahadevan case 2: 8249d66b22SVijay Mahadevan C = corner + 1 + genCtx.ystride; 8349d66b22SVijay Mahadevan D = corner + genCtx.ystride; 84*63d025dbSVijay Mahadevan M = corner + 0.5; /* technically -- need to modify vertex generation */ 85ce27a4eeSVijay Mahadevan subent_conn.resize(3); /* only linear TRI supported */ 86ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBTRI, 2, 0, subent_conn.data()); 87755f3dfbSVijay Mahadevan if (trigen_opts == 1) { 88ce27a4eeSVijay Mahadevan if (subelem) { /* 0 1 2 of a QUAD */ 8977a8c067SVijay Mahadevan connectivity[offset + subent_conn[0]] = B; 90755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[1]] = C; 9177a8c067SVijay Mahadevan connectivity[offset + subent_conn[2]] = A; 92c3dd00c7SVijay Mahadevan } 93ce27a4eeSVijay Mahadevan else { /* 2 3 0 of a QUAD */ 9477a8c067SVijay Mahadevan connectivity[offset + subent_conn[0]] = D; 9577a8c067SVijay Mahadevan connectivity[offset + subent_conn[1]] = A; 96755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = C; 97755f3dfbSVijay Mahadevan } 98755f3dfbSVijay Mahadevan } 99755f3dfbSVijay Mahadevan else if (trigen_opts == 2) { 100755f3dfbSVijay Mahadevan if (subelem) { /* 0 1 2 of a QUAD */ 101755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 10277a8c067SVijay Mahadevan connectivity[offset + subent_conn[1]] = B; 10377a8c067SVijay Mahadevan connectivity[offset + subent_conn[2]] = D; 104755f3dfbSVijay Mahadevan } 105755f3dfbSVijay Mahadevan else { /* 2 3 0 of a QUAD */ 10677a8c067SVijay Mahadevan connectivity[offset + subent_conn[0]] = C; 10777a8c067SVijay Mahadevan connectivity[offset + subent_conn[1]] = D; 108755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = B; 109755f3dfbSVijay Mahadevan } 110755f3dfbSVijay Mahadevan } 111755f3dfbSVijay Mahadevan else { 112755f3dfbSVijay Mahadevan switch (subelem) { /* 0 1 2 of a QUAD */ 113755f3dfbSVijay Mahadevan case 0: 114755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 115755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[1]] = B; 116755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = M; 117755f3dfbSVijay Mahadevan break; 118755f3dfbSVijay Mahadevan case 1: 119755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[0]] = B; 120755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[1]] = C; 121755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = M; 122755f3dfbSVijay Mahadevan break; 123755f3dfbSVijay Mahadevan case 2: 12449d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = C; 12549d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = D; 126755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = M; 127755f3dfbSVijay Mahadevan break; 128755f3dfbSVijay Mahadevan case 3: 129755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[0]] = D; 130755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[1]] = A; 131755f3dfbSVijay Mahadevan connectivity[offset + subent_conn[2]] = M; 132755f3dfbSVijay Mahadevan break; 133755f3dfbSVijay Mahadevan } 134c3dd00c7SVijay Mahadevan } 135c3dd00c7SVijay Mahadevan break; 136c3dd00c7SVijay Mahadevan case 3: 137c3dd00c7SVijay Mahadevan default: 13849d66b22SVijay Mahadevan C = corner + 1 + genCtx.ystride; 13949d66b22SVijay Mahadevan D = corner + genCtx.ystride; 14049d66b22SVijay Mahadevan E = corner + genCtx.zstride; 14149d66b22SVijay Mahadevan F = corner + 1 + genCtx.zstride; 14249d66b22SVijay Mahadevan G = corner + 1 + genCtx.ystride + genCtx.zstride; 14349d66b22SVijay Mahadevan H = corner + genCtx.ystride + genCtx.zstride; 144ce27a4eeSVijay Mahadevan subent_conn.resize(4); /* only linear TET supported */ 145ce27a4eeSVijay Mahadevan moab::CN::SubEntityVertexIndices(moab::MBTET, 3, 0, subent_conn.data()); 146c3dd00c7SVijay Mahadevan switch (subelem) { 14749d66b22SVijay Mahadevan case 0: /* 4 3 7 6 of a HEX */ 14849d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = E; 14949d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = D; 15049d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = H; 15149d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = G; 152c3dd00c7SVijay Mahadevan break; 15349d66b22SVijay Mahadevan case 1: /* 0 1 2 5 of a HEX */ 15449d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 15549d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = B; 15649d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = C; 15749d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 158c3dd00c7SVijay Mahadevan break; 15949d66b22SVijay Mahadevan case 2: /* 0 3 4 5 of a HEX */ 16049d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 16149d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = D; 16249d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = E; 16349d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 164c3dd00c7SVijay Mahadevan break; 16549d66b22SVijay Mahadevan case 3: /* 2 6 3 5 of a HEX */ 16649d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = C; 16749d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = G; 16849d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = D; 16949d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 170c3dd00c7SVijay Mahadevan break; 17149d66b22SVijay Mahadevan case 4: /* 0 2 3 5 of a HEX */ 17249d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = A; 17349d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = C; 17449d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = D; 17549d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 17649d66b22SVijay Mahadevan break; 17749d66b22SVijay Mahadevan case 5: /* 3 6 4 5 of a HEX */ 17849d66b22SVijay Mahadevan connectivity[offset + subent_conn[0]] = D; 17949d66b22SVijay Mahadevan connectivity[offset + subent_conn[1]] = G; 18049d66b22SVijay Mahadevan connectivity[offset + subent_conn[2]] = E; 18149d66b22SVijay Mahadevan connectivity[offset + subent_conn[3]] = F; 182c3dd00c7SVijay Mahadevan break; 183c3dd00c7SVijay Mahadevan } 184c3dd00c7SVijay Mahadevan break; 185c3dd00c7SVijay Mahadevan } 186ce27a4eeSVijay Mahadevan return subent_conn.size(); 187c3dd00c7SVijay Mahadevan } 188c3dd00c7SVijay Mahadevan 189ce27a4eeSVijay Mahadevan 190ce27a4eeSVijay Mahadevan #undef __FUNCT__ 191ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_SetElementConnectivity_Private" 19249d66b22SVijay Mahadevan std::pair<PetscInt, PetscInt> DMMoab_SetElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt offset, PetscInt corner, moab::EntityHandle *connectivity) 193c3dd00c7SVijay Mahadevan { 19449d66b22SVijay Mahadevan PetscInt vcount = 0; 19549d66b22SVijay Mahadevan PetscInt simplices_per_tensor[4] = {0, 1, 2, 6}; 19649d66b22SVijay Mahadevan std::vector<PetscInt> subent_conn; /* only linear edge, tri, tet supported now */ 197ce27a4eeSVijay Mahadevan subent_conn.reserve(27); 19849d66b22SVijay Mahadevan PetscInt m, subelem; 199ce27a4eeSVijay Mahadevan if (genCtx.simplex) { 20049d66b22SVijay Mahadevan subelem = simplices_per_tensor[genCtx.dim]; 201ce27a4eeSVijay Mahadevan for (m = 0; m < subelem; m++) { 202ce27a4eeSVijay Mahadevan vcount = DMMoab_SetSimplexElementConnectivity_Private(genCtx, m, offset, corner, subent_conn, connectivity); 203ce27a4eeSVijay Mahadevan offset += vcount; 204ce27a4eeSVijay Mahadevan } 205c3dd00c7SVijay Mahadevan } 206c3dd00c7SVijay Mahadevan else { 207c3dd00c7SVijay Mahadevan subelem = 1; 208ce27a4eeSVijay Mahadevan vcount = DMMoab_SetTensorElementConnectivity_Private(genCtx, offset, corner, subent_conn, connectivity); 209c3dd00c7SVijay Mahadevan } 21049d66b22SVijay Mahadevan return std::pair<PetscInt, PetscInt>(vcount * subelem, subelem); 211c3dd00c7SVijay Mahadevan } 212c3dd00c7SVijay Mahadevan 213c3dd00c7SVijay Mahadevan 214ce27a4eeSVijay Mahadevan #undef __FUNCT__ 215ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_GenerateVertices_Private" 21649d66b22SVijay Mahadevan PetscErrorCode DMMoab_GenerateVertices_Private(moab::Interface *mbImpl, moab::ReadUtilIface *iface, DMMoabMeshGeneratorCtx& genCtx, PetscInt m, PetscInt n, PetscInt k, 21749d66b22SVijay Mahadevan PetscInt a, PetscInt b, PetscInt c, moab::Tag& global_id_tag, moab::EntityHandle& startv, moab::Range& uverts) 218ce27a4eeSVijay Mahadevan { 21949d66b22SVijay Mahadevan PetscInt x, y, z, ix, nnodes; 220755f3dfbSVijay Mahadevan PetscErrorCode ierr; 22149d66b22SVijay Mahadevan PetscInt ii, jj, kk; 22249d66b22SVijay Mahadevan std::vector<PetscReal*> arrays; 223755f3dfbSVijay Mahadevan PetscInt* gids; 224ce27a4eeSVijay Mahadevan moab::ErrorCode merr; 225ce27a4eeSVijay Mahadevan 226ce27a4eeSVijay Mahadevan PetscFunctionBegin; 227ce27a4eeSVijay Mahadevan /* we will generate (q*block+1)^3 vertices, and block^3 hexas; q is 1 for linear, 2 for quadratic 228ce27a4eeSVijay Mahadevan * the global id of the vertices will come from m, n, k, a, b, c 229ce27a4eeSVijay Mahadevan * x will vary from m*A*q*block + a*q*block to m*A*q*block+(a+1)*q*block etc. 230ce27a4eeSVijay Mahadevan */ 231ce27a4eeSVijay Mahadevan nnodes = genCtx.blockSizeVertexXYZ[0] * (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] * (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1) : 1); 232755f3dfbSVijay Mahadevan ierr = PetscMalloc1(nnodes, &gids);CHKERRQ(ierr); 233ce27a4eeSVijay Mahadevan 234ce27a4eeSVijay Mahadevan merr = iface->get_node_coords(3, nnodes, 0, startv, arrays); MBERR("Can't get node coords.", merr); 235ce27a4eeSVijay Mahadevan 236ce27a4eeSVijay Mahadevan /* will start with the lower corner: */ 237*63d025dbSVijay Mahadevan /* x = ( m * genCtx.A + a ) * genCtx.q * genCtx.blockSizeElementXYZ[0]; */ 238*63d025dbSVijay Mahadevan /* y = ( n * genCtx.B + b ) * genCtx.q * genCtx.blockSizeElementXYZ[1]; */ 239*63d025dbSVijay Mahadevan /* z = ( k * genCtx.C + c ) * genCtx.q * genCtx.blockSizeElementXYZ[2]; */ 240755f3dfbSVijay Mahadevan 241755f3dfbSVijay Mahadevan x = ( m * genCtx.A + a ) * genCtx.q; 242755f3dfbSVijay Mahadevan y = ( n * genCtx.B + b ) * genCtx.q; 243755f3dfbSVijay Mahadevan z = ( k * genCtx.C + c ) * genCtx.q; 244755f3dfbSVijay Mahadevan PetscInfo3(NULL, "Starting offset for coordinates := %d, %d, %d\n", x, y, z); 245ce27a4eeSVijay Mahadevan ix = 0; 246ce27a4eeSVijay Mahadevan moab::Range verts(startv, startv + nnodes - 1); 247ce27a4eeSVijay Mahadevan for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1); kk++) { 248ce27a4eeSVijay Mahadevan for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] : 1); jj++) { 249ce27a4eeSVijay Mahadevan for (ii = 0; ii < genCtx.blockSizeVertexXYZ[0]; ii++, ix++) { 250ce27a4eeSVijay Mahadevan /* set coordinates for the vertices */ 251ce27a4eeSVijay Mahadevan arrays[0][ix] = (x + ii) * genCtx.dx + genCtx.xyzbounds[0]; 252ce27a4eeSVijay Mahadevan arrays[1][ix] = (y + jj) * genCtx.dy + genCtx.xyzbounds[2]; 253ce27a4eeSVijay Mahadevan arrays[2][ix] = (z + kk) * genCtx.dz + genCtx.xyzbounds[4]; 254755f3dfbSVijay Mahadevan PetscInfo3(NULL, "Creating vertex with coordinates := %f, %f, %f\n", arrays[0][ix], arrays[1][ix], arrays[2][ix]); 255ce27a4eeSVijay Mahadevan 256ce27a4eeSVijay Mahadevan /* If we want to set some tags on the vertices -> use the following entity handle definition: 257ce27a4eeSVijay Mahadevan moab::EntityHandle v = startv + ix; 258ce27a4eeSVijay Mahadevan */ 259ce27a4eeSVijay Mahadevan /* compute the global ID for vertex */ 260ce27a4eeSVijay Mahadevan gids[ix] = 1 + (x + ii) + (y + jj) * genCtx.NX + (z + kk) * (genCtx.NX * genCtx.NY); 261ce27a4eeSVijay Mahadevan } 262ce27a4eeSVijay Mahadevan } 263ce27a4eeSVijay Mahadevan } 264ce27a4eeSVijay Mahadevan /* set global ID data on vertices */ 265ce27a4eeSVijay Mahadevan mbImpl->tag_set_data(global_id_tag, verts, &gids[0]); 266ce27a4eeSVijay Mahadevan verts.swap(uverts); 267755f3dfbSVijay Mahadevan ierr = PetscFree(gids);CHKERRQ(ierr); 268ce27a4eeSVijay Mahadevan PetscFunctionReturn(0); 269ce27a4eeSVijay Mahadevan } 270ce27a4eeSVijay Mahadevan 271ce27a4eeSVijay Mahadevan #undef __FUNCT__ 272ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_GenerateElements_Private" 27349d66b22SVijay Mahadevan PetscErrorCode DMMoab_GenerateElements_Private(moab::Interface* mbImpl, moab::ReadUtilIface* iface, DMMoabMeshGeneratorCtx& genCtx, PetscInt m, PetscInt n, PetscInt k, 27449d66b22SVijay Mahadevan PetscInt a, PetscInt b, PetscInt c, moab::Tag& global_id_tag, moab::EntityHandle startv, moab::Range& cells) 275ce27a4eeSVijay Mahadevan { 276ce27a4eeSVijay Mahadevan moab::ErrorCode merr; 277ce27a4eeSVijay Mahadevan PetscInt ix, ie, xe, ye, ze; 278ce27a4eeSVijay Mahadevan PetscInt ii, jj, kk, nvperelem; 27949d66b22SVijay Mahadevan PetscInt simplices_per_tensor[4] = {0, 1, 2, 6}; 280ce27a4eeSVijay 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);*/ 281ce27a4eeSVijay Mahadevan PetscInt nelems = ntensorelems; 282*63d025dbSVijay Mahadevan moab::EntityHandle starte; /* connectivity */ 283ce27a4eeSVijay Mahadevan moab::EntityHandle* conn; 284ce27a4eeSVijay Mahadevan 285ce27a4eeSVijay Mahadevan PetscFunctionBegin; 286ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 287ce27a4eeSVijay Mahadevan case 1: 288ce27a4eeSVijay Mahadevan nvperelem = 2; 289ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 2, moab::MBEDGE, 0, starte, conn); MBERR("Can't get EDGE2 element connectivity.", merr); 290ce27a4eeSVijay Mahadevan break; 291ce27a4eeSVijay Mahadevan case 2: 292ce27a4eeSVijay Mahadevan if (genCtx.simplex) { 293ce27a4eeSVijay Mahadevan nvperelem = 3; 29449d66b22SVijay Mahadevan nelems = ntensorelems * simplices_per_tensor[genCtx.dim]; 295ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 3, moab::MBTRI, 0, starte, conn); MBERR("Can't get TRI3 element connectivity.", merr); 296ce27a4eeSVijay Mahadevan } 297ce27a4eeSVijay Mahadevan else { 298ce27a4eeSVijay Mahadevan nvperelem = 4; 299ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 4, moab::MBQUAD, 0, starte, conn); MBERR("Can't get QUAD4 element connectivity.", merr); 300ce27a4eeSVijay Mahadevan } 301ce27a4eeSVijay Mahadevan break; 302ce27a4eeSVijay Mahadevan case 3: 303ce27a4eeSVijay Mahadevan default: 304ce27a4eeSVijay Mahadevan if (genCtx.simplex) { 305ce27a4eeSVijay Mahadevan nvperelem = 4; 30649d66b22SVijay Mahadevan nelems = ntensorelems * simplices_per_tensor[genCtx.dim]; 307ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 4, moab::MBTET, 0, starte, conn); MBERR("Can't get TET4 element connectivity.", merr); 308ce27a4eeSVijay Mahadevan } 309ce27a4eeSVijay Mahadevan else { 310ce27a4eeSVijay Mahadevan nvperelem = 8; 311ce27a4eeSVijay Mahadevan merr = iface->get_element_connect(nelems, 8, moab::MBHEX, 0, starte, conn); MBERR("Can't get HEX8 element connectivity.", merr); 312ce27a4eeSVijay Mahadevan } 313ce27a4eeSVijay Mahadevan break; 314ce27a4eeSVijay Mahadevan } 315ce27a4eeSVijay Mahadevan 316*63d025dbSVijay Mahadevan ix = ie = 0; /* index now in the elements, for global ids */ 317ce27a4eeSVijay Mahadevan 318ce27a4eeSVijay Mahadevan /* create a temporary range to store local element handles */ 319ce27a4eeSVijay Mahadevan moab::Range tmp(starte, starte + nelems - 1); 32049d66b22SVijay Mahadevan std::vector<PetscInt> gids(nelems); 321ce27a4eeSVijay Mahadevan 322ce27a4eeSVijay Mahadevan /* identify the elements at the lower corner, for their global ids */ 323ce27a4eeSVijay Mahadevan xe = m * genCtx.A * genCtx.blockSizeElementXYZ[0] + a * genCtx.blockSizeElementXYZ[0]; 324ce27a4eeSVijay Mahadevan ye = (genCtx.dim > 1 ? n * genCtx.B * genCtx.blockSizeElementXYZ[1] + b * genCtx.blockSizeElementXYZ[1] : 0); 325ce27a4eeSVijay Mahadevan ze = (genCtx.dim > 2 ? k * genCtx.C * genCtx.blockSizeElementXYZ[2] + c * genCtx.blockSizeElementXYZ[2] : 0); 326ce27a4eeSVijay Mahadevan 327ce27a4eeSVijay Mahadevan /* create owned elements requested by genCtx */ 328ce27a4eeSVijay Mahadevan for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeElementXYZ[2] : 1); kk++) { 329ce27a4eeSVijay Mahadevan for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeElementXYZ[1] : 1); jj++) { 330ce27a4eeSVijay Mahadevan for (ii = 0; ii < genCtx.blockSizeElementXYZ[0]; ii++) { 331ce27a4eeSVijay Mahadevan 332ce27a4eeSVijay Mahadevan moab::EntityHandle corner = startv + genCtx.q * ii + genCtx.q * jj * genCtx.ystride + genCtx.q * kk * genCtx.zstride; 333ce27a4eeSVijay Mahadevan 33449d66b22SVijay Mahadevan std::pair<PetscInt, PetscInt> entoffset = DMMoab_SetElementConnectivity_Private(genCtx, ix, corner, conn); 335ce27a4eeSVijay Mahadevan 33649d66b22SVijay Mahadevan for (PetscInt j = 0; j < entoffset.second; j++) { 337ce27a4eeSVijay Mahadevan /* The entity handle for the particular element -> if we want to set some tags is 338ce27a4eeSVijay Mahadevan moab::EntityHandle eh = starte + ie + j; 339ce27a4eeSVijay Mahadevan */ 340ce27a4eeSVijay Mahadevan gids[ie + j] = 1 + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)); 341*63d025dbSVijay Mahadevan /* gids[ie+j] = ie + j + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)); */ 342*63d025dbSVijay Mahadevan /* gids[ie+j] = 1 + ie; */ 343*63d025dbSVijay Mahadevan /* ie++; */ 344ce27a4eeSVijay Mahadevan } 345ce27a4eeSVijay Mahadevan 346ce27a4eeSVijay Mahadevan ix += entoffset.first; 34749d66b22SVijay Mahadevan ie += entoffset.second; 348ce27a4eeSVijay Mahadevan } 349ce27a4eeSVijay Mahadevan } 350ce27a4eeSVijay Mahadevan } 351ce27a4eeSVijay Mahadevan if (genCtx.adjEnts) { /* we need to update adjacencies now, because some elements are new */ 352ce27a4eeSVijay Mahadevan merr = iface->update_adjacencies(starte, nelems, nvperelem, conn); MBERR("Can't update adjacencies", merr); 353ce27a4eeSVijay Mahadevan } 354ce27a4eeSVijay Mahadevan tmp.swap(cells); 355ce27a4eeSVijay Mahadevan merr = mbImpl->tag_set_data(global_id_tag, cells, &gids[0]); MBERR("Can't set global ids to elements.", merr); 356ce27a4eeSVijay Mahadevan PetscFunctionReturn(0); 357ce27a4eeSVijay Mahadevan } 358ce27a4eeSVijay Mahadevan 359ce27a4eeSVijay Mahadevan #undef __FUNCT__ 360ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMBUtil_InitializeOptions" 361755f3dfbSVijay Mahadevan PetscErrorCode DMMBUtil_InitializeOptions(DMMoabMeshGeneratorCtx& genCtx, PetscInt dim, PetscBool simplex, PetscInt rank, PetscInt nprocs, const PetscReal* bounds, PetscInt nelems) 362ce27a4eeSVijay Mahadevan { 363ce27a4eeSVijay Mahadevan PetscFunctionBegin; 364ce27a4eeSVijay Mahadevan /* Initialize all genCtx data */ 365ce27a4eeSVijay Mahadevan genCtx.dim = dim; 366ce27a4eeSVijay Mahadevan genCtx.simplex = simplex; 36749d66b22SVijay Mahadevan genCtx.newMergeMethod = genCtx.keep_skins = genCtx.adjEnts = true; 368ce27a4eeSVijay Mahadevan /* determine other global quantities for the mesh used for nodes increments */ 369ce27a4eeSVijay Mahadevan genCtx.q = 1; 370755f3dfbSVijay Mahadevan genCtx.fraction = genCtx.remainder = genCtx.cumfraction = 0; 371ce27a4eeSVijay Mahadevan 372ce27a4eeSVijay Mahadevan if (!genCtx.usrxyzgrid) { /* not overridden by genCtx - assume nele equally and that genCtx wants a uniform cube mesh */ 373755f3dfbSVijay Mahadevan 374755f3dfbSVijay Mahadevan genCtx.fraction = nelems / nprocs; /* partition only by the largest dimension */ 375755f3dfbSVijay Mahadevan genCtx.remainder = nelems % nprocs; /* remainder after partition which gets evenly distributed by round-robin */ 376755f3dfbSVijay Mahadevan genCtx.cumfraction = (rank > 0 ? (genCtx.fraction) * (rank) + (rank - 1 < genCtx.remainder ? rank : genCtx.remainder ) : 0); 377755f3dfbSVijay Mahadevan if (rank < genCtx.remainder) /* This process gets "fraction+1" elements */ 378755f3dfbSVijay Mahadevan genCtx.fraction++; 379755f3dfbSVijay Mahadevan 3809daf19fdSVijay Mahadevan PetscInfo3(NULL, "Fraction = %D, Remainder = %D, Cumulative fraction = %D\n", genCtx.fraction, genCtx.remainder, genCtx.cumfraction); 381755f3dfbSVijay Mahadevan switch (genCtx.dim) { 382755f3dfbSVijay Mahadevan case 1: 383755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = genCtx.fraction; 384755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = 1; 385755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = 1; 386755f3dfbSVijay Mahadevan break; 387755f3dfbSVijay Mahadevan case 2: 388755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = nelems; 389755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = genCtx.fraction; 390755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = 1; 391755f3dfbSVijay Mahadevan break; 392755f3dfbSVijay Mahadevan case 3: 393755f3dfbSVijay Mahadevan default: 394755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = nelems; 395755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = nelems; 396755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = genCtx.fraction; 397755f3dfbSVijay Mahadevan break; 39849d66b22SVijay Mahadevan } 399ce27a4eeSVijay Mahadevan } 400ce27a4eeSVijay Mahadevan 401ce27a4eeSVijay Mahadevan /* partition only by the largest dimension */ 402ce27a4eeSVijay Mahadevan /* Total number of local elements := genCtx.blockSizeElementXYZ[0]*(genCtx.dim>1? genCtx.blockSizeElementXYZ[1]*(genCtx.dim>2 ? genCtx.blockSizeElementXYZ[2]:1) :1); */ 403ce27a4eeSVijay Mahadevan if (bounds) { 40449d66b22SVijay Mahadevan for (PetscInt i = 0; i < 6; i++) 405ce27a4eeSVijay Mahadevan genCtx.xyzbounds[i] = bounds[i]; 406ce27a4eeSVijay Mahadevan } 407ce27a4eeSVijay Mahadevan else { 408ce27a4eeSVijay Mahadevan genCtx.xyzbounds[0] = genCtx.xyzbounds[2] = genCtx.xyzbounds[4] = 0.0; 409ce27a4eeSVijay Mahadevan genCtx.xyzbounds[1] = genCtx.xyzbounds[3] = genCtx.xyzbounds[5] = 1.0; 410ce27a4eeSVijay Mahadevan } 411ce27a4eeSVijay Mahadevan 412ce27a4eeSVijay Mahadevan if (!genCtx.usrprocgrid) { 413ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 414ce27a4eeSVijay Mahadevan case 1: 415ce27a4eeSVijay Mahadevan genCtx.M = nprocs; 416ce27a4eeSVijay Mahadevan genCtx.N = genCtx.K = 1; 417ce27a4eeSVijay Mahadevan break; 418ce27a4eeSVijay Mahadevan case 2: 419755f3dfbSVijay Mahadevan genCtx.N = nprocs; 420*63d025dbSVijay Mahadevan genCtx.M = genCtx.K = 1; 421ce27a4eeSVijay Mahadevan break; 422ce27a4eeSVijay Mahadevan default: 423755f3dfbSVijay Mahadevan genCtx.K = nprocs; 424*63d025dbSVijay Mahadevan genCtx.M = genCtx.N = 1; 425ce27a4eeSVijay Mahadevan break; 426ce27a4eeSVijay Mahadevan } 427ce27a4eeSVijay Mahadevan } 428ce27a4eeSVijay Mahadevan 429ce27a4eeSVijay Mahadevan if (!genCtx.usrrefgrid) { 430ce27a4eeSVijay Mahadevan genCtx.A = genCtx.B = genCtx.C = 1; 431ce27a4eeSVijay Mahadevan } 432ce27a4eeSVijay Mahadevan 433ce27a4eeSVijay Mahadevan /* more default values */ 434ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.ney = genCtx.nez = 0; 435ce27a4eeSVijay Mahadevan genCtx.xstride = genCtx.ystride = genCtx.zstride = 0; 436ce27a4eeSVijay Mahadevan genCtx.NX = genCtx.NY = genCtx.NZ = 0; 437ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.ney = genCtx.nez = 0; 438ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 1; 439ce27a4eeSVijay Mahadevan 440ce27a4eeSVijay Mahadevan switch (genCtx.dim) { 441ce27a4eeSVijay Mahadevan case 3: 442ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 443ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1; 444ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[2] = genCtx.q * genCtx.blockSizeElementXYZ[2] + 1; 445ce27a4eeSVijay Mahadevan 446ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 447755f3dfbSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */ 448ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1); 449ce27a4eeSVijay Mahadevan genCtx.xstride = 1; 450ce27a4eeSVijay Mahadevan genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1]; /* number of elements in y direction .... */ 451755f3dfbSVijay Mahadevan genCtx.dy = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */ 452ce27a4eeSVijay Mahadevan genCtx.NY = (genCtx.q * genCtx.ney + 1); 453ce27a4eeSVijay Mahadevan genCtx.ystride = genCtx.blockSizeVertexXYZ[0]; 454ce27a4eeSVijay Mahadevan genCtx.nez = genCtx.K * genCtx.C * genCtx.blockSizeElementXYZ[2]; /* number of elements in z direction .... */ 455755f3dfbSVijay Mahadevan genCtx.dz = (genCtx.xyzbounds[5] - genCtx.xyzbounds[4]) / (nelems * genCtx.q); /* distance between 2 nodes in z direction */ 456ce27a4eeSVijay Mahadevan genCtx.NZ = (genCtx.q * genCtx.nez + 1); 457ce27a4eeSVijay Mahadevan genCtx.zstride = genCtx.blockSizeVertexXYZ[0] * genCtx.blockSizeVertexXYZ[1]; 458ce27a4eeSVijay Mahadevan break; 459ce27a4eeSVijay Mahadevan case 2: 460ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 461ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1; 4629daf19fdSVijay Mahadevan genCtx.blockSizeVertexXYZ[2] = 0; 463ce27a4eeSVijay Mahadevan 464ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 465ce27a4eeSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (genCtx.nex * genCtx.q); /* distance between 2 nodes in x direction */ 466ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1); 467ce27a4eeSVijay Mahadevan genCtx.xstride = 1; 468ce27a4eeSVijay Mahadevan genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1]; /* number of elements in y direction .... */ 469755f3dfbSVijay Mahadevan genCtx.dy = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */ 470ce27a4eeSVijay Mahadevan genCtx.NY = (genCtx.q * genCtx.ney + 1); 471ce27a4eeSVijay Mahadevan genCtx.ystride = genCtx.blockSizeVertexXYZ[0]; 472ce27a4eeSVijay Mahadevan break; 473ce27a4eeSVijay Mahadevan case 1: 4749daf19fdSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 0; 475ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1; 476ce27a4eeSVijay Mahadevan 477ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */ 478755f3dfbSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */ 479ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1); 480ce27a4eeSVijay Mahadevan genCtx.xstride = 1; 481ce27a4eeSVijay Mahadevan break; 482ce27a4eeSVijay Mahadevan } 483ce27a4eeSVijay Mahadevan 484755f3dfbSVijay Mahadevan /* Lets check for some valid input */ 485755f3dfbSVijay Mahadevan if (genCtx.dim < 1 || genCtx.dim > 3) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid topological dimension specified: %d.\n", genCtx.dim); 486755f3dfbSVijay Mahadevan if (genCtx.M * genCtx.N * genCtx.K != nprocs) SETERRQ4(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid [m, n, k] data: %d, %d, %d. Product must be equal to global size = %d.\n", genCtx.M, genCtx.N, genCtx.K, nprocs); 487755f3dfbSVijay Mahadevan /* validate the bounds data */ 488755f3dfbSVijay Mahadevan if (genCtx.xyzbounds[0] >= genCtx.xyzbounds[1]) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "X-dim: Left boundary cannot be greater than right. [%G >= %G]\n", genCtx.xyzbounds[0], genCtx.xyzbounds[1]); 489755f3dfbSVijay Mahadevan if (genCtx.dim > 1 && (genCtx.xyzbounds[2] >= genCtx.xyzbounds[3])) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Y-dim: Left boundary cannot be greater than right. [%G >= %G]\n", genCtx.xyzbounds[2], genCtx.xyzbounds[3]); 490755f3dfbSVijay Mahadevan if (genCtx.dim > 2 && (genCtx.xyzbounds[4] >= genCtx.xyzbounds[5])) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Z-dim: Left boundary cannot be greater than right. [%G >= %G]\n", genCtx.xyzbounds[4], genCtx.xyzbounds[5]); 491755f3dfbSVijay Mahadevan 492ce27a4eeSVijay Mahadevan PetscInfo3(NULL, "Local elements:= %d, %d, %d\n", genCtx.blockSizeElementXYZ[0], genCtx.blockSizeElementXYZ[1], genCtx.blockSizeElementXYZ[2]); 493ce27a4eeSVijay Mahadevan PetscInfo3(NULL, "Local vertices:= %d, %d, %d\n", genCtx.blockSizeVertexXYZ[0], genCtx.blockSizeVertexXYZ[1], genCtx.blockSizeVertexXYZ[2]); 494755f3dfbSVijay Mahadevan PetscInfo3(NULL, "Local blocks/processors := %d, %d, %d\n", genCtx.A, genCtx.B, genCtx.C); 495755f3dfbSVijay Mahadevan PetscInfo3(NULL, "Local processors := %d, %d, %d\n", genCtx.M, genCtx.N, genCtx.K); 496ce27a4eeSVijay Mahadevan PetscInfo3(NULL, "Local nexyz:= %d, %d, %d\n", genCtx.nex, genCtx.ney, genCtx.nez); 497ce27a4eeSVijay Mahadevan PetscInfo3(NULL, "Local delxyz:= %g, %g, %g\n", genCtx.dx, genCtx.dy, genCtx.dz); 498ce27a4eeSVijay Mahadevan PetscInfo3(NULL, "Local strides:= %d, %d, %d\n", genCtx.xstride, genCtx.ystride, genCtx.zstride); 499ce27a4eeSVijay Mahadevan PetscFunctionReturn(0); 500ce27a4eeSVijay Mahadevan } 501ce27a4eeSVijay Mahadevan 502aa859218SVijay Mahadevan /*@ 503ce27a4eeSVijay Mahadevan DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with genCtx specified bounds. 504aa859218SVijay Mahadevan 505aa859218SVijay Mahadevan Collective on MPI_Comm 506aa859218SVijay Mahadevan 507aa859218SVijay Mahadevan Input Parameters: 508aa859218SVijay Mahadevan + comm - The communicator for the DM object 509aa859218SVijay Mahadevan . dim - The spatial dimension 510b8ecf6d3SVijay 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 511b8ecf6d3SVijay Mahadevan . nele - The number of discrete elements in each direction 512b8ecf6d3SVijay Mahadevan . user_nghost - The number of ghosted layers needed in the partitioned mesh 513aa859218SVijay Mahadevan 514aa859218SVijay Mahadevan Output Parameter: 515aa859218SVijay Mahadevan . dm - The DM object 516aa859218SVijay Mahadevan 517aa859218SVijay Mahadevan Level: beginner 518aa859218SVijay Mahadevan 519aa859218SVijay Mahadevan .keywords: DM, create 520aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile() 521aa859218SVijay Mahadevan @*/ 522ce27a4eeSVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal* bounds, PetscInt nele, PetscInt nghost, DM *dm) 52351d15aeeSVijay Mahadevan { 52451d15aeeSVijay Mahadevan PetscErrorCode ierr; 52551d15aeeSVijay Mahadevan moab::ErrorCode merr; 5269daf19fdSVijay Mahadevan PetscInt a, b, c, n, global_size, global_rank; 52751d15aeeSVijay Mahadevan DM_Moab *dmmoab; 528ce27a4eeSVijay Mahadevan moab::Interface *mbImpl; 5299daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 53051d15aeeSVijay Mahadevan moab::ParallelComm *pcomm; 5319daf19fdSVijay Mahadevan #endif 532a4717116SVijay Mahadevan moab::ReadUtilIface *readMeshIface; 533ce27a4eeSVijay Mahadevan moab::Range verts, cells, edges, faces, adj, dim3, dim2; 534ce27a4eeSVijay Mahadevan DMMoabMeshGeneratorCtx genCtx; 535a4717116SVijay Mahadevan const PetscInt npts = nele + 1; /* Number of points in every dimension */ 536ce27a4eeSVijay Mahadevan 537ce27a4eeSVijay Mahadevan moab::Tag global_id_tag, part_tag, geom_tag; 538ce27a4eeSVijay Mahadevan moab::Range ownedvtx, ownedelms, localvtxs, localelms; 539ce27a4eeSVijay Mahadevan moab::EntityHandle regionset; 540755f3dfbSVijay Mahadevan PetscInt ml = 0, nl = 0, kl = 0; 54151d15aeeSVijay Mahadevan 54251d15aeeSVijay Mahadevan PetscFunctionBegin; 543e5136372SVijay Mahadevan if (dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension argument for mesh: dim=[1,3].\n"); 544e5136372SVijay Mahadevan 545755f3dfbSVijay Mahadevan ierr = PetscLogEventRegister("GenerateMesh", DM_CLASSID, &genCtx.generateMesh);CHKERRQ(ierr); 546755f3dfbSVijay Mahadevan ierr = PetscLogEventRegister("AddVertices", DM_CLASSID, &genCtx.generateVertices);CHKERRQ(ierr); 547755f3dfbSVijay Mahadevan ierr = PetscLogEventRegister("AddElements", DM_CLASSID, &genCtx.generateElements);CHKERRQ(ierr); 548755f3dfbSVijay Mahadevan ierr = PetscLogEventRegister("ParResolve", DM_CLASSID, &genCtx.parResolve);CHKERRQ(ierr); 549755f3dfbSVijay Mahadevan ierr = PetscLogEventBegin(genCtx.generateMesh, 0, 0, 0, 0);CHKERRQ(ierr); 550ce27a4eeSVijay Mahadevan ierr = MPI_Comm_size(comm, &global_size);CHKERRQ(ierr); 551e5136372SVijay Mahadevan /* total number of vertices in all dimensions */ 552e5136372SVijay Mahadevan n = pow(npts, dim); 553e5136372SVijay Mahadevan 554e5136372SVijay Mahadevan /* do some error checking */ 555e5136372SVijay Mahadevan if (n < 2) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be >= 2.\n"); 556ce27a4eeSVijay Mahadevan if (global_size > n) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of processors must be less than or equal to number of elements.\n"); 557e5136372SVijay Mahadevan if (nghost < 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of ghost layers cannot be negative.\n"); 558e5136372SVijay Mahadevan 559a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 5609daf19fdSVijay Mahadevan ierr = DMMoabCreateMoab(comm, NULL, NULL, NULL, dm);CHKERRQ(ierr); 56151d15aeeSVijay Mahadevan 562a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 56351d15aeeSVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 564ce27a4eeSVijay Mahadevan mbImpl = dmmoab->mbiface; 5659daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 56651d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm; 567ce27a4eeSVijay Mahadevan global_rank = pcomm->rank(); 5689daf19fdSVijay Mahadevan #else 5699daf19fdSVijay Mahadevan global_rank = 0; 5709daf19fdSVijay Mahadevan global_size = 1; 5719daf19fdSVijay Mahadevan #endif 5729daf19fdSVijay Mahadevan global_id_tag = dmmoab->ltog_tag; 573c8d5365dSVijay Mahadevan dmmoab->dim = dim; 57449d66b22SVijay Mahadevan dmmoab->nghostrings = nghost; 575304006b3SVijay Mahadevan dmmoab->refct = 1; 57651d15aeeSVijay Mahadevan 577b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */ 578ce27a4eeSVijay Mahadevan merr = mbImpl->create_meshset(moab::MESHSET_SET, dmmoab->fileset); MBERR("Creating file set failed", merr); 579b5410836SVijay Mahadevan 580a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */ 581ce27a4eeSVijay Mahadevan merr = mbImpl->query_interface(readMeshIface); MBERRNM(merr); 58251d15aeeSVijay Mahadevan 583755f3dfbSVijay Mahadevan genCtx.M = genCtx.N = genCtx.K = 1; 584755f3dfbSVijay Mahadevan genCtx.A = genCtx.B = genCtx.C = 1; 585755f3dfbSVijay Mahadevan 586755f3dfbSVijay Mahadevan ierr = PetscOptionsBegin(comm, "", "DMMoab Creation Options", "DMMOAB");CHKERRQ(ierr); 587ce27a4eeSVijay Mahadevan /* Handle DMMoab spatial resolution */ 588ce27a4eeSVijay Mahadevan PetscOptionsInt("-dmb_grid_x", "Number of grid points in x direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[0], &genCtx.blockSizeElementXYZ[0], &genCtx.usrxyzgrid); 589ce27a4eeSVijay Mahadevan if (dim > 1) {PetscOptionsInt("-dmb_grid_y", "Number of grid points in y direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[1], &genCtx.blockSizeElementXYZ[1], &genCtx.usrxyzgrid);} 590ce27a4eeSVijay Mahadevan if (dim > 2) {PetscOptionsInt("-dmb_grid_z", "Number of grid points in z direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[2], &genCtx.blockSizeElementXYZ[2], &genCtx.usrxyzgrid);} 591a4717116SVijay Mahadevan 592ce27a4eeSVijay Mahadevan /* Handle DMMoab parallel distibution */ 593ce27a4eeSVijay Mahadevan PetscOptionsInt("-dmb_processors_x", "Number of processors in x direction", "DMMoabSetNumProcs", genCtx.M, &genCtx.M, &genCtx.usrprocgrid); 594ce27a4eeSVijay Mahadevan if (dim > 1) {PetscOptionsInt("-dmb_processors_y", "Number of processors in y direction", "DMMoabSetNumProcs", genCtx.N, &genCtx.N, &genCtx.usrprocgrid);} 595ce27a4eeSVijay Mahadevan if (dim > 2) {PetscOptionsInt("-dmb_processors_z", "Number of processors in z direction", "DMMoabSetNumProcs", genCtx.K, &genCtx.K, &genCtx.usrprocgrid);} 59651d15aeeSVijay Mahadevan 597ce27a4eeSVijay Mahadevan /* Handle DMMoab block refinement */ 598ce27a4eeSVijay Mahadevan PetscOptionsInt("-dmb_refine_x", "Number of refinement blocks in x direction", "DMMoabSetRefinement", genCtx.A, &genCtx.A, &genCtx.usrrefgrid); 599ce27a4eeSVijay Mahadevan if (dim > 1) {PetscOptionsInt("-dmb_refine_y", "Number of refinement blocks in y direction", "DMMoabSetRefinement", genCtx.B, &genCtx.B, &genCtx.usrrefgrid);} 600ce27a4eeSVijay Mahadevan if (dim > 2) {PetscOptionsInt("-dmb_refine_z", "Number of refinement blocks in z direction", "DMMoabSetRefinement", genCtx.C, &genCtx.C, &genCtx.usrrefgrid);} 601ce27a4eeSVijay Mahadevan PetscOptionsEnd(); 60251d15aeeSVijay Mahadevan 603ce27a4eeSVijay Mahadevan ierr = DMMBUtil_InitializeOptions(genCtx, dim, useSimplex, global_rank, global_size, bounds, nele);CHKERRQ(ierr); 60451d15aeeSVijay Mahadevan 60564e1c140SVijay Mahadevan //if(nele<nprocs) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The dimensional discretization size should be greater or equal to number of processors: %D < %D",nele,nprocs); 60664e1c140SVijay Mahadevan 607ce27a4eeSVijay Mahadevan if (genCtx.adjEnts) genCtx.keep_skins = true; /* do not delete anything - consumes more memory */ 60866f88a78SVijay Mahadevan 609ce27a4eeSVijay Mahadevan /* determine m, n, k for processor rank */ 610755f3dfbSVijay Mahadevan ml = nl = kl = 0; 611755f3dfbSVijay Mahadevan switch (genCtx.dim) { 612755f3dfbSVijay Mahadevan case 1: 613755f3dfbSVijay Mahadevan ml = (genCtx.cumfraction); 614755f3dfbSVijay Mahadevan break; 615755f3dfbSVijay Mahadevan case 2: 616755f3dfbSVijay Mahadevan nl = (genCtx.cumfraction); 617755f3dfbSVijay Mahadevan break; 618755f3dfbSVijay Mahadevan default: 619755f3dfbSVijay Mahadevan kl = (genCtx.cumfraction) / genCtx.q / genCtx.blockSizeElementXYZ[2] / genCtx.C; //genCtx.K 620755f3dfbSVijay Mahadevan break; 621755f3dfbSVijay Mahadevan } 622e5136372SVijay Mahadevan 623ce27a4eeSVijay Mahadevan /* 624ce27a4eeSVijay Mahadevan * so there are a total of M * A * blockSizeElement elements in x direction (so M * A * blockSizeElement + 1 verts in x direction) 625ce27a4eeSVijay Mahadevan * so there are a total of N * B * blockSizeElement elements in y direction (so N * B * blockSizeElement + 1 verts in y direction) 626ce27a4eeSVijay Mahadevan * so there are a total of K * C * blockSizeElement elements in z direction (so K * C * blockSizeElement + 1 verts in z direction) 627c8d5365dSVijay Mahadevan 628ce27a4eeSVijay Mahadevan * there are ( M * A blockSizeElement ) * ( N * B * blockSizeElement) * (K * C * blockSizeElement ) hexas 629ce27a4eeSVijay Mahadevan * there are ( M * A * blockSizeElement + 1) * ( N * B * blockSizeElement + 1 ) * (K * C * blockSizeElement + 1) vertices 630ce27a4eeSVijay Mahadevan * x is the first dimension that varies 631ce27a4eeSVijay Mahadevan */ 632e5136372SVijay Mahadevan 633ce27a4eeSVijay Mahadevan /* generate the block at (a, b, c); it will represent a partition , it will get a partition tag */ 63449d66b22SVijay Mahadevan PetscInt dum_id = -1; 635ce27a4eeSVijay Mahadevan merr = mbImpl->tag_get_handle("GLOBAL_ID", 1, moab::MB_TYPE_INTEGER, global_id_tag); MBERR("Getting Tag handle failed", merr); 63651d15aeeSVijay Mahadevan 637ce27a4eeSVijay 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 Tag handle failed", merr); 63851d15aeeSVijay Mahadevan 6393a4aeca1SVijay Mahadevan /* lets create some sets */ 640ce27a4eeSVijay 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); 64149d66b22SVijay Mahadevan merr = mbImpl->create_meshset(moab::MESHSET_SET, regionset); MBERRNM(merr); 642755f3dfbSVijay Mahadevan ierr = PetscLogEventEnd(genCtx.generateMesh, 0, 0, 0, 0);CHKERRQ(ierr); 6433a4aeca1SVijay Mahadevan 644ce27a4eeSVijay Mahadevan for (a = 0; a < (genCtx.dim > 0 ? genCtx.A : genCtx.A); a++) { 645ce27a4eeSVijay Mahadevan for (b = 0; b < (genCtx.dim > 1 ? genCtx.B : 1); b++) { 646ce27a4eeSVijay Mahadevan for (c = 0; c < (genCtx.dim > 2 ? genCtx.C : 1); c++) { 64749d66b22SVijay Mahadevan 64849d66b22SVijay Mahadevan moab::EntityHandle startv; 64949d66b22SVijay Mahadevan 650755f3dfbSVijay Mahadevan ierr = PetscLogEventBegin(genCtx.generateVertices, 0, 0, 0, 0);CHKERRQ(ierr); 651ce27a4eeSVijay Mahadevan ierr = DMMoab_GenerateVertices_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, verts);CHKERRQ(ierr); 652755f3dfbSVijay Mahadevan ierr = PetscLogEventEnd(genCtx.generateVertices, 0, 0, 0, 0);CHKERRQ(ierr); 6533a4aeca1SVijay Mahadevan 654755f3dfbSVijay Mahadevan ierr = PetscLogEventBegin(genCtx.generateElements, 0, 0, 0, 0);CHKERRQ(ierr); 655ce27a4eeSVijay Mahadevan ierr = DMMoab_GenerateElements_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, cells);CHKERRQ(ierr); 656755f3dfbSVijay Mahadevan ierr = PetscLogEventEnd(genCtx.generateElements, 0, 0, 0, 0);CHKERRQ(ierr); 6573a4aeca1SVijay Mahadevan 65849d66b22SVijay Mahadevan PetscInt part_num = 0; 65949d66b22SVijay Mahadevan switch (genCtx.dim) { 66049d66b22SVijay Mahadevan case 3: 66149d66b22SVijay Mahadevan part_num += (c + kl * genCtx.C) * (genCtx.M * genCtx.A * genCtx.N * genCtx.B); 66249d66b22SVijay Mahadevan case 2: 66349d66b22SVijay Mahadevan part_num += (b + nl * genCtx.B) * (genCtx.M * genCtx.A); 66449d66b22SVijay Mahadevan case 1: 66549d66b22SVijay Mahadevan part_num += (a + ml * genCtx.A); 66649d66b22SVijay Mahadevan break; 66749d66b22SVijay Mahadevan } 66849d66b22SVijay Mahadevan 66949d66b22SVijay Mahadevan moab::EntityHandle part_set; 670ce27a4eeSVijay Mahadevan merr = mbImpl->create_meshset(moab::MESHSET_SET, part_set); MBERR("Can't create mesh set.", merr); 6713a4aeca1SVijay Mahadevan 672ce27a4eeSVijay Mahadevan merr = mbImpl->add_entities(part_set, verts); MBERR("Can't add vertices to set.", merr); 67349d66b22SVijay Mahadevan merr = mbImpl->add_entities(part_set, cells); MBERR("Can't add entities to set.", merr); 67449d66b22SVijay Mahadevan // PetscInfo2(NULL, "Generated local mesh: %D vertices and %D elements.\n", verts.size(), cells.size()); 6753a4aeca1SVijay Mahadevan 676ce27a4eeSVijay Mahadevan merr = mbImpl->add_entities(regionset, cells); MBERR("Can't add entities to set.", merr); 677ce27a4eeSVijay Mahadevan 678ce27a4eeSVijay Mahadevan /* if needed, add all edges and faces */ 679ce27a4eeSVijay Mahadevan if (genCtx.adjEnts) 680ce27a4eeSVijay Mahadevan { 681ce27a4eeSVijay Mahadevan if (genCtx.dim > 1) { 682ce27a4eeSVijay Mahadevan merr = mbImpl->get_adjacencies(cells, 1, true, edges, moab::Interface::UNION); MBERR("Can't get edges", merr); 683ce27a4eeSVijay Mahadevan merr = mbImpl->add_entities(part_set, edges); MBERR("Can't add edges to partition set.", merr); 684ce27a4eeSVijay Mahadevan } 685ce27a4eeSVijay Mahadevan if (genCtx.dim > 2) { 686ce27a4eeSVijay Mahadevan merr = mbImpl->get_adjacencies(cells, 2, true, faces, moab::Interface::UNION); MBERR("Can't get faces", merr); 687ce27a4eeSVijay Mahadevan merr = mbImpl->add_entities(part_set, faces); MBERR("Can't add faces to partition set.", merr); 688ce27a4eeSVijay Mahadevan } 689ce27a4eeSVijay Mahadevan edges.clear(); 690ce27a4eeSVijay Mahadevan faces.clear(); 691ce27a4eeSVijay Mahadevan } 69249d66b22SVijay Mahadevan verts.clear(); cells.clear(); 693ce27a4eeSVijay Mahadevan 694ce27a4eeSVijay Mahadevan merr = mbImpl->tag_set_data(part_tag, &part_set, 1, &part_num); MBERR("Can't set part tag on set", merr); 69549d66b22SVijay Mahadevan if (dmmoab->fileset) { 696ce27a4eeSVijay Mahadevan merr = mbImpl->add_parent_child(dmmoab->fileset, part_set); MBERR("Can't add part set to file set.", merr); 697ce27a4eeSVijay Mahadevan merr = mbImpl->unite_meshset(dmmoab->fileset, part_set); MBERRNM(merr); 69849d66b22SVijay Mahadevan } 69949d66b22SVijay Mahadevan merr = mbImpl->add_entities(dmmoab->fileset, &part_set, 1); MBERRNM(merr); 700ce27a4eeSVijay Mahadevan } 701ce27a4eeSVijay Mahadevan } 702ce27a4eeSVijay Mahadevan } 703ce27a4eeSVijay Mahadevan 70449d66b22SVijay Mahadevan /* set geometric dimension tag for regions */ 70549d66b22SVijay Mahadevan merr = mbImpl->tag_set_data(geom_tag, ®ionset, 1, &dmmoab->dim); MBERRNM(merr); 70649d66b22SVijay Mahadevan merr = mbImpl->add_parent_child(dmmoab->fileset, regionset); MBERRNM(merr); 707ce27a4eeSVijay Mahadevan 70849d66b22SVijay Mahadevan /* Only in parallel: resolve shared entities between processors and exchange ghost layers */ 709ce27a4eeSVijay Mahadevan if (global_size > 1) { 710ce27a4eeSVijay Mahadevan 71177a8c067SVijay Mahadevan ierr = PetscLogEventBegin(genCtx.parResolve, 0, 0, 0, 0);CHKERRQ(ierr); 71277a8c067SVijay Mahadevan 713ce27a4eeSVijay Mahadevan merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, genCtx.dim, cells); MBERR("Can't get all d-dimensional elements.", merr); 714ce27a4eeSVijay Mahadevan merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 0, verts); MBERR("Can't get all vertices.", merr); 715ce27a4eeSVijay Mahadevan 716ce27a4eeSVijay Mahadevan if (genCtx.A * genCtx.B * genCtx.C != 1) { // merge needed 71749d66b22SVijay Mahadevan moab::MergeMesh mm(mbImpl); 718ce27a4eeSVijay Mahadevan if (genCtx.newMergeMethod) { 719ce27a4eeSVijay Mahadevan merr = mm.merge_using_integer_tag(verts, global_id_tag); MBERR("Can't merge with GLOBAL_ID tag", merr); 720ce27a4eeSVijay Mahadevan } 721ce27a4eeSVijay Mahadevan else { 722ce27a4eeSVijay Mahadevan merr = mm.merge_entities(cells, 0.0001); MBERR("Can't merge with coordinates", merr); 7233a4aeca1SVijay Mahadevan } 7243a4aeca1SVijay Mahadevan } 7253a4aeca1SVijay Mahadevan 7269daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 727a4717116SVijay Mahadevan /* check the handles */ 728ce27a4eeSVijay Mahadevan merr = pcomm->check_all_shared_handles(); MBERRV(mbImpl, merr); 72951d15aeeSVijay Mahadevan 730a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */ 731ce27a4eeSVijay Mahadevan merr = pcomm->resolve_shared_ents(dmmoab->fileset, cells, dim, dim - 1, NULL, &global_id_tag); MBERRV(mbImpl, merr); 73249d66b22SVijay Mahadevan if (dmmoab->fileset) { 733ce27a4eeSVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false, &dmmoab->fileset); MBERRV(mbImpl, merr); 73449d66b22SVijay Mahadevan } 73549d66b22SVijay Mahadevan else { 73649d66b22SVijay Mahadevan merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false); MBERRV(mbImpl, merr); 73749d66b22SVijay Mahadevan } 738c8d5365dSVijay Mahadevan 739e427d9c9SVijay Mahadevan /* Reassign global IDs on all entities. */ 740e427d9c9SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, false, true, false); MBERRNM(merr); 7419daf19fdSVijay Mahadevan #endif 74277a8c067SVijay Mahadevan 74377a8c067SVijay Mahadevan ierr = PetscLogEventEnd(genCtx.parResolve, 0, 0, 0, 0);CHKERRQ(ierr); 74449d66b22SVijay Mahadevan } 7453f1c6e43SVijay Mahadevan 746ce27a4eeSVijay Mahadevan if (!genCtx.keep_skins) { // default is to delete the 1- and 2-dimensional entities 747ce27a4eeSVijay Mahadevan // delete all quads and edges 748ce27a4eeSVijay Mahadevan moab::Range toDelete; 749ce27a4eeSVijay Mahadevan if (genCtx.dim > 1) { 750ce27a4eeSVijay Mahadevan merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 1, toDelete); MBERR("Can't get edges", merr); 751ce27a4eeSVijay Mahadevan } 7523f1c6e43SVijay Mahadevan 753ce27a4eeSVijay Mahadevan if (genCtx.dim > 2) { 754ce27a4eeSVijay Mahadevan merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 2, toDelete); MBERR("Can't get faces", merr); 755ce27a4eeSVijay Mahadevan } 756c8d5365dSVijay Mahadevan 7579daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 758ce27a4eeSVijay Mahadevan merr = dmmoab->pcomm->delete_entities(toDelete) ; MBERR("Can't delete entities", merr); 7599daf19fdSVijay Mahadevan #endif 760ce27a4eeSVijay Mahadevan } 76151d15aeeSVijay Mahadevan PetscFunctionReturn(0); 76251d15aeeSVijay Mahadevan } 76351d15aeeSVijay Mahadevan 76451d15aeeSVijay Mahadevan 76549d66b22SVijay Mahadevan #undef __FUNCT__ 76649d66b22SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private" 76749d66b22SVijay 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) 76851d15aeeSVijay Mahadevan { 7696acfe860SVijay Mahadevan char *ropts; 77049d66b22SVijay Mahadevan char ropts_par[PETSC_MAX_PATH_LEN], ropts_pargh[PETSC_MAX_PATH_LEN]; 77161eb6e27SVijay Mahadevan char ropts_dbg[PETSC_MAX_PATH_LEN]; 772f90c3b0eSVijay Mahadevan PetscErrorCode ierr; 77351d15aeeSVijay Mahadevan 774a4717116SVijay Mahadevan PetscFunctionBegin; 77595dccacaSBarry Smith ierr = PetscMalloc1(PETSC_MAX_PATH_LEN, &ropts);CHKERRQ(ierr); 77661eb6e27SVijay Mahadevan ierr = PetscMemzero(&ropts_par, PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 77749d66b22SVijay Mahadevan ierr = PetscMemzero(&ropts_pargh, PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 77861eb6e27SVijay Mahadevan ierr = PetscMemzero(&ropts_dbg, PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 77961eb6e27SVijay Mahadevan 780e23c60ebSVijay Mahadevan /* do parallel read unless using only one processor */ 781a4717116SVijay Mahadevan if (numproc > 1) { 78249d66b22SVijay Mahadevan // ierr = 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":""));CHKERRQ(ierr); 78349d66b22SVijay Mahadevan ierr = 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;" : ""));CHKERRQ(ierr); 78449d66b22SVijay Mahadevan if (nghost) { 78549d66b22SVijay Mahadevan ierr = PetscSNPrintf(ropts_pargh, PETSC_MAX_PATH_LEN, "PARALLEL_GHOSTS=%d.0.%d;", dim, nghost);CHKERRQ(ierr); 78649d66b22SVijay Mahadevan } 7872e4e7c01SVijay Mahadevan } 7882e4e7c01SVijay Mahadevan 789c8d5365dSVijay Mahadevan if (dbglevel) { 790f90c3b0eSVijay Mahadevan if (numproc > 1) { 79149d66b22SVijay Mahadevan ierr = PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%d;DEBUG_PIO=%d;", dbglevel, dbglevel);CHKERRQ(ierr); 79261eb6e27SVijay Mahadevan } 79361eb6e27SVijay Mahadevan else { 79449d66b22SVijay Mahadevan ierr = PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%d;", dbglevel);CHKERRQ(ierr); 795f90c3b0eSVijay Mahadevan } 796c8d5365dSVijay Mahadevan } 79751d15aeeSVijay Mahadevan 79849d66b22SVijay Mahadevan ierr = 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 : ""));CHKERRQ(ierr); 799f90c3b0eSVijay Mahadevan *read_opts = ropts; 80051d15aeeSVijay Mahadevan PetscFunctionReturn(0); 80151d15aeeSVijay Mahadevan } 80251d15aeeSVijay Mahadevan 80351d15aeeSVijay Mahadevan 804aa859218SVijay Mahadevan /*@ 805aa859218SVijay Mahadevan DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file. 806aa859218SVijay Mahadevan 807aa859218SVijay Mahadevan Collective on MPI_Comm 808aa859218SVijay Mahadevan 809aa859218SVijay Mahadevan Input Parameters: 810aa859218SVijay Mahadevan + comm - The communicator for the DM object 811aa859218SVijay Mahadevan . dim - The spatial dimension 812b8ecf6d3SVijay Mahadevan . filename - The name of the mesh file to be loaded 813b8ecf6d3SVijay Mahadevan . usrreadopts - The options string to read a MOAB mesh. 814aa859218SVijay Mahadevan Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo) 815aa859218SVijay Mahadevan 816aa859218SVijay Mahadevan Output Parameter: 817aa859218SVijay Mahadevan . dm - The DM object 818aa859218SVijay Mahadevan 819aa859218SVijay Mahadevan Level: beginner 820aa859218SVijay Mahadevan 821aa859218SVijay Mahadevan .keywords: DM, create 822b8ecf6d3SVijay Mahadevan 823aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh() 824aa859218SVijay Mahadevan @*/ 82549d66b22SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm, PetscInt dim, PetscInt nghost, const char* filename, const char* usrreadopts, DM *dm) 82651d15aeeSVijay Mahadevan { 827a4717116SVijay Mahadevan moab::ErrorCode merr; 8282e4e7c01SVijay Mahadevan PetscInt nprocs; 829a4717116SVijay Mahadevan DM_Moab *dmmoab; 830a4717116SVijay Mahadevan moab::Interface *mbiface; 8319daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 832a4717116SVijay Mahadevan moab::ParallelComm *pcomm; 8339daf19fdSVijay Mahadevan #endif 834a4717116SVijay Mahadevan moab::Range verts, elems; 8357023aa44SVijay Mahadevan const char *readopts; 836a4717116SVijay Mahadevan PetscErrorCode ierr; 83751d15aeeSVijay Mahadevan 83851d15aeeSVijay Mahadevan PetscFunctionBegin; 83949d66b22SVijay Mahadevan PetscValidPointer(dm, 6); 84051d15aeeSVijay Mahadevan 841a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */ 8429daf19fdSVijay Mahadevan ierr = DMMoabCreateMoab(comm, NULL, NULL, NULL, dm);CHKERRQ(ierr); 84351d15aeeSVijay Mahadevan 844a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */ 845a4717116SVijay Mahadevan dmmoab = (DM_Moab*)(*dm)->data; 846a4717116SVijay Mahadevan mbiface = dmmoab->mbiface; 8479daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 848a4717116SVijay Mahadevan pcomm = dmmoab->pcomm; 849a4717116SVijay Mahadevan nprocs = pcomm->size(); 8509daf19fdSVijay Mahadevan #else 8519daf19fdSVijay Mahadevan nprocs = 1; 8529daf19fdSVijay Mahadevan #endif 853aa859218SVijay Mahadevan /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */ 854c8d5365dSVijay Mahadevan dmmoab->dim = dim; 85549d66b22SVijay Mahadevan dmmoab->nghostrings = nghost; 856304006b3SVijay Mahadevan dmmoab->refct = 1; 857a4717116SVijay Mahadevan 858b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */ 859b5410836SVijay Mahadevan merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset); MBERR("Creating file set failed", merr); 860b5410836SVijay Mahadevan 861a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */ 86249d66b22SVijay Mahadevan ierr = DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, nghost, dmmoab->read_mode, 8632e4e7c01SVijay Mahadevan dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts);CHKERRQ(ierr); 8647023aa44SVijay Mahadevan 865e5136372SVijay Mahadevan PetscInfo2(*dm, "Reading file %s with options: %s\n", filename, readopts); 866a4717116SVijay Mahadevan 867a4717116SVijay Mahadevan /* Load the mesh from a file. */ 86849d66b22SVijay Mahadevan if (dmmoab->fileset) { 8697023aa44SVijay Mahadevan merr = mbiface->load_file(filename, &dmmoab->fileset, readopts); MBERRVM(mbiface, "Reading MOAB file failed.", merr); 87049d66b22SVijay Mahadevan } 87149d66b22SVijay Mahadevan else { 87249d66b22SVijay Mahadevan merr = mbiface->load_file(filename, 0, readopts); MBERRVM(mbiface, "Reading MOAB file failed.", merr); 87349d66b22SVijay Mahadevan } 874a4717116SVijay Mahadevan 8759daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 8766e40195eSVijay Mahadevan /* Reassign global IDs on all entities. */ 877e5136372SVijay Mahadevan merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, true, true, true); MBERRNM(merr); 8789daf19fdSVijay Mahadevan #endif 879e5136372SVijay Mahadevan 880e5136372SVijay Mahadevan /* load the local vertices */ 881e5136372SVijay Mahadevan merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true); MBERRNM(merr); 882e5136372SVijay Mahadevan /* load the local elements */ 883e5136372SVijay Mahadevan merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true); MBERRNM(merr); 884e5136372SVijay Mahadevan 8859daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 886e5136372SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags 887e5136372SVijay Mahadevan on all of the ghost vertexes */ 888e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag, verts); MBERRV(mbiface, merr); 889e5136372SVijay Mahadevan merr = pcomm->exchange_tags(dmmoab->ltog_tag, elems); MBERRV(mbiface, merr); 890a4717116SVijay Mahadevan merr = pcomm->collective_sync_partition(); MBERR("Collective sync failed", merr); 8919daf19fdSVijay Mahadevan #endif 892a4717116SVijay Mahadevan 893e5136372SVijay Mahadevan PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size()); 8946acfe860SVijay Mahadevan ierr = PetscFree(readopts);CHKERRQ(ierr); 89551d15aeeSVijay Mahadevan PetscFunctionReturn(0); 89651d15aeeSVijay Mahadevan } 89751d15aeeSVijay Mahadevan 898304006b3SVijay Mahadevan 899304006b3SVijay Mahadevan #undef __FUNCT__ 900304006b3SVijay Mahadevan #define __FUNCT__ "DMMoabRenumberMeshEntities" 901304006b3SVijay Mahadevan /*@ 902304006b3SVijay Mahadevan DMMoabRenumberMeshEntities - Order and number all entities (vertices->elements) to be contiguously ordered 903304006b3SVijay Mahadevan in parallel 904304006b3SVijay Mahadevan 905304006b3SVijay Mahadevan Collective on MPI_Comm 906304006b3SVijay Mahadevan 907304006b3SVijay Mahadevan Input Parameters: 908304006b3SVijay Mahadevan + dm - The DM object 909304006b3SVijay Mahadevan 910304006b3SVijay Mahadevan Level: advanced 911304006b3SVijay Mahadevan 912304006b3SVijay Mahadevan .keywords: DM, reorder 913304006b3SVijay Mahadevan 914304006b3SVijay Mahadevan .seealso: DMSetUp(), DMCreate() 915304006b3SVijay Mahadevan @*/ 916304006b3SVijay Mahadevan PetscErrorCode DMMoabRenumberMeshEntities(DM dm) 917304006b3SVijay Mahadevan { 918304006b3SVijay Mahadevan moab::Range verts; 919304006b3SVijay Mahadevan 920304006b3SVijay Mahadevan PetscFunctionBegin; 921304006b3SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 922304006b3SVijay Mahadevan 9239daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 924304006b3SVijay Mahadevan /* Insert new points */ 9259daf19fdSVijay Mahadevan moab::ErrorCode merr; 9269daf19fdSVijay Mahadevan merr = ((DM_Moab*) dm->data)->pcomm->assign_global_ids(((DM_Moab*) dm->data)->fileset, 3, 0, false, true, false); MBERRNM(merr); 9279daf19fdSVijay Mahadevan #endif 928304006b3SVijay Mahadevan PetscFunctionReturn(0); 929304006b3SVijay Mahadevan } 930304006b3SVijay Mahadevan 931