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
DMMoab_SetTensorElementConnectivity_Private(DMMoabMeshGeneratorCtx & genCtx,PetscInt offset,PetscInt corner,std::vector<PetscInt> & subent_conn,moab::EntityHandle * connectivity)2766976f2fSJacob Faibussowitsch static PetscInt DMMoab_SetTensorElementConnectivity_Private(DMMoabMeshGeneratorCtx &genCtx, PetscInt offset, PetscInt corner, std::vector<PetscInt> &subent_conn, moab::EntityHandle *connectivity)
28d71ae5a4SJacob Faibussowitsch {
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
DMMoab_SetSimplexElementConnectivity_Private(DMMoabMeshGeneratorCtx & genCtx,PetscInt subelem,PetscInt offset,PetscInt corner,std::vector<PetscInt> & subent_conn,moab::EntityHandle * connectivity)6166976f2fSJacob Faibussowitsch static PetscInt DMMoab_SetSimplexElementConnectivity_Private(DMMoabMeshGeneratorCtx &genCtx, PetscInt subelem, PetscInt offset, PetscInt corner, std::vector<PetscInt> &subent_conn, moab::EntityHandle *connectivity)
62d71ae5a4SJacob Faibussowitsch {
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;
851baa6e33SBarry 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 }
901baa6e33SBarry 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;
951baa6e33SBarry 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 }
1001baa6e33SBarry 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
DMMoab_SetElementConnectivity_Private(DMMoabMeshGeneratorCtx & genCtx,PetscInt offset,PetscInt corner,moab::EntityHandle * connectivity)17866976f2fSJacob Faibussowitsch static std::pair<PetscInt, PetscInt> DMMoab_SetElementConnectivity_Private(DMMoabMeshGeneratorCtx &genCtx, PetscInt offset, PetscInt corner, moab::EntityHandle *connectivity)
179d71ae5a4SJacob Faibussowitsch {
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 }
1911baa6e33SBarry 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
DMMoab_GenerateVertices_Private(moab::Interface * mbImpl,moab::ReadUtilIface * iface,DMMoabMeshGeneratorCtx & genCtx,PetscInt m,PetscInt n,PetscInt k,PetscInt a,PetscInt b,PetscInt c,moab::Tag & global_id_tag,moab::EntityHandle & startv,moab::Range & uverts)19866976f2fSJacob Faibussowitsch static PetscErrorCode DMMoab_GenerateVertices_Private(moab::Interface *mbImpl, moab::ReadUtilIface *iface, DMMoabMeshGeneratorCtx &genCtx, PetscInt m, PetscInt n, PetscInt k, PetscInt a, PetscInt b, PetscInt c, moab::Tag &global_id_tag, moab::EntityHandle &startv, moab::Range &uverts)
199d71ae5a4SJacob Faibussowitsch {
20049d66b22SVijay Mahadevan PetscInt x, y, z, ix, nnodes;
20149d66b22SVijay Mahadevan PetscInt ii, jj, kk;
20249d66b22SVijay Mahadevan std::vector<PetscReal *> arrays;
203755f3dfbSVijay Mahadevan PetscInt *gids;
204ce27a4eeSVijay Mahadevan moab::ErrorCode merr;
205ce27a4eeSVijay Mahadevan
206ce27a4eeSVijay Mahadevan PetscFunctionBegin;
207ce27a4eeSVijay Mahadevan /* we will generate (q*block+1)^3 vertices, and block^3 hexas; q is 1 for linear, 2 for quadratic
208ce27a4eeSVijay Mahadevan * the global id of the vertices will come from m, n, k, a, b, c
209ce27a4eeSVijay Mahadevan * x will vary from m*A*q*block + a*q*block to m*A*q*block+(a+1)*q*block etc.
210ce27a4eeSVijay Mahadevan */
211ce27a4eeSVijay Mahadevan nnodes = genCtx.blockSizeVertexXYZ[0] * (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] * (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1) : 1);
2129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnodes, &gids));
213ce27a4eeSVijay Mahadevan
2149371c9d4SSatish Balay merr = iface->get_node_coords(3, nnodes, 0, startv, arrays);
2159371c9d4SSatish Balay 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));
2493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
250ce27a4eeSVijay Mahadevan }
251ce27a4eeSVijay Mahadevan
DMMoab_GenerateElements_Private(moab::Interface * mbImpl,moab::ReadUtilIface * iface,DMMoabMeshGeneratorCtx & genCtx,PetscInt m,PetscInt n,PetscInt k,PetscInt a,PetscInt b,PetscInt c,moab::Tag & global_id_tag,moab::EntityHandle startv,moab::Range & cells)25266976f2fSJacob Faibussowitsch static PetscErrorCode DMMoab_GenerateElements_Private(moab::Interface *mbImpl, moab::ReadUtilIface *iface, DMMoabMeshGeneratorCtx &genCtx, PetscInt m, PetscInt n, PetscInt k, PetscInt a, PetscInt b, PetscInt c, moab::Tag &global_id_tag, moab::EntityHandle startv, moab::Range &cells)
253d71ae5a4SJacob Faibussowitsch {
254ce27a4eeSVijay Mahadevan moab::ErrorCode merr;
255ce27a4eeSVijay Mahadevan PetscInt ix, ie, xe, ye, ze;
256ce27a4eeSVijay Mahadevan PetscInt ii, jj, kk, nvperelem;
25749d66b22SVijay Mahadevan PetscInt simplices_per_tensor[4] = {0, 1, 2, 6};
258ce27a4eeSVijay 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);*/
259ce27a4eeSVijay Mahadevan PetscInt nelems = ntensorelems;
26063d025dbSVijay Mahadevan moab::EntityHandle starte; /* connectivity */
261ce27a4eeSVijay Mahadevan moab::EntityHandle *conn;
262ce27a4eeSVijay Mahadevan
263ce27a4eeSVijay Mahadevan PetscFunctionBegin;
264ce27a4eeSVijay Mahadevan switch (genCtx.dim) {
265ce27a4eeSVijay Mahadevan case 1:
266ce27a4eeSVijay Mahadevan nvperelem = 2;
2679371c9d4SSatish Balay merr = iface->get_element_connect(nelems, 2, moab::MBEDGE, 0, starte, conn);
2689371c9d4SSatish Balay 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];
2749371c9d4SSatish Balay merr = iface->get_element_connect(nelems, 3, moab::MBTRI, 0, starte, conn);
2759371c9d4SSatish Balay MBERR("Can't get TRI3 element connectivity.", merr);
2761baa6e33SBarry Smith } else {
277ce27a4eeSVijay Mahadevan nvperelem = 4;
2789371c9d4SSatish Balay merr = iface->get_element_connect(nelems, 4, moab::MBQUAD, 0, starte, conn);
2799371c9d4SSatish Balay MBERR("Can't get QUAD4 element connectivity.", merr);
280ce27a4eeSVijay Mahadevan }
281ce27a4eeSVijay Mahadevan break;
282ce27a4eeSVijay Mahadevan case 3:
283ce27a4eeSVijay Mahadevan default:
284ce27a4eeSVijay Mahadevan if (genCtx.simplex) {
285ce27a4eeSVijay Mahadevan nvperelem = 4;
28649d66b22SVijay Mahadevan nelems = ntensorelems * simplices_per_tensor[genCtx.dim];
2879371c9d4SSatish Balay merr = iface->get_element_connect(nelems, 4, moab::MBTET, 0, starte, conn);
2889371c9d4SSatish Balay MBERR("Can't get TET4 element connectivity.", merr);
2891baa6e33SBarry Smith } else {
290ce27a4eeSVijay Mahadevan nvperelem = 8;
2919371c9d4SSatish Balay merr = iface->get_element_connect(nelems, 8, moab::MBHEX, 0, starte, conn);
2929371c9d4SSatish Balay MBERR("Can't get HEX8 element connectivity.", merr);
293ce27a4eeSVijay Mahadevan }
294ce27a4eeSVijay Mahadevan break;
295ce27a4eeSVijay Mahadevan }
296ce27a4eeSVijay Mahadevan
29763d025dbSVijay Mahadevan ix = ie = 0; /* index now in the elements, for global ids */
298ce27a4eeSVijay Mahadevan
299ce27a4eeSVijay Mahadevan /* create a temporary range to store local element handles */
300ce27a4eeSVijay Mahadevan moab::Range tmp(starte, starte + nelems - 1);
30149d66b22SVijay Mahadevan std::vector<PetscInt> gids(nelems);
302ce27a4eeSVijay Mahadevan
303ce27a4eeSVijay Mahadevan /* identify the elements at the lower corner, for their global ids */
304ce27a4eeSVijay Mahadevan xe = m * genCtx.A * genCtx.blockSizeElementXYZ[0] + a * genCtx.blockSizeElementXYZ[0];
305ce27a4eeSVijay Mahadevan ye = (genCtx.dim > 1 ? n * genCtx.B * genCtx.blockSizeElementXYZ[1] + b * genCtx.blockSizeElementXYZ[1] : 0);
306ce27a4eeSVijay Mahadevan ze = (genCtx.dim > 2 ? k * genCtx.C * genCtx.blockSizeElementXYZ[2] + c * genCtx.blockSizeElementXYZ[2] : 0);
307ce27a4eeSVijay Mahadevan
308ce27a4eeSVijay Mahadevan /* create owned elements requested by genCtx */
309ce27a4eeSVijay Mahadevan for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeElementXYZ[2] : 1); kk++) {
310ce27a4eeSVijay Mahadevan for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeElementXYZ[1] : 1); jj++) {
311ce27a4eeSVijay Mahadevan for (ii = 0; ii < genCtx.blockSizeElementXYZ[0]; ii++) {
312ce27a4eeSVijay Mahadevan moab::EntityHandle corner = startv + genCtx.q * ii + genCtx.q * jj * genCtx.ystride + genCtx.q * kk * genCtx.zstride;
313ce27a4eeSVijay Mahadevan
31449d66b22SVijay Mahadevan std::pair<PetscInt, PetscInt> entoffset = DMMoab_SetElementConnectivity_Private(genCtx, ix, corner, conn);
315ce27a4eeSVijay Mahadevan
31649d66b22SVijay Mahadevan for (PetscInt j = 0; j < entoffset.second; j++) {
317ce27a4eeSVijay Mahadevan /* The entity handle for the particular element -> if we want to set some tags is
318ce27a4eeSVijay Mahadevan moab::EntityHandle eh = starte + ie + j;
319ce27a4eeSVijay Mahadevan */
320ce27a4eeSVijay Mahadevan gids[ie + j] = 1 + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney));
32163d025dbSVijay Mahadevan /* gids[ie+j] = ie + j + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)); */
32263d025dbSVijay Mahadevan /* gids[ie+j] = 1 + ie; */
32363d025dbSVijay Mahadevan /* ie++; */
324ce27a4eeSVijay Mahadevan }
325ce27a4eeSVijay Mahadevan
326ce27a4eeSVijay Mahadevan ix += entoffset.first;
32749d66b22SVijay Mahadevan ie += entoffset.second;
328ce27a4eeSVijay Mahadevan }
329ce27a4eeSVijay Mahadevan }
330ce27a4eeSVijay Mahadevan }
331ce27a4eeSVijay Mahadevan if (genCtx.adjEnts) { /* we need to update adjacencies now, because some elements are new */
3329371c9d4SSatish Balay merr = iface->update_adjacencies(starte, nelems, nvperelem, conn);
3339371c9d4SSatish Balay MBERR("Can't update adjacencies", merr);
334ce27a4eeSVijay Mahadevan }
335ce27a4eeSVijay Mahadevan tmp.swap(cells);
3369371c9d4SSatish Balay merr = mbImpl->tag_set_data(global_id_tag, cells, &gids[0]);
3379371c9d4SSatish Balay MBERR("Can't set global ids to elements.", merr);
3383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
339ce27a4eeSVijay Mahadevan }
340ce27a4eeSVijay Mahadevan
DMMBUtil_InitializeOptions(DMMoabMeshGeneratorCtx & genCtx,PetscInt dim,PetscBool simplex,PetscInt rank,PetscInt nprocs,const PetscReal * bounds,PetscInt nelems)34166976f2fSJacob Faibussowitsch static PetscErrorCode DMMBUtil_InitializeOptions(DMMoabMeshGeneratorCtx &genCtx, PetscInt dim, PetscBool simplex, PetscInt rank, PetscInt nprocs, const PetscReal *bounds, PetscInt nelems)
342d71ae5a4SJacob Faibussowitsch {
343ce27a4eeSVijay Mahadevan PetscFunctionBegin;
344ce27a4eeSVijay Mahadevan /* Initialize all genCtx data */
345ce27a4eeSVijay Mahadevan genCtx.dim = dim;
346ce27a4eeSVijay Mahadevan genCtx.simplex = simplex;
34749d66b22SVijay Mahadevan genCtx.newMergeMethod = genCtx.keep_skins = genCtx.adjEnts = true;
348ce27a4eeSVijay Mahadevan /* determine other global quantities for the mesh used for nodes increments */
349ce27a4eeSVijay Mahadevan genCtx.q = 1;
350755f3dfbSVijay Mahadevan genCtx.fraction = genCtx.remainder = genCtx.cumfraction = 0;
351ce27a4eeSVijay Mahadevan
352ce27a4eeSVijay Mahadevan if (!genCtx.usrxyzgrid) { /* not overridden by genCtx - assume nele equally and that genCtx wants a uniform cube mesh */
353755f3dfbSVijay Mahadevan
354755f3dfbSVijay Mahadevan genCtx.fraction = nelems / nprocs; /* partition only by the largest dimension */
355755f3dfbSVijay Mahadevan genCtx.remainder = nelems % nprocs; /* remainder after partition which gets evenly distributed by round-robin */
356755f3dfbSVijay Mahadevan genCtx.cumfraction = (rank > 0 ? (genCtx.fraction) * (rank) + (rank - 1 < genCtx.remainder ? rank : genCtx.remainder) : 0);
357755f3dfbSVijay Mahadevan if (rank < genCtx.remainder) /* This process gets "fraction+1" elements */
358755f3dfbSVijay Mahadevan genCtx.fraction++;
359755f3dfbSVijay Mahadevan
3609566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Fraction = %" PetscInt_FMT ", Remainder = %" PetscInt_FMT ", Cumulative fraction = %" PetscInt_FMT "\n", genCtx.fraction, genCtx.remainder, genCtx.cumfraction));
361755f3dfbSVijay Mahadevan switch (genCtx.dim) {
362755f3dfbSVijay Mahadevan case 1:
363755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = genCtx.fraction;
364755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = 1;
365755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = 1;
366755f3dfbSVijay Mahadevan break;
367755f3dfbSVijay Mahadevan case 2:
368755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = nelems;
369755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = genCtx.fraction;
370755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = 1;
371755f3dfbSVijay Mahadevan break;
372755f3dfbSVijay Mahadevan case 3:
373755f3dfbSVijay Mahadevan default:
374755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[0] = nelems;
375755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[1] = nelems;
376755f3dfbSVijay Mahadevan genCtx.blockSizeElementXYZ[2] = genCtx.fraction;
377755f3dfbSVijay Mahadevan break;
37849d66b22SVijay Mahadevan }
379ce27a4eeSVijay Mahadevan }
380ce27a4eeSVijay Mahadevan
381ce27a4eeSVijay Mahadevan /* partition only by the largest dimension */
382ce27a4eeSVijay Mahadevan /* Total number of local elements := genCtx.blockSizeElementXYZ[0]*(genCtx.dim>1? genCtx.blockSizeElementXYZ[1]*(genCtx.dim>2 ? genCtx.blockSizeElementXYZ[2]:1) :1); */
383ce27a4eeSVijay Mahadevan if (bounds) {
3845f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < 6; i++) genCtx.xyzbounds[i] = bounds[i];
3851baa6e33SBarry Smith } else {
386ce27a4eeSVijay Mahadevan genCtx.xyzbounds[0] = genCtx.xyzbounds[2] = genCtx.xyzbounds[4] = 0.0;
387ce27a4eeSVijay Mahadevan genCtx.xyzbounds[1] = genCtx.xyzbounds[3] = genCtx.xyzbounds[5] = 1.0;
388ce27a4eeSVijay Mahadevan }
389ce27a4eeSVijay Mahadevan
390ce27a4eeSVijay Mahadevan if (!genCtx.usrprocgrid) {
391ce27a4eeSVijay Mahadevan switch (genCtx.dim) {
392ce27a4eeSVijay Mahadevan case 1:
393ce27a4eeSVijay Mahadevan genCtx.M = nprocs;
394ce27a4eeSVijay Mahadevan genCtx.N = genCtx.K = 1;
395ce27a4eeSVijay Mahadevan break;
396ce27a4eeSVijay Mahadevan case 2:
397755f3dfbSVijay Mahadevan genCtx.N = nprocs;
39863d025dbSVijay Mahadevan genCtx.M = genCtx.K = 1;
399ce27a4eeSVijay Mahadevan break;
400ce27a4eeSVijay Mahadevan default:
401755f3dfbSVijay Mahadevan genCtx.K = nprocs;
40263d025dbSVijay Mahadevan genCtx.M = genCtx.N = 1;
403ce27a4eeSVijay Mahadevan break;
404ce27a4eeSVijay Mahadevan }
405ce27a4eeSVijay Mahadevan }
406ce27a4eeSVijay Mahadevan
407ad540459SPierre Jolivet if (!genCtx.usrrefgrid) genCtx.A = genCtx.B = genCtx.C = 1;
408ce27a4eeSVijay Mahadevan
409ce27a4eeSVijay Mahadevan /* more default values */
410ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.ney = genCtx.nez = 0;
411ce27a4eeSVijay Mahadevan genCtx.xstride = genCtx.ystride = genCtx.zstride = 0;
412ce27a4eeSVijay Mahadevan genCtx.NX = genCtx.NY = genCtx.NZ = 0;
413ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.ney = genCtx.nez = 0;
414ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 1;
415ce27a4eeSVijay Mahadevan
416ce27a4eeSVijay Mahadevan switch (genCtx.dim) {
417ce27a4eeSVijay Mahadevan case 3:
418ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1;
419ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1;
420ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[2] = genCtx.q * genCtx.blockSizeElementXYZ[2] + 1;
421ce27a4eeSVijay Mahadevan
422ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */
423755f3dfbSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */
424ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1);
425ce27a4eeSVijay Mahadevan genCtx.xstride = 1;
426ce27a4eeSVijay Mahadevan genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1]; /* number of elements in y direction .... */
427755f3dfbSVijay Mahadevan genCtx.dy = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */
428ce27a4eeSVijay Mahadevan genCtx.NY = (genCtx.q * genCtx.ney + 1);
429ce27a4eeSVijay Mahadevan genCtx.ystride = genCtx.blockSizeVertexXYZ[0];
430ce27a4eeSVijay Mahadevan genCtx.nez = genCtx.K * genCtx.C * genCtx.blockSizeElementXYZ[2]; /* number of elements in z direction .... */
431755f3dfbSVijay Mahadevan genCtx.dz = (genCtx.xyzbounds[5] - genCtx.xyzbounds[4]) / (nelems * genCtx.q); /* distance between 2 nodes in z direction */
432ce27a4eeSVijay Mahadevan genCtx.NZ = (genCtx.q * genCtx.nez + 1);
433ce27a4eeSVijay Mahadevan genCtx.zstride = genCtx.blockSizeVertexXYZ[0] * genCtx.blockSizeVertexXYZ[1];
434ce27a4eeSVijay Mahadevan break;
435ce27a4eeSVijay Mahadevan case 2:
436ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1;
437ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1;
4389daf19fdSVijay Mahadevan genCtx.blockSizeVertexXYZ[2] = 0;
439ce27a4eeSVijay Mahadevan
440ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */
441ce27a4eeSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (genCtx.nex * genCtx.q); /* distance between 2 nodes in x direction */
442ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1);
443ce27a4eeSVijay Mahadevan genCtx.xstride = 1;
444ce27a4eeSVijay Mahadevan genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1]; /* number of elements in y direction .... */
445755f3dfbSVijay Mahadevan genCtx.dy = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */
446ce27a4eeSVijay Mahadevan genCtx.NY = (genCtx.q * genCtx.ney + 1);
447ce27a4eeSVijay Mahadevan genCtx.ystride = genCtx.blockSizeVertexXYZ[0];
448ce27a4eeSVijay Mahadevan break;
449ce27a4eeSVijay Mahadevan case 1:
4509daf19fdSVijay Mahadevan genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 0;
451ce27a4eeSVijay Mahadevan genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1;
452ce27a4eeSVijay Mahadevan
453ce27a4eeSVijay Mahadevan genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0]; /* number of elements in x direction, used for global id on element */
454755f3dfbSVijay Mahadevan genCtx.dx = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */
455ce27a4eeSVijay Mahadevan genCtx.NX = (genCtx.q * genCtx.nex + 1);
456ce27a4eeSVijay Mahadevan genCtx.xstride = 1;
457ce27a4eeSVijay Mahadevan break;
458ce27a4eeSVijay Mahadevan }
459ce27a4eeSVijay Mahadevan
460755f3dfbSVijay Mahadevan /* Lets check for some valid input */
4615f80ce2aSJacob Faibussowitsch PetscCheck(genCtx.dim >= 1 && genCtx.dim <= 3, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid topological dimension specified: %" PetscInt_FMT ".", genCtx.dim);
4629371c9d4SSatish Balay 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,
4639371c9d4SSatish Balay genCtx.N, genCtx.K, nprocs);
464755f3dfbSVijay Mahadevan /* validate the bounds data */
4655f80ce2aSJacob 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]);
4661dca8a05SBarry 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]);
4671dca8a05SBarry 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]);
468755f3dfbSVijay Mahadevan
4699566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local elements:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.blockSizeElementXYZ[0], genCtx.blockSizeElementXYZ[1], genCtx.blockSizeElementXYZ[2]));
4709566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local vertices:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.blockSizeVertexXYZ[0], genCtx.blockSizeVertexXYZ[1], genCtx.blockSizeVertexXYZ[2]));
4719566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local blocks/processors := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.A, genCtx.B, genCtx.C));
4729566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local processors := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.M, genCtx.N, genCtx.K));
4739566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local nexyz:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.nex, genCtx.ney, genCtx.nez));
4749566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local delxyz:= %g, %g, %g\n", genCtx.dx, genCtx.dy, genCtx.dz));
4759566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Local strides:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.xstride, genCtx.ystride, genCtx.zstride));
4763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
477ce27a4eeSVijay Mahadevan }
478ce27a4eeSVijay Mahadevan
479cab5ea25SPierre Jolivet /*@C
480ce27a4eeSVijay Mahadevan DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with genCtx specified bounds.
481aa859218SVijay Mahadevan
482d083f849SBarry Smith Collective
483aa859218SVijay Mahadevan
484aa859218SVijay Mahadevan Input Parameters:
485aa859218SVijay Mahadevan + comm - The communicator for the DM object
486aa859218SVijay Mahadevan . dim - The spatial dimension
48720f4b53cSBarry Smith . useSimplex - use a simplex mesh
488b8ecf6d3SVijay 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
489b8ecf6d3SVijay Mahadevan . nele - The number of discrete elements in each direction
49020f4b53cSBarry Smith - nghost - The number of ghosted layers needed in the partitioned mesh
491aa859218SVijay Mahadevan
492aa859218SVijay Mahadevan Output Parameter:
49320f4b53cSBarry Smith . dm - The `DM` object
494aa859218SVijay Mahadevan
495aa859218SVijay Mahadevan Level: beginner
496aa859218SVijay Mahadevan
497db781477SPatrick Sanan .seealso: `DMSetType()`, `DMCreate()`, `DMMoabLoadFromFile()`
498aa859218SVijay Mahadevan @*/
DMMoabCreateBoxMesh(MPI_Comm comm,PetscInt dim,PetscBool useSimplex,const PetscReal * bounds,PetscInt nele,PetscInt nghost,DM * dm)499d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal *bounds, PetscInt nele, PetscInt nghost, DM *dm)
500d71ae5a4SJacob Faibussowitsch {
50151d15aeeSVijay Mahadevan moab::ErrorCode merr;
5029daf19fdSVijay Mahadevan PetscInt a, b, c, n, global_size, global_rank;
50351d15aeeSVijay Mahadevan DM_Moab *dmmoab;
504ce27a4eeSVijay Mahadevan moab::Interface *mbImpl;
5059daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
50651d15aeeSVijay Mahadevan moab::ParallelComm *pcomm;
5079daf19fdSVijay Mahadevan #endif
508a4717116SVijay Mahadevan moab::ReadUtilIface *readMeshIface;
509ce27a4eeSVijay Mahadevan moab::Range verts, cells, edges, faces, adj, dim3, dim2;
510ce27a4eeSVijay Mahadevan DMMoabMeshGeneratorCtx genCtx;
511a4717116SVijay Mahadevan const PetscInt npts = nele + 1; /* Number of points in every dimension */
512ce27a4eeSVijay Mahadevan
5136ea892caSVijay Mahadevan moab::Tag global_id_tag, part_tag, geom_tag, mat_tag, dir_tag, neu_tag;
514ce27a4eeSVijay Mahadevan moab::Range ownedvtx, ownedelms, localvtxs, localelms;
515ce27a4eeSVijay Mahadevan moab::EntityHandle regionset;
516755f3dfbSVijay Mahadevan PetscInt ml = 0, nl = 0, kl = 0;
51751d15aeeSVijay Mahadevan
51851d15aeeSVijay Mahadevan PetscFunctionBegin;
5191dca8a05SBarry Smith PetscCheck(dim >= 1 && dim <= 3, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension argument for mesh: dim=[1,3].");
520e5136372SVijay Mahadevan
5219566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("GenerateMesh", DM_CLASSID, &genCtx.generateMesh));
5229566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("AddVertices", DM_CLASSID, &genCtx.generateVertices));
5239566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("AddElements", DM_CLASSID, &genCtx.generateElements));
5249566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("ParResolve", DM_CLASSID, &genCtx.parResolve));
5259566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(genCtx.generateMesh, 0, 0, 0, 0));
5269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &global_size));
527e5136372SVijay Mahadevan /* total number of vertices in all dimensions */
528e5136372SVijay Mahadevan n = pow(npts, dim);
529e5136372SVijay Mahadevan
530e5136372SVijay Mahadevan /* do some error checking */
53108401ef6SPierre Jolivet PetscCheck(n >= 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be >= 2.");
53208401ef6SPierre 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.");
53308401ef6SPierre Jolivet PetscCheck(nghost >= 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of ghost layers cannot be negative.");
534e5136372SVijay Mahadevan
535a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */
5369566063dSJacob Faibussowitsch PetscCall(DMMoabCreateMoab(comm, NULL, NULL, NULL, dm));
53751d15aeeSVijay Mahadevan
538a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */
53951d15aeeSVijay Mahadevan dmmoab = (DM_Moab *)(*dm)->data;
540ce27a4eeSVijay Mahadevan mbImpl = dmmoab->mbiface;
5419daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
54251d15aeeSVijay Mahadevan pcomm = dmmoab->pcomm;
543ce27a4eeSVijay Mahadevan global_rank = pcomm->rank();
5449daf19fdSVijay Mahadevan #else
5459daf19fdSVijay Mahadevan global_rank = 0;
5469daf19fdSVijay Mahadevan global_size = 1;
5479daf19fdSVijay Mahadevan #endif
5489daf19fdSVijay Mahadevan global_id_tag = dmmoab->ltog_tag;
549c8d5365dSVijay Mahadevan dmmoab->dim = dim;
55049d66b22SVijay Mahadevan dmmoab->nghostrings = nghost;
551304006b3SVijay Mahadevan dmmoab->refct = 1;
55251d15aeeSVijay Mahadevan
553b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */
5549371c9d4SSatish Balay merr = mbImpl->create_meshset(moab::MESHSET_SET, dmmoab->fileset);
5559371c9d4SSatish Balay MBERR("Creating file set failed", merr);
556b5410836SVijay Mahadevan
557a4717116SVijay Mahadevan /* No errors yet; proceed with building the mesh */
5589371c9d4SSatish Balay merr = mbImpl->query_interface(readMeshIface);
5599371c9d4SSatish Balay MBERRNM(merr);
56051d15aeeSVijay Mahadevan
561755f3dfbSVijay Mahadevan genCtx.M = genCtx.N = genCtx.K = 1;
562755f3dfbSVijay Mahadevan genCtx.A = genCtx.B = genCtx.C = 1;
563da8b1a3eSBarry Smith genCtx.blockSizeElementXYZ[0] = 0;
564da8b1a3eSBarry Smith genCtx.blockSizeElementXYZ[1] = 0;
565da8b1a3eSBarry Smith genCtx.blockSizeElementXYZ[2] = 0;
566755f3dfbSVijay Mahadevan
567d0609cedSBarry Smith PetscOptionsBegin(comm, "", "DMMoab Creation Options", "DMMOAB");
568ce27a4eeSVijay Mahadevan /* Handle DMMoab spatial resolution */
5699566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dmb_grid_x", "Number of grid points in x direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[0], &genCtx.blockSizeElementXYZ[0], &genCtx.usrxyzgrid));
5709566063dSJacob 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));
5719566063dSJacob 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));
572a4717116SVijay Mahadevan
5736aad120cSJose E. Roman /* Handle DMMoab parallel distribution */
5749566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dmb_processors_x", "Number of processors in x direction", "DMMoabSetNumProcs", genCtx.M, &genCtx.M, &genCtx.usrprocgrid));
5759566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsInt("-dmb_processors_y", "Number of processors in y direction", "DMMoabSetNumProcs", genCtx.N, &genCtx.N, &genCtx.usrprocgrid));
5769566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsInt("-dmb_processors_z", "Number of processors in z direction", "DMMoabSetNumProcs", genCtx.K, &genCtx.K, &genCtx.usrprocgrid));
57751d15aeeSVijay Mahadevan
578ce27a4eeSVijay Mahadevan /* Handle DMMoab block refinement */
5799566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dmb_refine_x", "Number of refinement blocks in x direction", "DMMoabSetRefinement", genCtx.A, &genCtx.A, &genCtx.usrrefgrid));
5809566063dSJacob Faibussowitsch if (dim > 1) PetscCall(PetscOptionsInt("-dmb_refine_y", "Number of refinement blocks in y direction", "DMMoabSetRefinement", genCtx.B, &genCtx.B, &genCtx.usrrefgrid));
5819566063dSJacob Faibussowitsch if (dim > 2) PetscCall(PetscOptionsInt("-dmb_refine_z", "Number of refinement blocks in z direction", "DMMoabSetRefinement", genCtx.C, &genCtx.C, &genCtx.usrrefgrid));
582d0609cedSBarry Smith PetscOptionsEnd();
58351d15aeeSVijay Mahadevan
5849566063dSJacob Faibussowitsch PetscCall(DMMBUtil_InitializeOptions(genCtx, dim, useSimplex, global_rank, global_size, bounds, nele));
58551d15aeeSVijay Mahadevan
58663a3b9bcSJacob 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);
58764e1c140SVijay Mahadevan
588ce27a4eeSVijay Mahadevan if (genCtx.adjEnts) genCtx.keep_skins = true; /* do not delete anything - consumes more memory */
58966f88a78SVijay Mahadevan
590ce27a4eeSVijay Mahadevan /* determine m, n, k for processor rank */
591755f3dfbSVijay Mahadevan ml = nl = kl = 0;
592755f3dfbSVijay Mahadevan switch (genCtx.dim) {
593d71ae5a4SJacob Faibussowitsch case 1:
594d71ae5a4SJacob Faibussowitsch ml = (genCtx.cumfraction);
595d71ae5a4SJacob Faibussowitsch break;
596d71ae5a4SJacob Faibussowitsch case 2:
597d71ae5a4SJacob Faibussowitsch nl = (genCtx.cumfraction);
598d71ae5a4SJacob Faibussowitsch break;
599755f3dfbSVijay Mahadevan default:
600755f3dfbSVijay Mahadevan kl = (genCtx.cumfraction) / genCtx.q / genCtx.blockSizeElementXYZ[2] / genCtx.C; //genCtx.K
601755f3dfbSVijay Mahadevan break;
602755f3dfbSVijay Mahadevan }
603e5136372SVijay Mahadevan
604ce27a4eeSVijay Mahadevan /*
605ce27a4eeSVijay Mahadevan * so there are a total of M * A * blockSizeElement elements in x direction (so M * A * blockSizeElement + 1 verts in x direction)
606ce27a4eeSVijay Mahadevan * so there are a total of N * B * blockSizeElement elements in y direction (so N * B * blockSizeElement + 1 verts in y direction)
607ce27a4eeSVijay Mahadevan * so there are a total of K * C * blockSizeElement elements in z direction (so K * C * blockSizeElement + 1 verts in z direction)
608c8d5365dSVijay Mahadevan
609ce27a4eeSVijay Mahadevan * there are ( M * A blockSizeElement) * ( N * B * blockSizeElement) * (K * C * blockSizeElement) hexas
610ce27a4eeSVijay Mahadevan * there are ( M * A * blockSizeElement + 1) * ( N * B * blockSizeElement + 1) * (K * C * blockSizeElement + 1) vertices
611ce27a4eeSVijay Mahadevan * x is the first dimension that varies
612ce27a4eeSVijay Mahadevan */
613e5136372SVijay Mahadevan
614ce27a4eeSVijay Mahadevan /* generate the block at (a, b, c); it will represent a partition , it will get a partition tag */
61549d66b22SVijay Mahadevan PetscInt dum_id = -1;
6169371c9d4SSatish Balay merr = mbImpl->tag_get_handle("GLOBAL_ID", 1, moab::MB_TYPE_INTEGER, global_id_tag);
6179371c9d4SSatish Balay MBERR("Getting Global_ID Tag handle failed", merr);
61851d15aeeSVijay Mahadevan
6199371c9d4SSatish Balay merr = mbImpl->tag_get_handle(MATERIAL_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, mat_tag);
6209371c9d4SSatish Balay MBERR("Getting Material set Tag handle failed", merr);
6219371c9d4SSatish Balay merr = mbImpl->tag_get_handle(DIRICHLET_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, dir_tag);
6229371c9d4SSatish Balay MBERR("Getting Dirichlet set Tag handle failed", merr);
6239371c9d4SSatish Balay merr = mbImpl->tag_get_handle(NEUMANN_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, neu_tag);
6249371c9d4SSatish Balay MBERR("Getting Neumann set Tag handle failed", merr);
6256ea892caSVijay Mahadevan
6269371c9d4SSatish Balay merr = mbImpl->tag_get_handle("PARALLEL_PARTITION", 1, moab::MB_TYPE_INTEGER, part_tag, moab::MB_TAG_CREAT | moab::MB_TAG_SPARSE, &dum_id);
6279371c9d4SSatish Balay MBERR("Getting Partition Tag handle failed", merr);
62851d15aeeSVijay Mahadevan
6293a4aeca1SVijay Mahadevan /* lets create some sets */
6309371c9d4SSatish Balay 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);
6319371c9d4SSatish Balay MBERRNM(merr);
6329371c9d4SSatish Balay merr = mbImpl->create_meshset(moab::MESHSET_SET, regionset);
6339371c9d4SSatish Balay MBERRNM(merr);
6349566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(genCtx.generateMesh, 0, 0, 0, 0));
6353a4aeca1SVijay Mahadevan
636ce27a4eeSVijay Mahadevan for (a = 0; a < (genCtx.dim > 0 ? genCtx.A : genCtx.A); a++) {
637ce27a4eeSVijay Mahadevan for (b = 0; b < (genCtx.dim > 1 ? genCtx.B : 1); b++) {
638ce27a4eeSVijay Mahadevan for (c = 0; c < (genCtx.dim > 2 ? genCtx.C : 1); c++) {
63949d66b22SVijay Mahadevan moab::EntityHandle startv;
64049d66b22SVijay Mahadevan
6419566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(genCtx.generateVertices, 0, 0, 0, 0));
6429566063dSJacob Faibussowitsch PetscCall(DMMoab_GenerateVertices_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, verts));
6439566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(genCtx.generateVertices, 0, 0, 0, 0));
6443a4aeca1SVijay Mahadevan
6459566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(genCtx.generateElements, 0, 0, 0, 0));
6469566063dSJacob Faibussowitsch PetscCall(DMMoab_GenerateElements_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, cells));
6479566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(genCtx.generateElements, 0, 0, 0, 0));
6483a4aeca1SVijay Mahadevan
64949d66b22SVijay Mahadevan PetscInt part_num = 0;
65049d66b22SVijay Mahadevan switch (genCtx.dim) {
651d71ae5a4SJacob Faibussowitsch case 3:
652d71ae5a4SJacob Faibussowitsch part_num += (c + kl * genCtx.C) * (genCtx.M * genCtx.A * genCtx.N * genCtx.B);
653d71ae5a4SJacob Faibussowitsch case 2:
654d71ae5a4SJacob Faibussowitsch part_num += (b + nl * genCtx.B) * (genCtx.M * genCtx.A);
655d71ae5a4SJacob Faibussowitsch case 1:
656d71ae5a4SJacob Faibussowitsch part_num += (a + ml * genCtx.A);
657d71ae5a4SJacob Faibussowitsch break;
65849d66b22SVijay Mahadevan }
65949d66b22SVijay Mahadevan
66049d66b22SVijay Mahadevan moab::EntityHandle part_set;
6619371c9d4SSatish Balay merr = mbImpl->create_meshset(moab::MESHSET_SET, part_set);
6629371c9d4SSatish Balay MBERR("Can't create mesh set.", merr);
6633a4aeca1SVijay Mahadevan
6649371c9d4SSatish Balay merr = mbImpl->add_entities(part_set, verts);
6659371c9d4SSatish Balay MBERR("Can't add vertices to set.", merr);
6669371c9d4SSatish Balay merr = mbImpl->add_entities(part_set, cells);
6679371c9d4SSatish Balay MBERR("Can't add entities to set.", merr);
6689371c9d4SSatish Balay merr = mbImpl->add_entities(regionset, cells);
6699371c9d4SSatish Balay MBERR("Can't add entities to set.", merr);
670ce27a4eeSVijay Mahadevan
671ce27a4eeSVijay Mahadevan /* if needed, add all edges and faces */
6729371c9d4SSatish Balay if (genCtx.adjEnts) {
673ce27a4eeSVijay Mahadevan if (genCtx.dim > 1) {
6749371c9d4SSatish Balay merr = mbImpl->get_adjacencies(cells, 1, true, edges, moab::Interface::UNION);
6759371c9d4SSatish Balay MBERR("Can't get edges", merr);
6769371c9d4SSatish Balay merr = mbImpl->add_entities(part_set, edges);
6779371c9d4SSatish Balay MBERR("Can't add edges to partition set.", merr);
678ce27a4eeSVijay Mahadevan }
679ce27a4eeSVijay Mahadevan if (genCtx.dim > 2) {
6809371c9d4SSatish Balay merr = mbImpl->get_adjacencies(cells, 2, true, faces, moab::Interface::UNION);
6819371c9d4SSatish Balay MBERR("Can't get faces", merr);
6829371c9d4SSatish Balay merr = mbImpl->add_entities(part_set, faces);
6839371c9d4SSatish Balay MBERR("Can't add faces to partition set.", merr);
684ce27a4eeSVijay Mahadevan }
685ce27a4eeSVijay Mahadevan edges.clear();
686ce27a4eeSVijay Mahadevan faces.clear();
687ce27a4eeSVijay Mahadevan }
6889371c9d4SSatish Balay verts.clear();
6899371c9d4SSatish Balay cells.clear();
690ce27a4eeSVijay Mahadevan
6919371c9d4SSatish Balay merr = mbImpl->tag_set_data(part_tag, &part_set, 1, &part_num);
6929371c9d4SSatish Balay MBERR("Can't set part tag on set", merr);
69349d66b22SVijay Mahadevan if (dmmoab->fileset) {
6949371c9d4SSatish Balay merr = mbImpl->add_parent_child(dmmoab->fileset, part_set);
6959371c9d4SSatish Balay MBERR("Can't add part set to file set.", merr);
6969371c9d4SSatish Balay merr = mbImpl->unite_meshset(dmmoab->fileset, part_set);
6979371c9d4SSatish Balay MBERRNM(merr);
69849d66b22SVijay Mahadevan }
6999371c9d4SSatish Balay merr = mbImpl->add_entities(dmmoab->fileset, &part_set, 1);
7009371c9d4SSatish Balay MBERRNM(merr);
701ce27a4eeSVijay Mahadevan }
702ce27a4eeSVijay Mahadevan }
703ce27a4eeSVijay Mahadevan }
704ce27a4eeSVijay Mahadevan
7059371c9d4SSatish Balay merr = mbImpl->add_parent_child(dmmoab->fileset, regionset);
7069371c9d4SSatish Balay MBERRNM(merr);
707ce27a4eeSVijay Mahadevan
70849d66b22SVijay Mahadevan /* Only in parallel: resolve shared entities between processors and exchange ghost layers */
709ce27a4eeSVijay Mahadevan if (global_size > 1) {
7109566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(genCtx.parResolve, 0, 0, 0, 0));
71177a8c067SVijay Mahadevan
7129371c9d4SSatish Balay merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, genCtx.dim, cells);
7139371c9d4SSatish Balay MBERR("Can't get all d-dimensional elements.", merr);
7149371c9d4SSatish Balay merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 0, verts);
7159371c9d4SSatish Balay MBERR("Can't get all vertices.", merr);
716ce27a4eeSVijay Mahadevan
717ce27a4eeSVijay Mahadevan if (genCtx.A * genCtx.B * genCtx.C != 1) { // merge needed
71849d66b22SVijay Mahadevan moab::MergeMesh mm(mbImpl);
719ce27a4eeSVijay Mahadevan if (genCtx.newMergeMethod) {
7209371c9d4SSatish Balay merr = mm.merge_using_integer_tag(verts, global_id_tag);
7219371c9d4SSatish Balay MBERR("Can't merge with GLOBAL_ID tag", merr);
7229371c9d4SSatish Balay } else {
7239371c9d4SSatish Balay merr = mm.merge_entities(cells, 0.0001);
7249371c9d4SSatish Balay MBERR("Can't merge with coordinates", merr);
7253a4aeca1SVijay Mahadevan }
7263a4aeca1SVijay Mahadevan }
7273a4aeca1SVijay Mahadevan
7289daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
729a4717116SVijay Mahadevan /* check the handles */
7309371c9d4SSatish Balay merr = pcomm->check_all_shared_handles();
7319371c9d4SSatish Balay MBERRV(mbImpl, merr);
73251d15aeeSVijay Mahadevan
733a4717116SVijay Mahadevan /* resolve the shared entities by exchanging information to adjacent processors */
7349371c9d4SSatish Balay merr = pcomm->resolve_shared_ents(dmmoab->fileset, cells, dim, dim - 1, NULL, &global_id_tag);
7359371c9d4SSatish Balay MBERRV(mbImpl, merr);
73649d66b22SVijay Mahadevan if (dmmoab->fileset) {
7379371c9d4SSatish Balay merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false, &dmmoab->fileset);
7389371c9d4SSatish Balay MBERRV(mbImpl, merr);
7399371c9d4SSatish Balay } else {
7409371c9d4SSatish Balay merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false);
7419371c9d4SSatish Balay MBERRV(mbImpl, merr);
74249d66b22SVijay Mahadevan }
743c8d5365dSVijay Mahadevan
744e427d9c9SVijay Mahadevan /* Reassign global IDs on all entities. */
7459371c9d4SSatish Balay merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, false, true, false);
7469371c9d4SSatish Balay MBERRNM(merr);
7479daf19fdSVijay Mahadevan #endif
74877a8c067SVijay Mahadevan
7499566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(genCtx.parResolve, 0, 0, 0, 0));
75049d66b22SVijay Mahadevan }
7513f1c6e43SVijay Mahadevan
752ce27a4eeSVijay Mahadevan if (!genCtx.keep_skins) { // default is to delete the 1- and 2-dimensional entities
753ce27a4eeSVijay Mahadevan // delete all quads and edges
754ce27a4eeSVijay Mahadevan moab::Range toDelete;
755ce27a4eeSVijay Mahadevan if (genCtx.dim > 1) {
7569371c9d4SSatish Balay merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 1, toDelete);
7579371c9d4SSatish Balay MBERR("Can't get edges", merr);
758ce27a4eeSVijay Mahadevan }
7593f1c6e43SVijay Mahadevan
760ce27a4eeSVijay Mahadevan if (genCtx.dim > 2) {
7619371c9d4SSatish Balay merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 2, toDelete);
7629371c9d4SSatish Balay MBERR("Can't get faces", merr);
763ce27a4eeSVijay Mahadevan }
764c8d5365dSVijay Mahadevan
7659daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
7669371c9d4SSatish Balay merr = dmmoab->pcomm->delete_entities(toDelete);
7679371c9d4SSatish Balay MBERR("Can't delete entities", merr);
7689daf19fdSVijay Mahadevan #endif
769ce27a4eeSVijay Mahadevan }
7706ea892caSVijay Mahadevan
7716ea892caSVijay Mahadevan /* set geometric dimension tag for regions */
7729371c9d4SSatish Balay merr = mbImpl->tag_set_data(geom_tag, ®ionset, 1, &dmmoab->dim);
7739371c9d4SSatish Balay MBERRNM(merr);
7746ea892caSVijay Mahadevan /* set default material ID for regions */
7756ea892caSVijay Mahadevan int default_material = 1;
7769371c9d4SSatish Balay merr = mbImpl->tag_set_data(mat_tag, ®ionset, 1, &default_material);
7779371c9d4SSatish Balay MBERRNM(merr);
7786ea892caSVijay Mahadevan /*
7796ea892caSVijay Mahadevan int default_dbc = 0;
7806ea892caSVijay Mahadevan merr = mbImpl->tag_set_data(dir_tag, &vertexset, 1, &default_dbc);MBERRNM(merr);
7816ea892caSVijay Mahadevan */
7823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
78351d15aeeSVijay Mahadevan }
78451d15aeeSVijay Mahadevan
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)785a4e35b19SJacob Faibussowitsch static 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)
786d71ae5a4SJacob Faibussowitsch {
7876acfe860SVijay Mahadevan char *ropts;
78849d66b22SVijay Mahadevan char ropts_par[PETSC_MAX_PATH_LEN], ropts_pargh[PETSC_MAX_PATH_LEN];
78961eb6e27SVijay Mahadevan char ropts_dbg[PETSC_MAX_PATH_LEN];
79051d15aeeSVijay Mahadevan
791a4717116SVijay Mahadevan PetscFunctionBegin;
7929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PETSC_MAX_PATH_LEN, &ropts));
7939566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ropts_par, PETSC_MAX_PATH_LEN));
7949566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ropts_pargh, PETSC_MAX_PATH_LEN));
7959566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ropts_dbg, PETSC_MAX_PATH_LEN));
79661eb6e27SVijay Mahadevan
797e23c60ebSVijay Mahadevan /* do parallel read unless using only one processor */
798a4717116SVijay Mahadevan if (numproc > 1) {
7999566063dSJacob 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":"")));
80057508eceSPierre Jolivet 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;" : ""));
80148a46eb9SPierre Jolivet if (nghost) PetscCall(PetscSNPrintf(ropts_pargh, PETSC_MAX_PATH_LEN, "PARALLEL_GHOSTS=%" PetscInt_FMT ".0.%" PetscInt_FMT ";", dim, nghost));
8022e4e7c01SVijay Mahadevan }
8032e4e7c01SVijay Mahadevan
804c8d5365dSVijay Mahadevan if (dbglevel) {
805*ac530a7eSPierre Jolivet if (numproc > 1) PetscCall(PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%" PetscInt_FMT ";DEBUG_PIO=%" PetscInt_FMT ";", dbglevel, dbglevel));
806*ac530a7eSPierre Jolivet else PetscCall(PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%" PetscInt_FMT ";", dbglevel));
807c8d5365dSVijay Mahadevan }
80851d15aeeSVijay Mahadevan
80957508eceSPierre Jolivet 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 : ""));
810f90c3b0eSVijay Mahadevan *read_opts = ropts;
8113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
81251d15aeeSVijay Mahadevan }
81351d15aeeSVijay Mahadevan
814cab5ea25SPierre Jolivet /*@C
8151d27aa22SBarry Smith DMMoabLoadFromFile - Creates a `DMMOAB` object by loading the mesh from a user specified file
8161d27aa22SBarry Smith <https://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo>
817aa859218SVijay Mahadevan
818d083f849SBarry Smith Collective
819aa859218SVijay Mahadevan
820aa859218SVijay Mahadevan Input Parameters:
8211d27aa22SBarry Smith + comm - The communicator for the `DMOAB` object
822aa859218SVijay Mahadevan . dim - The spatial dimension
82320f4b53cSBarry Smith . nghost - The number of ghosted layers needed in the partitioned mesh
824b8ecf6d3SVijay Mahadevan . filename - The name of the mesh file to be loaded
825a2b725a8SWilliam Gropp - usrreadopts - The options string to read a MOAB mesh.
826a2b725a8SWilliam Gropp
827aa859218SVijay Mahadevan Output Parameter:
8281d27aa22SBarry Smith . dm - The `DM` object
829aa859218SVijay Mahadevan
830aa859218SVijay Mahadevan Level: beginner
831aa859218SVijay Mahadevan
832db781477SPatrick Sanan .seealso: `DMSetType()`, `DMCreate()`, `DMMoabCreateBoxMesh()`
833aa859218SVijay Mahadevan @*/
DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,PetscInt nghost,const char * filename,const char * usrreadopts,DM * dm)834d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm, PetscInt dim, PetscInt nghost, const char *filename, const char *usrreadopts, DM *dm)
835d71ae5a4SJacob Faibussowitsch {
836a4717116SVijay Mahadevan moab::ErrorCode merr;
8372e4e7c01SVijay Mahadevan PetscInt nprocs;
838a4717116SVijay Mahadevan DM_Moab *dmmoab;
839a4717116SVijay Mahadevan moab::Interface *mbiface;
8409daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
841a4717116SVijay Mahadevan moab::ParallelComm *pcomm;
8429daf19fdSVijay Mahadevan #endif
843a4717116SVijay Mahadevan moab::Range verts, elems;
8447023aa44SVijay Mahadevan const char *readopts;
84551d15aeeSVijay Mahadevan
84651d15aeeSVijay Mahadevan PetscFunctionBegin;
8474f572ea9SToby Isaac PetscAssertPointer(dm, 6);
84851d15aeeSVijay Mahadevan
849a4717116SVijay Mahadevan /* Create the basic DMMoab object and keep the default parameters created by DM impls */
8509566063dSJacob Faibussowitsch PetscCall(DMMoabCreateMoab(comm, NULL, NULL, NULL, dm));
85151d15aeeSVijay Mahadevan
852a4717116SVijay Mahadevan /* get all the necessary handles from the private DM object */
853a4717116SVijay Mahadevan dmmoab = (DM_Moab *)(*dm)->data;
854a4717116SVijay Mahadevan mbiface = dmmoab->mbiface;
8559daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
856a4717116SVijay Mahadevan pcomm = dmmoab->pcomm;
857a4717116SVijay Mahadevan nprocs = pcomm->size();
8589daf19fdSVijay Mahadevan #else
8599daf19fdSVijay Mahadevan nprocs = 1;
8609daf19fdSVijay Mahadevan #endif
861aa859218SVijay Mahadevan /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */
862c8d5365dSVijay Mahadevan dmmoab->dim = dim;
86349d66b22SVijay Mahadevan dmmoab->nghostrings = nghost;
864304006b3SVijay Mahadevan dmmoab->refct = 1;
865a4717116SVijay Mahadevan
866b5410836SVijay Mahadevan /* create a file set to associate all entities in current mesh */
8679371c9d4SSatish Balay merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);
8689371c9d4SSatish Balay MBERR("Creating file set failed", merr);
869b5410836SVijay Mahadevan
870a4717116SVijay Mahadevan /* add mesh loading options specific to the DM */
8719371c9d4SSatish Balay PetscCall(DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, nghost, dmmoab->read_mode, dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts));
8727023aa44SVijay Mahadevan
8739566063dSJacob Faibussowitsch PetscCall(PetscInfo(*dm, "Reading file %s with options: %s\n", filename, readopts));
874a4717116SVijay Mahadevan
875a4717116SVijay Mahadevan /* Load the mesh from a file. */
87649d66b22SVijay Mahadevan if (dmmoab->fileset) {
8779371c9d4SSatish Balay merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);
8789371c9d4SSatish Balay MBERRVM(mbiface, "Reading MOAB file failed.", merr);
8799371c9d4SSatish Balay } else {
8809371c9d4SSatish Balay merr = mbiface->load_file(filename, 0, readopts);
8819371c9d4SSatish Balay MBERRVM(mbiface, "Reading MOAB file failed.", merr);
88249d66b22SVijay Mahadevan }
883a4717116SVijay Mahadevan
8849daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
8856e40195eSVijay Mahadevan /* Reassign global IDs on all entities. */
8866ea892caSVijay Mahadevan /* merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, true, true, true);MBERRNM(merr); */
8879daf19fdSVijay Mahadevan #endif
888e5136372SVijay Mahadevan
889e5136372SVijay Mahadevan /* load the local vertices */
8909371c9d4SSatish Balay merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);
8919371c9d4SSatish Balay MBERRNM(merr);
892e5136372SVijay Mahadevan /* load the local elements */
8939371c9d4SSatish Balay merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);
8949371c9d4SSatish Balay MBERRNM(merr);
895e5136372SVijay Mahadevan
8969daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
897e5136372SVijay Mahadevan /* Everything is set up, now just do a tag exchange to update tags
898e5136372SVijay Mahadevan on all of the ghost vertexes */
8999371c9d4SSatish Balay merr = pcomm->exchange_tags(dmmoab->ltog_tag, verts);
9009371c9d4SSatish Balay MBERRV(mbiface, merr);
9019371c9d4SSatish Balay merr = pcomm->exchange_tags(dmmoab->ltog_tag, elems);
9029371c9d4SSatish Balay MBERRV(mbiface, merr);
9039371c9d4SSatish Balay merr = pcomm->collective_sync_partition();
9049371c9d4SSatish Balay MBERR("Collective sync failed", merr);
9059daf19fdSVijay Mahadevan #endif
906a4717116SVijay Mahadevan
90763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(*dm, "MOAB file '%s' was successfully loaded. Found %zu vertices and %zu elements.\n", filename, verts.size(), elems.size()));
9089566063dSJacob Faibussowitsch PetscCall(PetscFree(readopts));
9093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
91051d15aeeSVijay Mahadevan }
91151d15aeeSVijay Mahadevan
912cab5ea25SPierre Jolivet /*@C
913304006b3SVijay Mahadevan DMMoabRenumberMeshEntities - Order and number all entities (vertices->elements) to be contiguously ordered
914304006b3SVijay Mahadevan in parallel
915304006b3SVijay Mahadevan
916d083f849SBarry Smith Collective
917304006b3SVijay Mahadevan
918304006b3SVijay Mahadevan Input Parameters:
919a2b725a8SWilliam Gropp . dm - The DM object
920304006b3SVijay Mahadevan
921304006b3SVijay Mahadevan Level: advanced
922304006b3SVijay Mahadevan
923db781477SPatrick Sanan .seealso: `DMSetUp()`, `DMCreate()`
924304006b3SVijay Mahadevan @*/
DMMoabRenumberMeshEntities(DM dm)925d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabRenumberMeshEntities(DM dm)
926d71ae5a4SJacob Faibussowitsch {
927304006b3SVijay Mahadevan moab::Range verts;
928304006b3SVijay Mahadevan
929304006b3SVijay Mahadevan PetscFunctionBegin;
930304006b3SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
931304006b3SVijay Mahadevan
9329daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
933304006b3SVijay Mahadevan /* Insert new points */
9349daf19fdSVijay Mahadevan moab::ErrorCode merr;
9359371c9d4SSatish Balay merr = ((DM_Moab *)dm->data)->pcomm->assign_global_ids(((DM_Moab *)dm->data)->fileset, 3, 0, false, true, false);
9369371c9d4SSatish Balay MBERRNM(merr);
9379daf19fdSVijay Mahadevan #endif
9383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
939304006b3SVijay Mahadevan }
940