xref: /petsc/src/dm/impls/moab/dmmbutil.cxx (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
279371c9d4SSatish Balay PetscInt DMMoab_SetTensorElementConnectivity_Private(DMMoabMeshGeneratorCtx &genCtx, PetscInt offset, PetscInt corner, std::vector<PetscInt> &subent_conn, moab::EntityHandle *connectivity) {
28ce27a4eeSVijay Mahadevan   switch (genCtx.dim) {
29e5136372SVijay Mahadevan   case 1:
30ce27a4eeSVijay Mahadevan     subent_conn.resize(2);
31ce27a4eeSVijay Mahadevan     moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data());
32ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[0]] = corner;
33ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[1]] = corner + 1;
34e5136372SVijay Mahadevan     break;
35e5136372SVijay Mahadevan   case 2:
36ce27a4eeSVijay Mahadevan     subent_conn.resize(4);
37ce27a4eeSVijay Mahadevan     moab::CN::SubEntityVertexIndices(moab::MBQUAD, 2, 0, subent_conn.data());
38ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[0]] = corner;
39ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[1]] = corner + 1;
40ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[2]] = corner + 1 + genCtx.ystride;
41ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[3]] = corner + genCtx.ystride;
42e5136372SVijay Mahadevan     break;
43e5136372SVijay Mahadevan   case 3:
44e5136372SVijay Mahadevan   default:
45ce27a4eeSVijay Mahadevan     subent_conn.resize(8);
46ce27a4eeSVijay Mahadevan     moab::CN::SubEntityVertexIndices(moab::MBHEX, 3, 0, subent_conn.data());
47ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[0]] = corner;
48ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[1]] = corner + 1;
49ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[2]] = corner + 1 + genCtx.ystride;
50ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[3]] = corner + genCtx.ystride;
51ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[4]] = corner + genCtx.zstride;
52ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[5]] = corner + 1 + genCtx.zstride;
53ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[6]] = corner + 1 + genCtx.ystride + genCtx.zstride;
54ce27a4eeSVijay Mahadevan     connectivity[offset + subent_conn[7]] = corner + genCtx.ystride + genCtx.zstride;
55e5136372SVijay Mahadevan     break;
56e5136372SVijay Mahadevan   }
57ce27a4eeSVijay Mahadevan   return subent_conn.size();
58e5136372SVijay Mahadevan }
5951d15aeeSVijay Mahadevan 
609371c9d4SSatish Balay PetscInt DMMoab_SetSimplexElementConnectivity_Private(DMMoabMeshGeneratorCtx &genCtx, PetscInt subelem, PetscInt offset, PetscInt corner, std::vector<PetscInt> &subent_conn, moab::EntityHandle *connectivity) {
61755f3dfbSVijay Mahadevan   PetscInt       A, B, C, D, E, F, G, H, M;
62755f3dfbSVijay Mahadevan   const PetscInt trigen_opts = 1; /* 1 - Aligned diagonally to right, 2 - Aligned diagonally to left, 3 - 4 elements per quad */
6349d66b22SVijay Mahadevan   A                          = corner;
6449d66b22SVijay Mahadevan   B                          = corner + 1;
65ce27a4eeSVijay Mahadevan   switch (genCtx.dim) {
66c3dd00c7SVijay Mahadevan   case 1:
67ce27a4eeSVijay Mahadevan     subent_conn.resize(2); /* only linear EDGE supported now */
68ce27a4eeSVijay Mahadevan     moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data());
6949d66b22SVijay Mahadevan     connectivity[offset + subent_conn[0]] = A;
7049d66b22SVijay Mahadevan     connectivity[offset + subent_conn[1]] = B;
71c3dd00c7SVijay Mahadevan     break;
72c3dd00c7SVijay Mahadevan   case 2:
7349d66b22SVijay Mahadevan     C = corner + 1 + genCtx.ystride;
7449d66b22SVijay Mahadevan     D = corner + genCtx.ystride;
7563d025dbSVijay Mahadevan     M = corner + 0.5;      /* technically -- need to modify vertex generation */
76ce27a4eeSVijay Mahadevan     subent_conn.resize(3); /* only linear TRI supported */
77ce27a4eeSVijay Mahadevan     moab::CN::SubEntityVertexIndices(moab::MBTRI, 2, 0, subent_conn.data());
78755f3dfbSVijay Mahadevan     if (trigen_opts == 1) {
79ce27a4eeSVijay Mahadevan       if (subelem) { /* 0 1 2 of a QUAD */
8077a8c067SVijay Mahadevan         connectivity[offset + subent_conn[0]] = B;
81755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[1]] = C;
8277a8c067SVijay Mahadevan         connectivity[offset + subent_conn[2]] = A;
831baa6e33SBarry Smith       } else { /* 2 3 0 of a QUAD */
8477a8c067SVijay Mahadevan         connectivity[offset + subent_conn[0]] = D;
8577a8c067SVijay Mahadevan         connectivity[offset + subent_conn[1]] = A;
86755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[2]] = C;
87755f3dfbSVijay Mahadevan       }
881baa6e33SBarry Smith     } else if (trigen_opts == 2) {
89755f3dfbSVijay Mahadevan       if (subelem) { /* 0 1 2 of a QUAD */
90755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[0]] = A;
9177a8c067SVijay Mahadevan         connectivity[offset + subent_conn[1]] = B;
9277a8c067SVijay Mahadevan         connectivity[offset + subent_conn[2]] = D;
931baa6e33SBarry Smith       } else { /* 2 3 0 of a QUAD */
9477a8c067SVijay Mahadevan         connectivity[offset + subent_conn[0]] = C;
9577a8c067SVijay Mahadevan         connectivity[offset + subent_conn[1]] = D;
96755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[2]] = B;
97755f3dfbSVijay Mahadevan       }
981baa6e33SBarry Smith     } else {
99755f3dfbSVijay Mahadevan       switch (subelem) { /* 0 1 2 of a QUAD */
100755f3dfbSVijay Mahadevan       case 0:
101755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[0]] = A;
102755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[1]] = B;
103755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[2]] = M;
104755f3dfbSVijay Mahadevan         break;
105755f3dfbSVijay Mahadevan       case 1:
106755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[0]] = B;
107755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[1]] = C;
108755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[2]] = M;
109755f3dfbSVijay Mahadevan         break;
110755f3dfbSVijay Mahadevan       case 2:
11149d66b22SVijay Mahadevan         connectivity[offset + subent_conn[0]] = C;
11249d66b22SVijay Mahadevan         connectivity[offset + subent_conn[1]] = D;
113755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[2]] = M;
114755f3dfbSVijay Mahadevan         break;
115755f3dfbSVijay Mahadevan       case 3:
116755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[0]] = D;
117755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[1]] = A;
118755f3dfbSVijay Mahadevan         connectivity[offset + subent_conn[2]] = M;
119755f3dfbSVijay Mahadevan         break;
120755f3dfbSVijay Mahadevan       }
121c3dd00c7SVijay Mahadevan     }
122c3dd00c7SVijay Mahadevan     break;
123c3dd00c7SVijay Mahadevan   case 3:
124c3dd00c7SVijay Mahadevan   default:
12549d66b22SVijay Mahadevan     C = corner + 1 + genCtx.ystride;
12649d66b22SVijay Mahadevan     D = corner + genCtx.ystride;
12749d66b22SVijay Mahadevan     E = corner + genCtx.zstride;
12849d66b22SVijay Mahadevan     F = corner + 1 + genCtx.zstride;
12949d66b22SVijay Mahadevan     G = corner + 1 + genCtx.ystride + genCtx.zstride;
13049d66b22SVijay Mahadevan     H = corner + genCtx.ystride + genCtx.zstride;
131ce27a4eeSVijay Mahadevan     subent_conn.resize(4); /* only linear TET supported */
132ce27a4eeSVijay Mahadevan     moab::CN::SubEntityVertexIndices(moab::MBTET, 3, 0, subent_conn.data());
133c3dd00c7SVijay Mahadevan     switch (subelem) {
13449d66b22SVijay Mahadevan     case 0: /* 4 3 7 6 of a HEX */
13549d66b22SVijay Mahadevan       connectivity[offset + subent_conn[0]] = E;
13649d66b22SVijay Mahadevan       connectivity[offset + subent_conn[1]] = D;
13749d66b22SVijay Mahadevan       connectivity[offset + subent_conn[2]] = H;
13849d66b22SVijay Mahadevan       connectivity[offset + subent_conn[3]] = G;
139c3dd00c7SVijay Mahadevan       break;
14049d66b22SVijay Mahadevan     case 1: /* 0 1 2 5 of a HEX */
14149d66b22SVijay Mahadevan       connectivity[offset + subent_conn[0]] = A;
14249d66b22SVijay Mahadevan       connectivity[offset + subent_conn[1]] = B;
14349d66b22SVijay Mahadevan       connectivity[offset + subent_conn[2]] = C;
14449d66b22SVijay Mahadevan       connectivity[offset + subent_conn[3]] = F;
145c3dd00c7SVijay Mahadevan       break;
14649d66b22SVijay Mahadevan     case 2: /* 0 3 4 5 of a HEX */
14749d66b22SVijay Mahadevan       connectivity[offset + subent_conn[0]] = A;
14849d66b22SVijay Mahadevan       connectivity[offset + subent_conn[1]] = D;
14949d66b22SVijay Mahadevan       connectivity[offset + subent_conn[2]] = E;
15049d66b22SVijay Mahadevan       connectivity[offset + subent_conn[3]] = F;
151c3dd00c7SVijay Mahadevan       break;
15249d66b22SVijay Mahadevan     case 3: /* 2 6 3 5 of a HEX */
15349d66b22SVijay Mahadevan       connectivity[offset + subent_conn[0]] = C;
15449d66b22SVijay Mahadevan       connectivity[offset + subent_conn[1]] = G;
15549d66b22SVijay Mahadevan       connectivity[offset + subent_conn[2]] = D;
15649d66b22SVijay Mahadevan       connectivity[offset + subent_conn[3]] = F;
157c3dd00c7SVijay Mahadevan       break;
15849d66b22SVijay Mahadevan     case 4: /* 0 2 3 5 of a HEX */
15949d66b22SVijay Mahadevan       connectivity[offset + subent_conn[0]] = A;
16049d66b22SVijay Mahadevan       connectivity[offset + subent_conn[1]] = C;
16149d66b22SVijay Mahadevan       connectivity[offset + subent_conn[2]] = D;
16249d66b22SVijay Mahadevan       connectivity[offset + subent_conn[3]] = F;
16349d66b22SVijay Mahadevan       break;
16449d66b22SVijay Mahadevan     case 5: /* 3 6 4 5 of a HEX */
16549d66b22SVijay Mahadevan       connectivity[offset + subent_conn[0]] = D;
16649d66b22SVijay Mahadevan       connectivity[offset + subent_conn[1]] = G;
16749d66b22SVijay Mahadevan       connectivity[offset + subent_conn[2]] = E;
16849d66b22SVijay Mahadevan       connectivity[offset + subent_conn[3]] = F;
169c3dd00c7SVijay Mahadevan       break;
170c3dd00c7SVijay Mahadevan     }
171c3dd00c7SVijay Mahadevan     break;
172c3dd00c7SVijay Mahadevan   }
173ce27a4eeSVijay Mahadevan   return subent_conn.size();
174c3dd00c7SVijay Mahadevan }
175c3dd00c7SVijay Mahadevan 
1769371c9d4SSatish Balay std::pair<PetscInt, PetscInt> DMMoab_SetElementConnectivity_Private(DMMoabMeshGeneratorCtx &genCtx, PetscInt offset, PetscInt corner, moab::EntityHandle *connectivity) {
17749d66b22SVijay Mahadevan   PetscInt              vcount                  = 0;
17849d66b22SVijay Mahadevan   PetscInt              simplices_per_tensor[4] = {0, 1, 2, 6};
17949d66b22SVijay Mahadevan   std::vector<PetscInt> subent_conn; /* only linear edge, tri, tet supported now */
180ce27a4eeSVijay Mahadevan   subent_conn.reserve(27);
18149d66b22SVijay Mahadevan   PetscInt m, subelem;
182ce27a4eeSVijay Mahadevan   if (genCtx.simplex) {
18349d66b22SVijay Mahadevan     subelem = simplices_per_tensor[genCtx.dim];
184ce27a4eeSVijay Mahadevan     for (m = 0; m < subelem; m++) {
185ce27a4eeSVijay Mahadevan       vcount = DMMoab_SetSimplexElementConnectivity_Private(genCtx, m, offset, corner, subent_conn, connectivity);
186ce27a4eeSVijay Mahadevan       offset += vcount;
187ce27a4eeSVijay Mahadevan     }
1881baa6e33SBarry Smith   } else {
189c3dd00c7SVijay Mahadevan     subelem = 1;
190ce27a4eeSVijay Mahadevan     vcount  = DMMoab_SetTensorElementConnectivity_Private(genCtx, offset, corner, subent_conn, connectivity);
191c3dd00c7SVijay Mahadevan   }
19249d66b22SVijay Mahadevan   return std::pair<PetscInt, PetscInt>(vcount * subelem, subelem);
193c3dd00c7SVijay Mahadevan }
194c3dd00c7SVijay Mahadevan 
1959371c9d4SSatish Balay 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) {
19649d66b22SVijay Mahadevan   PetscInt                 x, y, z, ix, nnodes;
19749d66b22SVijay Mahadevan   PetscInt                 ii, jj, kk;
19849d66b22SVijay Mahadevan   std::vector<PetscReal *> arrays;
199755f3dfbSVijay Mahadevan   PetscInt                *gids;
200ce27a4eeSVijay Mahadevan   moab::ErrorCode          merr;
201ce27a4eeSVijay Mahadevan 
202ce27a4eeSVijay Mahadevan   PetscFunctionBegin;
203ce27a4eeSVijay Mahadevan   /* we will generate (q*block+1)^3 vertices, and block^3 hexas; q is 1 for linear, 2 for quadratic
204ce27a4eeSVijay Mahadevan    * the global id of the vertices will come from m, n, k, a, b, c
205ce27a4eeSVijay Mahadevan    * x will vary from  m*A*q*block + a*q*block to m*A*q*block+(a+1)*q*block etc.
206ce27a4eeSVijay Mahadevan    */
207ce27a4eeSVijay Mahadevan   nnodes = genCtx.blockSizeVertexXYZ[0] * (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] * (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1) : 1);
2089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nnodes, &gids));
209ce27a4eeSVijay Mahadevan 
2109371c9d4SSatish Balay   merr = iface->get_node_coords(3, nnodes, 0, startv, arrays);
2119371c9d4SSatish Balay   MBERR("Can't get node coords.", merr);
212ce27a4eeSVijay Mahadevan 
213ce27a4eeSVijay Mahadevan   /* will start with the lower corner: */
21463d025dbSVijay Mahadevan   /* x = ( m * genCtx.A + a) * genCtx.q * genCtx.blockSizeElementXYZ[0]; */
21563d025dbSVijay Mahadevan   /* y = ( n * genCtx.B + b) * genCtx.q * genCtx.blockSizeElementXYZ[1]; */
21663d025dbSVijay Mahadevan   /* z = ( k * genCtx.C + c) * genCtx.q * genCtx.blockSizeElementXYZ[2]; */
217755f3dfbSVijay Mahadevan 
218755f3dfbSVijay Mahadevan   x = (m * genCtx.A + a) * genCtx.q;
219755f3dfbSVijay Mahadevan   y = (n * genCtx.B + b) * genCtx.q;
220755f3dfbSVijay Mahadevan   z = (k * genCtx.C + c) * genCtx.q;
2219566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Starting offset for coordinates := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", x, y, z));
222ce27a4eeSVijay Mahadevan   ix = 0;
223ce27a4eeSVijay Mahadevan   moab::Range verts(startv, startv + nnodes - 1);
224ce27a4eeSVijay Mahadevan   for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeVertexXYZ[2] : 1); kk++) {
225ce27a4eeSVijay Mahadevan     for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeVertexXYZ[1] : 1); jj++) {
226ce27a4eeSVijay Mahadevan       for (ii = 0; ii < genCtx.blockSizeVertexXYZ[0]; ii++, ix++) {
227ce27a4eeSVijay Mahadevan         /* set coordinates for the vertices */
228ce27a4eeSVijay Mahadevan         arrays[0][ix] = (x + ii) * genCtx.dx + genCtx.xyzbounds[0];
229ce27a4eeSVijay Mahadevan         arrays[1][ix] = (y + jj) * genCtx.dy + genCtx.xyzbounds[2];
230ce27a4eeSVijay Mahadevan         arrays[2][ix] = (z + kk) * genCtx.dz + genCtx.xyzbounds[4];
2319566063dSJacob Faibussowitsch         PetscCall(PetscInfo(NULL, "Creating vertex with coordinates := %f, %f, %f\n", arrays[0][ix], arrays[1][ix], arrays[2][ix]));
232ce27a4eeSVijay Mahadevan 
233ce27a4eeSVijay Mahadevan         /* If we want to set some tags on the vertices -> use the following entity handle definition:
234ce27a4eeSVijay Mahadevan            moab::EntityHandle v = startv + ix;
235ce27a4eeSVijay Mahadevan         */
236ce27a4eeSVijay Mahadevan         /* compute the global ID for vertex */
237ce27a4eeSVijay Mahadevan         gids[ix] = 1 + (x + ii) + (y + jj) * genCtx.NX + (z + kk) * (genCtx.NX * genCtx.NY);
238ce27a4eeSVijay Mahadevan       }
239ce27a4eeSVijay Mahadevan     }
240ce27a4eeSVijay Mahadevan   }
241ce27a4eeSVijay Mahadevan   /* set global ID data on vertices */
242ce27a4eeSVijay Mahadevan   mbImpl->tag_set_data(global_id_tag, verts, &gids[0]);
243ce27a4eeSVijay Mahadevan   verts.swap(uverts);
2449566063dSJacob Faibussowitsch   PetscCall(PetscFree(gids));
245ce27a4eeSVijay Mahadevan   PetscFunctionReturn(0);
246ce27a4eeSVijay Mahadevan }
247ce27a4eeSVijay Mahadevan 
2489371c9d4SSatish Balay 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) {
249ce27a4eeSVijay Mahadevan   moab::ErrorCode     merr;
250ce27a4eeSVijay Mahadevan   PetscInt            ix, ie, xe, ye, ze;
251ce27a4eeSVijay Mahadevan   PetscInt            ii, jj, kk, nvperelem;
25249d66b22SVijay Mahadevan   PetscInt            simplices_per_tensor[4] = {0, 1, 2, 6};
253ce27a4eeSVijay 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);*/
254ce27a4eeSVijay Mahadevan   PetscInt            nelems                  = ntensorelems;
25563d025dbSVijay Mahadevan   moab::EntityHandle  starte; /* connectivity */
256ce27a4eeSVijay Mahadevan   moab::EntityHandle *conn;
257ce27a4eeSVijay Mahadevan 
258ce27a4eeSVijay Mahadevan   PetscFunctionBegin;
259ce27a4eeSVijay Mahadevan   switch (genCtx.dim) {
260ce27a4eeSVijay Mahadevan   case 1:
261ce27a4eeSVijay Mahadevan     nvperelem = 2;
2629371c9d4SSatish Balay     merr      = iface->get_element_connect(nelems, 2, moab::MBEDGE, 0, starte, conn);
2639371c9d4SSatish Balay     MBERR("Can't get EDGE2 element connectivity.", merr);
264ce27a4eeSVijay Mahadevan     break;
265ce27a4eeSVijay Mahadevan   case 2:
266ce27a4eeSVijay Mahadevan     if (genCtx.simplex) {
267ce27a4eeSVijay Mahadevan       nvperelem = 3;
26849d66b22SVijay Mahadevan       nelems    = ntensorelems * simplices_per_tensor[genCtx.dim];
2699371c9d4SSatish Balay       merr      = iface->get_element_connect(nelems, 3, moab::MBTRI, 0, starte, conn);
2709371c9d4SSatish Balay       MBERR("Can't get TRI3 element connectivity.", merr);
2711baa6e33SBarry Smith     } else {
272ce27a4eeSVijay Mahadevan       nvperelem = 4;
2739371c9d4SSatish Balay       merr      = iface->get_element_connect(nelems, 4, moab::MBQUAD, 0, starte, conn);
2749371c9d4SSatish Balay       MBERR("Can't get QUAD4 element connectivity.", merr);
275ce27a4eeSVijay Mahadevan     }
276ce27a4eeSVijay Mahadevan     break;
277ce27a4eeSVijay Mahadevan   case 3:
278ce27a4eeSVijay Mahadevan   default:
279ce27a4eeSVijay Mahadevan     if (genCtx.simplex) {
280ce27a4eeSVijay Mahadevan       nvperelem = 4;
28149d66b22SVijay Mahadevan       nelems    = ntensorelems * simplices_per_tensor[genCtx.dim];
2829371c9d4SSatish Balay       merr      = iface->get_element_connect(nelems, 4, moab::MBTET, 0, starte, conn);
2839371c9d4SSatish Balay       MBERR("Can't get TET4 element connectivity.", merr);
2841baa6e33SBarry Smith     } else {
285ce27a4eeSVijay Mahadevan       nvperelem = 8;
2869371c9d4SSatish Balay       merr      = iface->get_element_connect(nelems, 8, moab::MBHEX, 0, starte, conn);
2879371c9d4SSatish Balay       MBERR("Can't get HEX8 element connectivity.", merr);
288ce27a4eeSVijay Mahadevan     }
289ce27a4eeSVijay Mahadevan     break;
290ce27a4eeSVijay Mahadevan   }
291ce27a4eeSVijay Mahadevan 
29263d025dbSVijay Mahadevan   ix = ie = 0; /* index now in the elements, for global ids */
293ce27a4eeSVijay Mahadevan 
294ce27a4eeSVijay Mahadevan   /* create a temporary range to store local element handles */
295ce27a4eeSVijay Mahadevan   moab::Range           tmp(starte, starte + nelems - 1);
29649d66b22SVijay Mahadevan   std::vector<PetscInt> gids(nelems);
297ce27a4eeSVijay Mahadevan 
298ce27a4eeSVijay Mahadevan   /* identify the elements at the lower corner, for their global ids */
299ce27a4eeSVijay Mahadevan   xe = m * genCtx.A * genCtx.blockSizeElementXYZ[0] + a * genCtx.blockSizeElementXYZ[0];
300ce27a4eeSVijay Mahadevan   ye = (genCtx.dim > 1 ? n * genCtx.B * genCtx.blockSizeElementXYZ[1] + b * genCtx.blockSizeElementXYZ[1] : 0);
301ce27a4eeSVijay Mahadevan   ze = (genCtx.dim > 2 ? k * genCtx.C * genCtx.blockSizeElementXYZ[2] + c * genCtx.blockSizeElementXYZ[2] : 0);
302ce27a4eeSVijay Mahadevan 
303ce27a4eeSVijay Mahadevan   /* create owned elements requested by genCtx */
304ce27a4eeSVijay Mahadevan   for (kk = 0; kk < (genCtx.dim > 2 ? genCtx.blockSizeElementXYZ[2] : 1); kk++) {
305ce27a4eeSVijay Mahadevan     for (jj = 0; jj < (genCtx.dim > 1 ? genCtx.blockSizeElementXYZ[1] : 1); jj++) {
306ce27a4eeSVijay Mahadevan       for (ii = 0; ii < genCtx.blockSizeElementXYZ[0]; ii++) {
307ce27a4eeSVijay Mahadevan         moab::EntityHandle corner = startv + genCtx.q * ii + genCtx.q * jj * genCtx.ystride + genCtx.q * kk * genCtx.zstride;
308ce27a4eeSVijay Mahadevan 
30949d66b22SVijay Mahadevan         std::pair<PetscInt, PetscInt> entoffset = DMMoab_SetElementConnectivity_Private(genCtx, ix, corner, conn);
310ce27a4eeSVijay Mahadevan 
31149d66b22SVijay Mahadevan         for (PetscInt j = 0; j < entoffset.second; j++) {
312ce27a4eeSVijay Mahadevan           /* The entity handle for the particular element -> if we want to set some tags is
313ce27a4eeSVijay Mahadevan              moab::EntityHandle eh = starte + ie + j;
314ce27a4eeSVijay Mahadevan           */
315ce27a4eeSVijay Mahadevan           gids[ie + j] = 1 + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney));
31663d025dbSVijay Mahadevan           /* gids[ie+j] = ie + j + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)); */
31763d025dbSVijay Mahadevan           /* gids[ie+j] = 1 + ie; */
31863d025dbSVijay Mahadevan           /* ie++; */
319ce27a4eeSVijay Mahadevan         }
320ce27a4eeSVijay Mahadevan 
321ce27a4eeSVijay Mahadevan         ix += entoffset.first;
32249d66b22SVijay Mahadevan         ie += entoffset.second;
323ce27a4eeSVijay Mahadevan       }
324ce27a4eeSVijay Mahadevan     }
325ce27a4eeSVijay Mahadevan   }
326ce27a4eeSVijay Mahadevan   if (genCtx.adjEnts) { /* we need to update adjacencies now, because some elements are new */
3279371c9d4SSatish Balay     merr = iface->update_adjacencies(starte, nelems, nvperelem, conn);
3289371c9d4SSatish Balay     MBERR("Can't update adjacencies", merr);
329ce27a4eeSVijay Mahadevan   }
330ce27a4eeSVijay Mahadevan   tmp.swap(cells);
3319371c9d4SSatish Balay   merr = mbImpl->tag_set_data(global_id_tag, cells, &gids[0]);
3329371c9d4SSatish Balay   MBERR("Can't set global ids to elements.", merr);
333ce27a4eeSVijay Mahadevan   PetscFunctionReturn(0);
334ce27a4eeSVijay Mahadevan }
335ce27a4eeSVijay Mahadevan 
3369371c9d4SSatish Balay PetscErrorCode DMMBUtil_InitializeOptions(DMMoabMeshGeneratorCtx &genCtx, PetscInt dim, PetscBool simplex, PetscInt rank, PetscInt nprocs, const PetscReal *bounds, PetscInt nelems) {
337ce27a4eeSVijay Mahadevan   PetscFunctionBegin;
338ce27a4eeSVijay Mahadevan   /* Initialize all genCtx data */
339ce27a4eeSVijay Mahadevan   genCtx.dim            = dim;
340ce27a4eeSVijay Mahadevan   genCtx.simplex        = simplex;
34149d66b22SVijay Mahadevan   genCtx.newMergeMethod = genCtx.keep_skins = genCtx.adjEnts = true;
342ce27a4eeSVijay Mahadevan   /* determine other global quantities for the mesh used for nodes increments */
343ce27a4eeSVijay Mahadevan   genCtx.q                                                   = 1;
344755f3dfbSVijay Mahadevan   genCtx.fraction = genCtx.remainder = genCtx.cumfraction = 0;
345ce27a4eeSVijay Mahadevan 
346ce27a4eeSVijay Mahadevan   if (!genCtx.usrxyzgrid) { /* not overridden by genCtx - assume nele equally and that genCtx wants a uniform cube mesh */
347755f3dfbSVijay Mahadevan 
348755f3dfbSVijay Mahadevan     genCtx.fraction    = nelems / nprocs; /* partition only by the largest dimension */
349755f3dfbSVijay Mahadevan     genCtx.remainder   = nelems % nprocs; /* remainder after partition which gets evenly distributed by round-robin */
350755f3dfbSVijay Mahadevan     genCtx.cumfraction = (rank > 0 ? (genCtx.fraction) * (rank) + (rank - 1 < genCtx.remainder ? rank : genCtx.remainder) : 0);
351755f3dfbSVijay Mahadevan     if (rank < genCtx.remainder) /* This process gets "fraction+1" elements */
352755f3dfbSVijay Mahadevan       genCtx.fraction++;
353755f3dfbSVijay Mahadevan 
3549566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "Fraction = %" PetscInt_FMT ", Remainder = %" PetscInt_FMT ", Cumulative fraction = %" PetscInt_FMT "\n", genCtx.fraction, genCtx.remainder, genCtx.cumfraction));
355755f3dfbSVijay Mahadevan     switch (genCtx.dim) {
356755f3dfbSVijay Mahadevan     case 1:
357755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[0] = genCtx.fraction;
358755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[1] = 1;
359755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[2] = 1;
360755f3dfbSVijay Mahadevan       break;
361755f3dfbSVijay Mahadevan     case 2:
362755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[0] = nelems;
363755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[1] = genCtx.fraction;
364755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[2] = 1;
365755f3dfbSVijay Mahadevan       break;
366755f3dfbSVijay Mahadevan     case 3:
367755f3dfbSVijay Mahadevan     default:
368755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[0] = nelems;
369755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[1] = nelems;
370755f3dfbSVijay Mahadevan       genCtx.blockSizeElementXYZ[2] = genCtx.fraction;
371755f3dfbSVijay Mahadevan       break;
37249d66b22SVijay Mahadevan     }
373ce27a4eeSVijay Mahadevan   }
374ce27a4eeSVijay Mahadevan 
375ce27a4eeSVijay Mahadevan   /* partition only by the largest dimension */
376ce27a4eeSVijay Mahadevan   /* Total number of local elements := genCtx.blockSizeElementXYZ[0]*(genCtx.dim>1? genCtx.blockSizeElementXYZ[1]*(genCtx.dim>2 ? genCtx.blockSizeElementXYZ[2]:1) :1); */
377ce27a4eeSVijay Mahadevan   if (bounds) {
3785f80ce2aSJacob Faibussowitsch     for (PetscInt i = 0; i < 6; i++) genCtx.xyzbounds[i] = bounds[i];
3791baa6e33SBarry Smith   } else {
380ce27a4eeSVijay Mahadevan     genCtx.xyzbounds[0] = genCtx.xyzbounds[2] = genCtx.xyzbounds[4] = 0.0;
381ce27a4eeSVijay Mahadevan     genCtx.xyzbounds[1] = genCtx.xyzbounds[3] = genCtx.xyzbounds[5] = 1.0;
382ce27a4eeSVijay Mahadevan   }
383ce27a4eeSVijay Mahadevan 
384ce27a4eeSVijay Mahadevan   if (!genCtx.usrprocgrid) {
385ce27a4eeSVijay Mahadevan     switch (genCtx.dim) {
386ce27a4eeSVijay Mahadevan     case 1:
387ce27a4eeSVijay Mahadevan       genCtx.M = nprocs;
388ce27a4eeSVijay Mahadevan       genCtx.N = genCtx.K = 1;
389ce27a4eeSVijay Mahadevan       break;
390ce27a4eeSVijay Mahadevan     case 2:
391755f3dfbSVijay Mahadevan       genCtx.N = nprocs;
39263d025dbSVijay Mahadevan       genCtx.M = genCtx.K = 1;
393ce27a4eeSVijay Mahadevan       break;
394ce27a4eeSVijay Mahadevan     default:
395755f3dfbSVijay Mahadevan       genCtx.K = nprocs;
39663d025dbSVijay Mahadevan       genCtx.M = genCtx.N = 1;
397ce27a4eeSVijay Mahadevan       break;
398ce27a4eeSVijay Mahadevan     }
399ce27a4eeSVijay Mahadevan   }
400ce27a4eeSVijay Mahadevan 
4019371c9d4SSatish Balay   if (!genCtx.usrrefgrid) { genCtx.A = genCtx.B = genCtx.C = 1; }
402ce27a4eeSVijay Mahadevan 
403ce27a4eeSVijay Mahadevan   /* more default values */
404ce27a4eeSVijay Mahadevan   genCtx.nex = genCtx.ney = genCtx.nez = 0;
405ce27a4eeSVijay Mahadevan   genCtx.xstride = genCtx.ystride = genCtx.zstride = 0;
406ce27a4eeSVijay Mahadevan   genCtx.NX = genCtx.NY = genCtx.NZ = 0;
407ce27a4eeSVijay Mahadevan   genCtx.nex = genCtx.ney = genCtx.nez = 0;
408ce27a4eeSVijay Mahadevan   genCtx.blockSizeVertexXYZ[0] = genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 1;
409ce27a4eeSVijay Mahadevan 
410ce27a4eeSVijay Mahadevan   switch (genCtx.dim) {
411ce27a4eeSVijay Mahadevan   case 3:
412ce27a4eeSVijay Mahadevan     genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1;
413ce27a4eeSVijay Mahadevan     genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1;
414ce27a4eeSVijay Mahadevan     genCtx.blockSizeVertexXYZ[2] = genCtx.q * genCtx.blockSizeElementXYZ[2] + 1;
415ce27a4eeSVijay Mahadevan 
416ce27a4eeSVijay Mahadevan     genCtx.nex     = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0];               /* number of elements in x direction, used for global id on element */
417755f3dfbSVijay Mahadevan     genCtx.dx      = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */
418ce27a4eeSVijay Mahadevan     genCtx.NX      = (genCtx.q * genCtx.nex + 1);
419ce27a4eeSVijay Mahadevan     genCtx.xstride = 1;
420ce27a4eeSVijay Mahadevan     genCtx.ney     = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1];               /* number of elements in y direction  .... */
421755f3dfbSVijay Mahadevan     genCtx.dy      = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */
422ce27a4eeSVijay Mahadevan     genCtx.NY      = (genCtx.q * genCtx.ney + 1);
423ce27a4eeSVijay Mahadevan     genCtx.ystride = genCtx.blockSizeVertexXYZ[0];
424ce27a4eeSVijay Mahadevan     genCtx.nez     = genCtx.K * genCtx.C * genCtx.blockSizeElementXYZ[2];               /* number of elements in z direction  .... */
425755f3dfbSVijay Mahadevan     genCtx.dz      = (genCtx.xyzbounds[5] - genCtx.xyzbounds[4]) / (nelems * genCtx.q); /* distance between 2 nodes in z direction */
426ce27a4eeSVijay Mahadevan     genCtx.NZ      = (genCtx.q * genCtx.nez + 1);
427ce27a4eeSVijay Mahadevan     genCtx.zstride = genCtx.blockSizeVertexXYZ[0] * genCtx.blockSizeVertexXYZ[1];
428ce27a4eeSVijay Mahadevan     break;
429ce27a4eeSVijay Mahadevan   case 2:
430ce27a4eeSVijay Mahadevan     genCtx.blockSizeVertexXYZ[0] = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1;
431ce27a4eeSVijay Mahadevan     genCtx.blockSizeVertexXYZ[1] = genCtx.q * genCtx.blockSizeElementXYZ[1] + 1;
4329daf19fdSVijay Mahadevan     genCtx.blockSizeVertexXYZ[2] = 0;
433ce27a4eeSVijay Mahadevan 
434ce27a4eeSVijay Mahadevan     genCtx.nex     = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0];                   /* number of elements in x direction, used for global id on element */
435ce27a4eeSVijay Mahadevan     genCtx.dx      = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (genCtx.nex * genCtx.q); /* distance between 2 nodes in x direction */
436ce27a4eeSVijay Mahadevan     genCtx.NX      = (genCtx.q * genCtx.nex + 1);
437ce27a4eeSVijay Mahadevan     genCtx.xstride = 1;
438ce27a4eeSVijay Mahadevan     genCtx.ney     = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1];               /* number of elements in y direction  .... */
439755f3dfbSVijay Mahadevan     genCtx.dy      = (genCtx.xyzbounds[3] - genCtx.xyzbounds[2]) / (nelems * genCtx.q); /* distance between 2 nodes in y direction */
440ce27a4eeSVijay Mahadevan     genCtx.NY      = (genCtx.q * genCtx.ney + 1);
441ce27a4eeSVijay Mahadevan     genCtx.ystride = genCtx.blockSizeVertexXYZ[0];
442ce27a4eeSVijay Mahadevan     break;
443ce27a4eeSVijay Mahadevan   case 1:
4449daf19fdSVijay Mahadevan     genCtx.blockSizeVertexXYZ[1] = genCtx.blockSizeVertexXYZ[2] = 0;
445ce27a4eeSVijay Mahadevan     genCtx.blockSizeVertexXYZ[0]                                = genCtx.q * genCtx.blockSizeElementXYZ[0] + 1;
446ce27a4eeSVijay Mahadevan 
447ce27a4eeSVijay Mahadevan     genCtx.nex     = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0];               /* number of elements in x direction, used for global id on element */
448755f3dfbSVijay Mahadevan     genCtx.dx      = (genCtx.xyzbounds[1] - genCtx.xyzbounds[0]) / (nelems * genCtx.q); /* distance between 2 nodes in x direction */
449ce27a4eeSVijay Mahadevan     genCtx.NX      = (genCtx.q * genCtx.nex + 1);
450ce27a4eeSVijay Mahadevan     genCtx.xstride = 1;
451ce27a4eeSVijay Mahadevan     break;
452ce27a4eeSVijay Mahadevan   }
453ce27a4eeSVijay Mahadevan 
454755f3dfbSVijay Mahadevan   /* Lets check for some valid input */
4555f80ce2aSJacob Faibussowitsch   PetscCheck(genCtx.dim >= 1 && genCtx.dim <= 3, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid topological dimension specified: %" PetscInt_FMT ".", genCtx.dim);
4569371c9d4SSatish 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,
4579371c9d4SSatish Balay              genCtx.N, genCtx.K, nprocs);
458755f3dfbSVijay Mahadevan   /* validate the bounds data */
4595f80ce2aSJacob 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]);
4601dca8a05SBarry 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]);
4611dca8a05SBarry 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]);
462755f3dfbSVijay Mahadevan 
4639566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Local elements:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.blockSizeElementXYZ[0], genCtx.blockSizeElementXYZ[1], genCtx.blockSizeElementXYZ[2]));
4649566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Local vertices:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.blockSizeVertexXYZ[0], genCtx.blockSizeVertexXYZ[1], genCtx.blockSizeVertexXYZ[2]));
4659566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Local blocks/processors := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.A, genCtx.B, genCtx.C));
4669566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Local processors := %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.M, genCtx.N, genCtx.K));
4679566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Local nexyz:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.nex, genCtx.ney, genCtx.nez));
4689566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Local delxyz:= %g, %g, %g\n", genCtx.dx, genCtx.dy, genCtx.dz));
4699566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "Local strides:= %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "\n", genCtx.xstride, genCtx.ystride, genCtx.zstride));
470ce27a4eeSVijay Mahadevan   PetscFunctionReturn(0);
471ce27a4eeSVijay Mahadevan }
472ce27a4eeSVijay Mahadevan 
473cab5ea25SPierre Jolivet /*@C
474ce27a4eeSVijay Mahadevan   DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with genCtx specified bounds.
475aa859218SVijay Mahadevan 
476d083f849SBarry Smith   Collective
477aa859218SVijay Mahadevan 
478aa859218SVijay Mahadevan   Input Parameters:
479aa859218SVijay Mahadevan + comm - The communicator for the DM object
480aa859218SVijay Mahadevan . dim - The spatial dimension
481b8ecf6d3SVijay 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
482b8ecf6d3SVijay Mahadevan . nele - The number of discrete elements in each direction
483a2b725a8SWilliam Gropp - user_nghost - The number of ghosted layers needed in the partitioned mesh
484aa859218SVijay Mahadevan 
485aa859218SVijay Mahadevan   Output Parameter:
486aa859218SVijay Mahadevan . dm  - The DM object
487aa859218SVijay Mahadevan 
488aa859218SVijay Mahadevan   Level: beginner
489aa859218SVijay Mahadevan 
490db781477SPatrick Sanan .seealso: `DMSetType()`, `DMCreate()`, `DMMoabLoadFromFile()`
491aa859218SVijay Mahadevan @*/
4929371c9d4SSatish Balay PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal *bounds, PetscInt nele, PetscInt nghost, DM *dm) {
49351d15aeeSVijay Mahadevan   moab::ErrorCode  merr;
4949daf19fdSVijay Mahadevan   PetscInt         a, b, c, n, global_size, global_rank;
49551d15aeeSVijay Mahadevan   DM_Moab         *dmmoab;
496ce27a4eeSVijay Mahadevan   moab::Interface *mbImpl;
4979daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
49851d15aeeSVijay Mahadevan   moab::ParallelComm *pcomm;
4999daf19fdSVijay Mahadevan #endif
500a4717116SVijay Mahadevan   moab::ReadUtilIface   *readMeshIface;
501ce27a4eeSVijay Mahadevan   moab::Range            verts, cells, edges, faces, adj, dim3, dim2;
502ce27a4eeSVijay Mahadevan   DMMoabMeshGeneratorCtx genCtx;
503a4717116SVijay Mahadevan   const PetscInt         npts = nele + 1; /* Number of points in every dimension */
504ce27a4eeSVijay Mahadevan 
5056ea892caSVijay Mahadevan   moab::Tag          global_id_tag, part_tag, geom_tag, mat_tag, dir_tag, neu_tag;
506ce27a4eeSVijay Mahadevan   moab::Range        ownedvtx, ownedelms, localvtxs, localelms;
507ce27a4eeSVijay Mahadevan   moab::EntityHandle regionset;
508755f3dfbSVijay Mahadevan   PetscInt           ml = 0, nl = 0, kl = 0;
50951d15aeeSVijay Mahadevan 
51051d15aeeSVijay Mahadevan   PetscFunctionBegin;
5111dca8a05SBarry Smith   PetscCheck(dim >= 1 && dim <= 3, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension argument for mesh: dim=[1,3].");
512e5136372SVijay Mahadevan 
5139566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("GenerateMesh", DM_CLASSID, &genCtx.generateMesh));
5149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("AddVertices", DM_CLASSID, &genCtx.generateVertices));
5159566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("AddElements", DM_CLASSID, &genCtx.generateElements));
5169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("ParResolve", DM_CLASSID, &genCtx.parResolve));
5179566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(genCtx.generateMesh, 0, 0, 0, 0));
5189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &global_size));
519e5136372SVijay Mahadevan   /* total number of vertices in all dimensions */
520e5136372SVijay Mahadevan   n = pow(npts, dim);
521e5136372SVijay Mahadevan 
522e5136372SVijay Mahadevan   /* do some error checking */
52308401ef6SPierre Jolivet   PetscCheck(n >= 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of points must be >= 2.");
52408401ef6SPierre 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.");
52508401ef6SPierre Jolivet   PetscCheck(nghost >= 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of ghost layers cannot be negative.");
526e5136372SVijay Mahadevan 
527a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
5289566063dSJacob Faibussowitsch   PetscCall(DMMoabCreateMoab(comm, NULL, NULL, NULL, dm));
52951d15aeeSVijay Mahadevan 
530a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
53151d15aeeSVijay Mahadevan   dmmoab = (DM_Moab *)(*dm)->data;
532ce27a4eeSVijay Mahadevan   mbImpl = dmmoab->mbiface;
5339daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
53451d15aeeSVijay Mahadevan   pcomm       = dmmoab->pcomm;
535ce27a4eeSVijay Mahadevan   global_rank = pcomm->rank();
5369daf19fdSVijay Mahadevan #else
5379daf19fdSVijay Mahadevan   global_rank = 0;
5389daf19fdSVijay Mahadevan   global_size = 1;
5399daf19fdSVijay Mahadevan #endif
5409daf19fdSVijay Mahadevan   global_id_tag       = dmmoab->ltog_tag;
541c8d5365dSVijay Mahadevan   dmmoab->dim         = dim;
54249d66b22SVijay Mahadevan   dmmoab->nghostrings = nghost;
543304006b3SVijay Mahadevan   dmmoab->refct       = 1;
54451d15aeeSVijay Mahadevan 
545b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
5469371c9d4SSatish Balay   merr = mbImpl->create_meshset(moab::MESHSET_SET, dmmoab->fileset);
5479371c9d4SSatish Balay   MBERR("Creating file set failed", merr);
548b5410836SVijay Mahadevan 
549a4717116SVijay Mahadevan   /* No errors yet; proceed with building the mesh */
5509371c9d4SSatish Balay   merr = mbImpl->query_interface(readMeshIface);
5519371c9d4SSatish Balay   MBERRNM(merr);
55251d15aeeSVijay Mahadevan 
553755f3dfbSVijay Mahadevan   genCtx.M = genCtx.N = genCtx.K = 1;
554755f3dfbSVijay Mahadevan   genCtx.A = genCtx.B = genCtx.C = 1;
555da8b1a3eSBarry Smith   genCtx.blockSizeElementXYZ[0]  = 0;
556da8b1a3eSBarry Smith   genCtx.blockSizeElementXYZ[1]  = 0;
557da8b1a3eSBarry Smith   genCtx.blockSizeElementXYZ[2]  = 0;
558755f3dfbSVijay Mahadevan 
559d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "DMMoab Creation Options", "DMMOAB");
560ce27a4eeSVijay Mahadevan   /* Handle DMMoab spatial resolution */
5619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dmb_grid_x", "Number of grid points in x direction", "DMMoabSetSizes", genCtx.blockSizeElementXYZ[0], &genCtx.blockSizeElementXYZ[0], &genCtx.usrxyzgrid));
5629566063dSJacob 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));
5639566063dSJacob 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));
564a4717116SVijay Mahadevan 
5656aad120cSJose E. Roman   /* Handle DMMoab parallel distribution */
5669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dmb_processors_x", "Number of processors in x direction", "DMMoabSetNumProcs", genCtx.M, &genCtx.M, &genCtx.usrprocgrid));
5679566063dSJacob Faibussowitsch   if (dim > 1) PetscCall(PetscOptionsInt("-dmb_processors_y", "Number of processors in y direction", "DMMoabSetNumProcs", genCtx.N, &genCtx.N, &genCtx.usrprocgrid));
5689566063dSJacob Faibussowitsch   if (dim > 2) PetscCall(PetscOptionsInt("-dmb_processors_z", "Number of processors in z direction", "DMMoabSetNumProcs", genCtx.K, &genCtx.K, &genCtx.usrprocgrid));
56951d15aeeSVijay Mahadevan 
570ce27a4eeSVijay Mahadevan   /* Handle DMMoab block refinement */
5719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dmb_refine_x", "Number of refinement blocks in x direction", "DMMoabSetRefinement", genCtx.A, &genCtx.A, &genCtx.usrrefgrid));
5729566063dSJacob Faibussowitsch   if (dim > 1) PetscCall(PetscOptionsInt("-dmb_refine_y", "Number of refinement blocks in y direction", "DMMoabSetRefinement", genCtx.B, &genCtx.B, &genCtx.usrrefgrid));
5739566063dSJacob Faibussowitsch   if (dim > 2) PetscCall(PetscOptionsInt("-dmb_refine_z", "Number of refinement blocks in z direction", "DMMoabSetRefinement", genCtx.C, &genCtx.C, &genCtx.usrrefgrid));
574d0609cedSBarry Smith   PetscOptionsEnd();
57551d15aeeSVijay Mahadevan 
5769566063dSJacob Faibussowitsch   PetscCall(DMMBUtil_InitializeOptions(genCtx, dim, useSimplex, global_rank, global_size, bounds, nele));
57751d15aeeSVijay Mahadevan 
57863a3b9bcSJacob 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);
57964e1c140SVijay Mahadevan 
580ce27a4eeSVijay Mahadevan   if (genCtx.adjEnts) genCtx.keep_skins = true; /* do not delete anything - consumes more memory */
58166f88a78SVijay Mahadevan 
582ce27a4eeSVijay Mahadevan   /* determine m, n, k for processor rank */
583755f3dfbSVijay Mahadevan   ml = nl = kl = 0;
584755f3dfbSVijay Mahadevan   switch (genCtx.dim) {
5859371c9d4SSatish Balay   case 1: ml = (genCtx.cumfraction); break;
5869371c9d4SSatish Balay   case 2: nl = (genCtx.cumfraction); break;
587755f3dfbSVijay Mahadevan   default:
588755f3dfbSVijay Mahadevan     kl = (genCtx.cumfraction) / genCtx.q / genCtx.blockSizeElementXYZ[2] / genCtx.C; //genCtx.K
589755f3dfbSVijay Mahadevan     break;
590755f3dfbSVijay Mahadevan   }
591e5136372SVijay Mahadevan 
592ce27a4eeSVijay Mahadevan   /*
593ce27a4eeSVijay Mahadevan    * so there are a total of M * A * blockSizeElement elements in x direction (so M * A * blockSizeElement + 1 verts in x direction)
594ce27a4eeSVijay Mahadevan    * so there are a total of N * B * blockSizeElement elements in y direction (so N * B * blockSizeElement + 1 verts in y direction)
595ce27a4eeSVijay Mahadevan    * so there are a total of K * C * blockSizeElement elements in z direction (so K * C * blockSizeElement + 1 verts in z direction)
596c8d5365dSVijay Mahadevan 
597ce27a4eeSVijay Mahadevan    * there are ( M * A blockSizeElement)      *  ( N * B * blockSizeElement)      * (K * C * blockSizeElement)    hexas
598ce27a4eeSVijay Mahadevan    * there are ( M * A * blockSizeElement + 1) *  ( N * B * blockSizeElement + 1) * (K * C * blockSizeElement + 1) vertices
599ce27a4eeSVijay Mahadevan    * x is the first dimension that varies
600ce27a4eeSVijay Mahadevan    */
601e5136372SVijay Mahadevan 
602ce27a4eeSVijay Mahadevan   /* generate the block at (a, b, c); it will represent a partition , it will get a partition tag */
60349d66b22SVijay Mahadevan   PetscInt dum_id = -1;
6049371c9d4SSatish Balay   merr            = mbImpl->tag_get_handle("GLOBAL_ID", 1, moab::MB_TYPE_INTEGER, global_id_tag);
6059371c9d4SSatish Balay   MBERR("Getting Global_ID Tag handle failed", merr);
60651d15aeeSVijay Mahadevan 
6079371c9d4SSatish Balay   merr = mbImpl->tag_get_handle(MATERIAL_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, mat_tag);
6089371c9d4SSatish Balay   MBERR("Getting Material set Tag handle failed", merr);
6099371c9d4SSatish Balay   merr = mbImpl->tag_get_handle(DIRICHLET_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, dir_tag);
6109371c9d4SSatish Balay   MBERR("Getting Dirichlet set Tag handle failed", merr);
6119371c9d4SSatish Balay   merr = mbImpl->tag_get_handle(NEUMANN_SET_TAG_NAME, 1, moab::MB_TYPE_INTEGER, neu_tag);
6129371c9d4SSatish Balay   MBERR("Getting Neumann set Tag handle failed", merr);
6136ea892caSVijay Mahadevan 
6149371c9d4SSatish 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);
6159371c9d4SSatish Balay   MBERR("Getting Partition Tag handle failed", merr);
61651d15aeeSVijay Mahadevan 
6173a4aeca1SVijay Mahadevan   /* lets create some sets */
6189371c9d4SSatish 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);
6199371c9d4SSatish Balay   MBERRNM(merr);
6209371c9d4SSatish Balay   merr = mbImpl->create_meshset(moab::MESHSET_SET, regionset);
6219371c9d4SSatish Balay   MBERRNM(merr);
6229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(genCtx.generateMesh, 0, 0, 0, 0));
6233a4aeca1SVijay Mahadevan 
624ce27a4eeSVijay Mahadevan   for (a = 0; a < (genCtx.dim > 0 ? genCtx.A : genCtx.A); a++) {
625ce27a4eeSVijay Mahadevan     for (b = 0; b < (genCtx.dim > 1 ? genCtx.B : 1); b++) {
626ce27a4eeSVijay Mahadevan       for (c = 0; c < (genCtx.dim > 2 ? genCtx.C : 1); c++) {
62749d66b22SVijay Mahadevan         moab::EntityHandle startv;
62849d66b22SVijay Mahadevan 
6299566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(genCtx.generateVertices, 0, 0, 0, 0));
6309566063dSJacob Faibussowitsch         PetscCall(DMMoab_GenerateVertices_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, verts));
6319566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(genCtx.generateVertices, 0, 0, 0, 0));
6323a4aeca1SVijay Mahadevan 
6339566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(genCtx.generateElements, 0, 0, 0, 0));
6349566063dSJacob Faibussowitsch         PetscCall(DMMoab_GenerateElements_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, cells));
6359566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(genCtx.generateElements, 0, 0, 0, 0));
6363a4aeca1SVijay Mahadevan 
63749d66b22SVijay Mahadevan         PetscInt part_num = 0;
63849d66b22SVijay Mahadevan         switch (genCtx.dim) {
6399371c9d4SSatish Balay         case 3: part_num += (c + kl * genCtx.C) * (genCtx.M * genCtx.A * genCtx.N * genCtx.B);
6409371c9d4SSatish Balay         case 2: part_num += (b + nl * genCtx.B) * (genCtx.M * genCtx.A);
6419371c9d4SSatish Balay         case 1: part_num += (a + ml * genCtx.A); break;
64249d66b22SVijay Mahadevan         }
64349d66b22SVijay Mahadevan 
64449d66b22SVijay Mahadevan         moab::EntityHandle part_set;
6459371c9d4SSatish Balay         merr = mbImpl->create_meshset(moab::MESHSET_SET, part_set);
6469371c9d4SSatish Balay         MBERR("Can't create mesh set.", merr);
6473a4aeca1SVijay Mahadevan 
6489371c9d4SSatish Balay         merr = mbImpl->add_entities(part_set, verts);
6499371c9d4SSatish Balay         MBERR("Can't add vertices to set.", merr);
6509371c9d4SSatish Balay         merr = mbImpl->add_entities(part_set, cells);
6519371c9d4SSatish Balay         MBERR("Can't add entities to set.", merr);
6529371c9d4SSatish Balay         merr = mbImpl->add_entities(regionset, cells);
6539371c9d4SSatish Balay         MBERR("Can't add entities to set.", merr);
654ce27a4eeSVijay Mahadevan 
655ce27a4eeSVijay Mahadevan         /* if needed, add all edges and faces */
6569371c9d4SSatish Balay         if (genCtx.adjEnts) {
657ce27a4eeSVijay Mahadevan           if (genCtx.dim > 1) {
6589371c9d4SSatish Balay             merr = mbImpl->get_adjacencies(cells, 1, true, edges, moab::Interface::UNION);
6599371c9d4SSatish Balay             MBERR("Can't get edges", merr);
6609371c9d4SSatish Balay             merr = mbImpl->add_entities(part_set, edges);
6619371c9d4SSatish Balay             MBERR("Can't add edges to partition set.", merr);
662ce27a4eeSVijay Mahadevan           }
663ce27a4eeSVijay Mahadevan           if (genCtx.dim > 2) {
6649371c9d4SSatish Balay             merr = mbImpl->get_adjacencies(cells, 2, true, faces, moab::Interface::UNION);
6659371c9d4SSatish Balay             MBERR("Can't get faces", merr);
6669371c9d4SSatish Balay             merr = mbImpl->add_entities(part_set, faces);
6679371c9d4SSatish Balay             MBERR("Can't add faces to partition set.", merr);
668ce27a4eeSVijay Mahadevan           }
669ce27a4eeSVijay Mahadevan           edges.clear();
670ce27a4eeSVijay Mahadevan           faces.clear();
671ce27a4eeSVijay Mahadevan         }
6729371c9d4SSatish Balay         verts.clear();
6739371c9d4SSatish Balay         cells.clear();
674ce27a4eeSVijay Mahadevan 
6759371c9d4SSatish Balay         merr = mbImpl->tag_set_data(part_tag, &part_set, 1, &part_num);
6769371c9d4SSatish Balay         MBERR("Can't set part tag on set", merr);
67749d66b22SVijay Mahadevan         if (dmmoab->fileset) {
6789371c9d4SSatish Balay           merr = mbImpl->add_parent_child(dmmoab->fileset, part_set);
6799371c9d4SSatish Balay           MBERR("Can't add part set to file set.", merr);
6809371c9d4SSatish Balay           merr = mbImpl->unite_meshset(dmmoab->fileset, part_set);
6819371c9d4SSatish Balay           MBERRNM(merr);
68249d66b22SVijay Mahadevan         }
6839371c9d4SSatish Balay         merr = mbImpl->add_entities(dmmoab->fileset, &part_set, 1);
6849371c9d4SSatish Balay         MBERRNM(merr);
685ce27a4eeSVijay Mahadevan       }
686ce27a4eeSVijay Mahadevan     }
687ce27a4eeSVijay Mahadevan   }
688ce27a4eeSVijay Mahadevan 
6899371c9d4SSatish Balay   merr = mbImpl->add_parent_child(dmmoab->fileset, regionset);
6909371c9d4SSatish Balay   MBERRNM(merr);
691ce27a4eeSVijay Mahadevan 
69249d66b22SVijay Mahadevan   /* Only in parallel: resolve shared entities between processors and exchange ghost layers */
693ce27a4eeSVijay Mahadevan   if (global_size > 1) {
6949566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(genCtx.parResolve, 0, 0, 0, 0));
69577a8c067SVijay Mahadevan 
6969371c9d4SSatish Balay     merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, genCtx.dim, cells);
6979371c9d4SSatish Balay     MBERR("Can't get all d-dimensional elements.", merr);
6989371c9d4SSatish Balay     merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 0, verts);
6999371c9d4SSatish Balay     MBERR("Can't get all vertices.", merr);
700ce27a4eeSVijay Mahadevan 
701ce27a4eeSVijay Mahadevan     if (genCtx.A * genCtx.B * genCtx.C != 1) { //  merge needed
70249d66b22SVijay Mahadevan       moab::MergeMesh mm(mbImpl);
703ce27a4eeSVijay Mahadevan       if (genCtx.newMergeMethod) {
7049371c9d4SSatish Balay         merr = mm.merge_using_integer_tag(verts, global_id_tag);
7059371c9d4SSatish Balay         MBERR("Can't merge with GLOBAL_ID tag", merr);
7069371c9d4SSatish Balay       } else {
7079371c9d4SSatish Balay         merr = mm.merge_entities(cells, 0.0001);
7089371c9d4SSatish Balay         MBERR("Can't merge with coordinates", merr);
7093a4aeca1SVijay Mahadevan       }
7103a4aeca1SVijay Mahadevan     }
7113a4aeca1SVijay Mahadevan 
7129daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
713a4717116SVijay Mahadevan     /* check the handles */
7149371c9d4SSatish Balay     merr = pcomm->check_all_shared_handles();
7159371c9d4SSatish Balay     MBERRV(mbImpl, merr);
71651d15aeeSVijay Mahadevan 
717a4717116SVijay Mahadevan     /* resolve the shared entities by exchanging information to adjacent processors */
7189371c9d4SSatish Balay     merr = pcomm->resolve_shared_ents(dmmoab->fileset, cells, dim, dim - 1, NULL, &global_id_tag);
7199371c9d4SSatish Balay     MBERRV(mbImpl, merr);
72049d66b22SVijay Mahadevan     if (dmmoab->fileset) {
7219371c9d4SSatish Balay       merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false, &dmmoab->fileset);
7229371c9d4SSatish Balay       MBERRV(mbImpl, merr);
7239371c9d4SSatish Balay     } else {
7249371c9d4SSatish Balay       merr = pcomm->exchange_ghost_cells(dim, 0, nghost, dim, true, false);
7259371c9d4SSatish Balay       MBERRV(mbImpl, merr);
72649d66b22SVijay Mahadevan     }
727c8d5365dSVijay Mahadevan 
728e427d9c9SVijay Mahadevan     /* Reassign global IDs on all entities. */
7299371c9d4SSatish Balay     merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, false, true, false);
7309371c9d4SSatish Balay     MBERRNM(merr);
7319daf19fdSVijay Mahadevan #endif
73277a8c067SVijay Mahadevan 
7339566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(genCtx.parResolve, 0, 0, 0, 0));
73449d66b22SVijay Mahadevan   }
7353f1c6e43SVijay Mahadevan 
736ce27a4eeSVijay Mahadevan   if (!genCtx.keep_skins) { // default is to delete the 1- and 2-dimensional entities
737ce27a4eeSVijay Mahadevan     // delete all quads and edges
738ce27a4eeSVijay Mahadevan     moab::Range toDelete;
739ce27a4eeSVijay Mahadevan     if (genCtx.dim > 1) {
7409371c9d4SSatish Balay       merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 1, toDelete);
7419371c9d4SSatish Balay       MBERR("Can't get edges", merr);
742ce27a4eeSVijay Mahadevan     }
7433f1c6e43SVijay Mahadevan 
744ce27a4eeSVijay Mahadevan     if (genCtx.dim > 2) {
7459371c9d4SSatish Balay       merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 2, toDelete);
7469371c9d4SSatish Balay       MBERR("Can't get faces", merr);
747ce27a4eeSVijay Mahadevan     }
748c8d5365dSVijay Mahadevan 
7499daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
7509371c9d4SSatish Balay     merr = dmmoab->pcomm->delete_entities(toDelete);
7519371c9d4SSatish Balay     MBERR("Can't delete entities", merr);
7529daf19fdSVijay Mahadevan #endif
753ce27a4eeSVijay Mahadevan   }
7546ea892caSVijay Mahadevan 
7556ea892caSVijay Mahadevan   /* set geometric dimension tag for regions */
7569371c9d4SSatish Balay   merr = mbImpl->tag_set_data(geom_tag, &regionset, 1, &dmmoab->dim);
7579371c9d4SSatish Balay   MBERRNM(merr);
7586ea892caSVijay Mahadevan   /* set default material ID for regions */
7596ea892caSVijay Mahadevan   int default_material = 1;
7609371c9d4SSatish Balay   merr                 = mbImpl->tag_set_data(mat_tag, &regionset, 1, &default_material);
7619371c9d4SSatish Balay   MBERRNM(merr);
7626ea892caSVijay Mahadevan   /*
7636ea892caSVijay Mahadevan     int default_dbc = 0;
7646ea892caSVijay Mahadevan     merr = mbImpl->tag_set_data(dir_tag, &vertexset, 1, &default_dbc);MBERRNM(merr);
7656ea892caSVijay Mahadevan   */
76651d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
76751d15aeeSVijay Mahadevan }
76851d15aeeSVijay Mahadevan 
7699371c9d4SSatish Balay 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) {
7706acfe860SVijay Mahadevan   char *ropts;
77149d66b22SVijay Mahadevan   char  ropts_par[PETSC_MAX_PATH_LEN], ropts_pargh[PETSC_MAX_PATH_LEN];
77261eb6e27SVijay Mahadevan   char  ropts_dbg[PETSC_MAX_PATH_LEN];
77351d15aeeSVijay Mahadevan 
774a4717116SVijay Mahadevan   PetscFunctionBegin;
7759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(PETSC_MAX_PATH_LEN, &ropts));
7769566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&ropts_par, PETSC_MAX_PATH_LEN));
7779566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&ropts_pargh, PETSC_MAX_PATH_LEN));
7789566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&ropts_dbg, PETSC_MAX_PATH_LEN));
77961eb6e27SVijay Mahadevan 
780e23c60ebSVijay Mahadevan   /* do parallel read unless using only one processor */
781a4717116SVijay Mahadevan   if (numproc > 1) {
7829566063dSJacob 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":"")));
7839566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(ropts_par, PETSC_MAX_PATH_LEN, "PARALLEL=%s;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS;%s", MoabReadModes[mode], (by_rank ? "PARTITION_BY_RANK;" : "")));
784*48a46eb9SPierre Jolivet     if (nghost) PetscCall(PetscSNPrintf(ropts_pargh, PETSC_MAX_PATH_LEN, "PARALLEL_GHOSTS=%" PetscInt_FMT ".0.%" PetscInt_FMT ";", dim, nghost));
7852e4e7c01SVijay Mahadevan   }
7862e4e7c01SVijay Mahadevan 
787c8d5365dSVijay Mahadevan   if (dbglevel) {
788f90c3b0eSVijay Mahadevan     if (numproc > 1) {
7899566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%" PetscInt_FMT ";DEBUG_PIO=%" PetscInt_FMT ";", dbglevel, dbglevel));
7909371c9d4SSatish Balay     } else PetscCall(PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%" PetscInt_FMT ";", dbglevel));
791c8d5365dSVijay Mahadevan   }
79251d15aeeSVijay Mahadevan 
7939566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(ropts, PETSC_MAX_PATH_LEN, "%s%s%s%s%s", ropts_par, (nghost ? ropts_pargh : ""), ropts_dbg, (extra_opts ? extra_opts : ""), (dm_opts ? dm_opts : "")));
794f90c3b0eSVijay Mahadevan   *read_opts = ropts;
79551d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
79651d15aeeSVijay Mahadevan }
79751d15aeeSVijay Mahadevan 
798cab5ea25SPierre Jolivet /*@C
799aa859218SVijay Mahadevan   DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file.
800aa859218SVijay Mahadevan 
801d083f849SBarry Smith   Collective
802aa859218SVijay Mahadevan 
803aa859218SVijay Mahadevan   Input Parameters:
804aa859218SVijay Mahadevan + comm - The communicator for the DM object
805aa859218SVijay Mahadevan . dim - The spatial dimension
806b8ecf6d3SVijay Mahadevan . filename - The name of the mesh file to be loaded
807a2b725a8SWilliam Gropp - usrreadopts - The options string to read a MOAB mesh.
808a2b725a8SWilliam Gropp 
809a8d69d7bSBarry Smith   Reference (Parallel Mesh Initialization: https://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo)
810aa859218SVijay Mahadevan 
811aa859218SVijay Mahadevan   Output Parameter:
812aa859218SVijay Mahadevan . dm  - The DM object
813aa859218SVijay Mahadevan 
814aa859218SVijay Mahadevan   Level: beginner
815aa859218SVijay Mahadevan 
816db781477SPatrick Sanan .seealso: `DMSetType()`, `DMCreate()`, `DMMoabCreateBoxMesh()`
817aa859218SVijay Mahadevan @*/
8189371c9d4SSatish Balay PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm, PetscInt dim, PetscInt nghost, const char *filename, const char *usrreadopts, DM *dm) {
819a4717116SVijay Mahadevan   moab::ErrorCode  merr;
8202e4e7c01SVijay Mahadevan   PetscInt         nprocs;
821a4717116SVijay Mahadevan   DM_Moab         *dmmoab;
822a4717116SVijay Mahadevan   moab::Interface *mbiface;
8239daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
824a4717116SVijay Mahadevan   moab::ParallelComm *pcomm;
8259daf19fdSVijay Mahadevan #endif
826a4717116SVijay Mahadevan   moab::Range verts, elems;
8277023aa44SVijay Mahadevan   const char *readopts;
82851d15aeeSVijay Mahadevan 
82951d15aeeSVijay Mahadevan   PetscFunctionBegin;
83049d66b22SVijay Mahadevan   PetscValidPointer(dm, 6);
83151d15aeeSVijay Mahadevan 
832a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
8339566063dSJacob Faibussowitsch   PetscCall(DMMoabCreateMoab(comm, NULL, NULL, NULL, dm));
83451d15aeeSVijay Mahadevan 
835a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
836a4717116SVijay Mahadevan   dmmoab  = (DM_Moab *)(*dm)->data;
837a4717116SVijay Mahadevan   mbiface = dmmoab->mbiface;
8389daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
839a4717116SVijay Mahadevan   pcomm  = dmmoab->pcomm;
840a4717116SVijay Mahadevan   nprocs = pcomm->size();
8419daf19fdSVijay Mahadevan #else
8429daf19fdSVijay Mahadevan   nprocs      = 1;
8439daf19fdSVijay Mahadevan #endif
844aa859218SVijay Mahadevan   /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */
845c8d5365dSVijay Mahadevan   dmmoab->dim         = dim;
84649d66b22SVijay Mahadevan   dmmoab->nghostrings = nghost;
847304006b3SVijay Mahadevan   dmmoab->refct       = 1;
848a4717116SVijay Mahadevan 
849b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
8509371c9d4SSatish Balay   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);
8519371c9d4SSatish Balay   MBERR("Creating file set failed", merr);
852b5410836SVijay Mahadevan 
853a4717116SVijay Mahadevan   /* add mesh loading options specific to the DM */
8549371c9d4SSatish Balay   PetscCall(DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, nghost, dmmoab->read_mode, dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts));
8557023aa44SVijay Mahadevan 
8569566063dSJacob Faibussowitsch   PetscCall(PetscInfo(*dm, "Reading file %s with options: %s\n", filename, readopts));
857a4717116SVijay Mahadevan 
858a4717116SVijay Mahadevan   /* Load the mesh from a file. */
85949d66b22SVijay Mahadevan   if (dmmoab->fileset) {
8609371c9d4SSatish Balay     merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);
8619371c9d4SSatish Balay     MBERRVM(mbiface, "Reading MOAB file failed.", merr);
8629371c9d4SSatish Balay   } else {
8639371c9d4SSatish Balay     merr = mbiface->load_file(filename, 0, readopts);
8649371c9d4SSatish Balay     MBERRVM(mbiface, "Reading MOAB file failed.", merr);
86549d66b22SVijay Mahadevan   }
866a4717116SVijay Mahadevan 
8679daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
8686e40195eSVijay Mahadevan   /* Reassign global IDs on all entities. */
8696ea892caSVijay Mahadevan   /* merr = pcomm->assign_global_ids(dmmoab->fileset, dim, 1, true, true, true);MBERRNM(merr); */
8709daf19fdSVijay Mahadevan #endif
871e5136372SVijay Mahadevan 
872e5136372SVijay Mahadevan   /* load the local vertices */
8739371c9d4SSatish Balay   merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);
8749371c9d4SSatish Balay   MBERRNM(merr);
875e5136372SVijay Mahadevan   /* load the local elements */
8769371c9d4SSatish Balay   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);
8779371c9d4SSatish Balay   MBERRNM(merr);
878e5136372SVijay Mahadevan 
8799daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
880e5136372SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
881e5136372SVijay Mahadevan      on all of the ghost vertexes */
8829371c9d4SSatish Balay   merr = pcomm->exchange_tags(dmmoab->ltog_tag, verts);
8839371c9d4SSatish Balay   MBERRV(mbiface, merr);
8849371c9d4SSatish Balay   merr = pcomm->exchange_tags(dmmoab->ltog_tag, elems);
8859371c9d4SSatish Balay   MBERRV(mbiface, merr);
8869371c9d4SSatish Balay   merr = pcomm->collective_sync_partition();
8879371c9d4SSatish Balay   MBERR("Collective sync failed", merr);
8889daf19fdSVijay Mahadevan #endif
889a4717116SVijay Mahadevan 
89063a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(*dm, "MOAB file '%s' was successfully loaded. Found %zu vertices and %zu elements.\n", filename, verts.size(), elems.size()));
8919566063dSJacob Faibussowitsch   PetscCall(PetscFree(readopts));
89251d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
89351d15aeeSVijay Mahadevan }
89451d15aeeSVijay Mahadevan 
895cab5ea25SPierre Jolivet /*@C
896304006b3SVijay Mahadevan   DMMoabRenumberMeshEntities - Order and number all entities (vertices->elements) to be contiguously ordered
897304006b3SVijay Mahadevan   in parallel
898304006b3SVijay Mahadevan 
899d083f849SBarry Smith   Collective
900304006b3SVijay Mahadevan 
901304006b3SVijay Mahadevan   Input Parameters:
902a2b725a8SWilliam Gropp . dm  - The DM object
903304006b3SVijay Mahadevan 
904304006b3SVijay Mahadevan   Level: advanced
905304006b3SVijay Mahadevan 
906db781477SPatrick Sanan .seealso: `DMSetUp()`, `DMCreate()`
907304006b3SVijay Mahadevan @*/
9089371c9d4SSatish Balay PetscErrorCode DMMoabRenumberMeshEntities(DM dm) {
909304006b3SVijay Mahadevan   moab::Range verts;
910304006b3SVijay Mahadevan 
911304006b3SVijay Mahadevan   PetscFunctionBegin;
912304006b3SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
913304006b3SVijay Mahadevan 
9149daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
915304006b3SVijay Mahadevan   /* Insert new points */
9169daf19fdSVijay Mahadevan   moab::ErrorCode merr;
9179371c9d4SSatish Balay   merr = ((DM_Moab *)dm->data)->pcomm->assign_global_ids(((DM_Moab *)dm->data)->fileset, 3, 0, false, true, false);
9189371c9d4SSatish Balay   MBERRNM(merr);
9199daf19fdSVijay Mahadevan #endif
920304006b3SVijay Mahadevan   PetscFunctionReturn(0);
921304006b3SVijay Mahadevan }
922