xref: /petsc/src/dm/impls/moab/dmmbutil.cxx (revision 49d66b22a367b93e538da1d54516fdaf6756fee1)
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
12*49d66b22SVijay 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;
23ce27a4eeSVijay Mahadevan } DMMoabMeshGeneratorCtx;
24ce27a4eeSVijay Mahadevan 
25ce27a4eeSVijay Mahadevan 
26ce27a4eeSVijay Mahadevan #undef __FUNCT__
27ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_SetTensorElementConnectivity_Private"
28*49d66b22SVijay Mahadevan PetscInt DMMoab_SetTensorElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt offset, PetscInt corner, std::vector<PetscInt>& subent_conn, moab::EntityHandle *connectivity)
2951d15aeeSVijay Mahadevan {
30ce27a4eeSVijay Mahadevan   switch(genCtx.dim) {
31e5136372SVijay Mahadevan     case 1:
32ce27a4eeSVijay Mahadevan       subent_conn.resize(2);
33ce27a4eeSVijay Mahadevan       moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data());
34ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[0]] = corner;
35ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[1]] = corner + 1;
36e5136372SVijay Mahadevan       break;
37e5136372SVijay Mahadevan     case 2:
38ce27a4eeSVijay Mahadevan       subent_conn.resize(4);
39ce27a4eeSVijay Mahadevan       moab::CN::SubEntityVertexIndices(moab::MBQUAD, 2, 0, subent_conn.data());
40ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[0]] = corner;
41ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[1]] = corner + 1;
42ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[2]] = corner + 1 + genCtx.ystride;
43ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[3]] = corner + genCtx.ystride;
44e5136372SVijay Mahadevan       break;
45e5136372SVijay Mahadevan     case 3:
46e5136372SVijay Mahadevan     default:
47ce27a4eeSVijay Mahadevan       subent_conn.resize(8);
48ce27a4eeSVijay Mahadevan       moab::CN::SubEntityVertexIndices(moab::MBHEX, 3, 0, subent_conn.data());
49ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[0]] = corner;
50ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[1]] = corner + 1;
51ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[2]] = corner + 1 + genCtx.ystride;
52ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[3]] = corner + genCtx.ystride;
53ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[4]] = corner + genCtx.zstride;
54ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[5]] = corner + 1 + genCtx.zstride;
55ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[6]] = corner + 1 + genCtx.ystride + genCtx.zstride;
56ce27a4eeSVijay Mahadevan       connectivity[offset+subent_conn[7]] = corner + genCtx.ystride + genCtx.zstride;
57e5136372SVijay Mahadevan       break;
58e5136372SVijay Mahadevan   }
59ce27a4eeSVijay Mahadevan   return subent_conn.size();
60e5136372SVijay Mahadevan }
6151d15aeeSVijay Mahadevan 
62ce27a4eeSVijay Mahadevan 
63ce27a4eeSVijay Mahadevan #undef __FUNCT__
64ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_SetSimplexElementConnectivity_Private"
65*49d66b22SVijay Mahadevan PetscInt DMMoab_SetSimplexElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt subelem, PetscInt offset, PetscInt corner, std::vector<PetscInt>& subent_conn, moab::EntityHandle *connectivity)
66c3dd00c7SVijay Mahadevan {
67*49d66b22SVijay Mahadevan   PetscInt A, B, C, D, E, F, G, H;
68*49d66b22SVijay Mahadevan   A = corner;
69*49d66b22SVijay Mahadevan   B = corner + 1;
70ce27a4eeSVijay Mahadevan   switch(genCtx.dim) {
71c3dd00c7SVijay Mahadevan     case 1:
72ce27a4eeSVijay Mahadevan       subent_conn.resize(2);  /* only linear EDGE supported now */
73ce27a4eeSVijay Mahadevan       moab::CN::SubEntityVertexIndices(moab::MBEDGE, 1, 0, subent_conn.data());
74*49d66b22SVijay Mahadevan       connectivity[offset+subent_conn[0]] = A;
75*49d66b22SVijay Mahadevan       connectivity[offset+subent_conn[1]] = B;
76c3dd00c7SVijay Mahadevan       break;
77c3dd00c7SVijay Mahadevan     case 2:
78*49d66b22SVijay Mahadevan       C = corner + 1 + genCtx.ystride;
79*49d66b22SVijay Mahadevan       D = corner +     genCtx.ystride;
80ce27a4eeSVijay Mahadevan       subent_conn.resize(3);  /* only linear TRI supported */
81ce27a4eeSVijay Mahadevan       moab::CN::SubEntityVertexIndices(moab::MBTRI, 2, 0, subent_conn.data());
82ce27a4eeSVijay Mahadevan       if (subelem) { /* 0 1 2 of a QUAD */
83*49d66b22SVijay Mahadevan         connectivity[offset+subent_conn[0]] = A;
84*49d66b22SVijay Mahadevan         connectivity[offset+subent_conn[1]] = B;
85*49d66b22SVijay Mahadevan         connectivity[offset+subent_conn[2]] = C;
86c3dd00c7SVijay Mahadevan       }
87ce27a4eeSVijay Mahadevan       else {        /* 2 3 0 of a QUAD */
88*49d66b22SVijay Mahadevan         connectivity[offset+subent_conn[0]] = C;
89*49d66b22SVijay Mahadevan         connectivity[offset+subent_conn[1]] = D;
90*49d66b22SVijay Mahadevan         connectivity[offset+subent_conn[2]] = A;
91c3dd00c7SVijay Mahadevan       }
92c3dd00c7SVijay Mahadevan       break;
93c3dd00c7SVijay Mahadevan     case 3:
94c3dd00c7SVijay Mahadevan     default:
95*49d66b22SVijay Mahadevan       C = corner + 1 + genCtx.ystride;
96*49d66b22SVijay Mahadevan       D = corner +     genCtx.ystride;
97*49d66b22SVijay Mahadevan       E = corner +                      genCtx.zstride;
98*49d66b22SVijay Mahadevan       F = corner + 1 +                  genCtx.zstride;
99*49d66b22SVijay Mahadevan       G = corner + 1 + genCtx.ystride + genCtx.zstride;
100*49d66b22SVijay Mahadevan       H = corner +     genCtx.ystride + genCtx.zstride;
101ce27a4eeSVijay Mahadevan       subent_conn.resize(4);  /* only linear TET supported */
102ce27a4eeSVijay Mahadevan       moab::CN::SubEntityVertexIndices(moab::MBTET, 3, 0, subent_conn.data());
103c3dd00c7SVijay Mahadevan       switch(subelem) {
104*49d66b22SVijay Mahadevan         case 0: /* 4 3 7 6 of a HEX */
105*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[0]] = E;
106*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[1]] = D;
107*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[2]] = H;
108*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[3]] = G;
109c3dd00c7SVijay Mahadevan           break;
110*49d66b22SVijay Mahadevan         case 1: /* 0 1 2 5 of a HEX */
111*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[0]] = A;
112*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[1]] = B;
113*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[2]] = C;
114*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[3]] = F;
115c3dd00c7SVijay Mahadevan           break;
116*49d66b22SVijay Mahadevan         case 2: /* 0 3 4 5 of a HEX */
117*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[0]] = A;
118*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[1]] = D;
119*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[2]] = E;
120*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[3]] = F;
121c3dd00c7SVijay Mahadevan           break;
122*49d66b22SVijay Mahadevan         case 3: /* 2 6 3 5 of a HEX */
123*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[0]] = C;
124*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[1]] = G;
125*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[2]] = D;
126*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[3]] = F;
127c3dd00c7SVijay Mahadevan           break;
128*49d66b22SVijay Mahadevan         case 4: /* 0 2 3 5 of a HEX */
129*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[0]] = A;
130*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[1]] = C;
131*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[2]] = D;
132*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[3]] = F;
133*49d66b22SVijay Mahadevan           break;
134*49d66b22SVijay Mahadevan         case 5: /* 3 6 4 5 of a HEX */
135*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[0]] = D;
136*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[1]] = G;
137*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[2]] = E;
138*49d66b22SVijay Mahadevan           connectivity[offset+subent_conn[3]] = F;
139c3dd00c7SVijay Mahadevan           break;
140c3dd00c7SVijay Mahadevan       }
141c3dd00c7SVijay Mahadevan       break;
142c3dd00c7SVijay Mahadevan   }
143ce27a4eeSVijay Mahadevan   return subent_conn.size();
144c3dd00c7SVijay Mahadevan }
145c3dd00c7SVijay Mahadevan 
146ce27a4eeSVijay Mahadevan 
147ce27a4eeSVijay Mahadevan #undef __FUNCT__
148ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_SetElementConnectivity_Private"
149*49d66b22SVijay Mahadevan std::pair<PetscInt,PetscInt> DMMoab_SetElementConnectivity_Private(DMMoabMeshGeneratorCtx& genCtx, PetscInt offset, PetscInt corner, moab::EntityHandle *connectivity)
150c3dd00c7SVijay Mahadevan {
151*49d66b22SVijay Mahadevan   PetscInt vcount=0;
152*49d66b22SVijay Mahadevan   PetscInt simplices_per_tensor[4] = {0,1,2,6};
153*49d66b22SVijay Mahadevan   std::vector<PetscInt> subent_conn;  /* only linear edge, tri, tet supported now */
154ce27a4eeSVijay Mahadevan   subent_conn.reserve(27);
155*49d66b22SVijay Mahadevan   PetscInt m,subelem;
156ce27a4eeSVijay Mahadevan   if (genCtx.simplex) {
157*49d66b22SVijay Mahadevan     subelem=simplices_per_tensor[genCtx.dim];
158ce27a4eeSVijay Mahadevan     for (m=0; m<subelem; m++) {
159ce27a4eeSVijay Mahadevan       vcount=DMMoab_SetSimplexElementConnectivity_Private(genCtx, m, offset, corner, subent_conn, connectivity);
160ce27a4eeSVijay Mahadevan       offset+=vcount;
161ce27a4eeSVijay Mahadevan     }
162c3dd00c7SVijay Mahadevan   }
163c3dd00c7SVijay Mahadevan   else {
164c3dd00c7SVijay Mahadevan     subelem=1;
165ce27a4eeSVijay Mahadevan     vcount=DMMoab_SetTensorElementConnectivity_Private(genCtx, offset, corner, subent_conn, connectivity);
166c3dd00c7SVijay Mahadevan   }
167*49d66b22SVijay Mahadevan   return std::pair<PetscInt,PetscInt>(vcount*subelem,subelem);
168c3dd00c7SVijay Mahadevan }
169c3dd00c7SVijay Mahadevan 
170c3dd00c7SVijay Mahadevan 
171ce27a4eeSVijay Mahadevan #undef __FUNCT__
172ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_GenerateVertices_Private"
173*49d66b22SVijay Mahadevan PetscErrorCode DMMoab_GenerateVertices_Private(moab::Interface *mbImpl, moab::ReadUtilIface *iface, DMMoabMeshGeneratorCtx& genCtx, PetscInt m, PetscInt n, PetscInt k,
174*49d66b22SVijay Mahadevan     PetscInt a, PetscInt b, PetscInt c, moab::Tag& global_id_tag, moab::EntityHandle& startv, moab::Range& uverts)
175ce27a4eeSVijay Mahadevan {
176*49d66b22SVijay Mahadevan   PetscInt x,y,z,ix,nnodes;
177*49d66b22SVijay Mahadevan   PetscInt ii,jj,kk;
178*49d66b22SVijay Mahadevan   std::vector<PetscReal*> arrays;
179*49d66b22SVijay Mahadevan   std::vector<PetscInt> gids;
180ce27a4eeSVijay Mahadevan   moab::ErrorCode merr;
181ce27a4eeSVijay Mahadevan 
182ce27a4eeSVijay Mahadevan   PetscFunctionBegin;
183ce27a4eeSVijay Mahadevan   /* we will generate (q*block+1)^3 vertices, and block^3 hexas; q is 1 for linear, 2 for quadratic
184ce27a4eeSVijay Mahadevan    * the global id of the vertices will come from m, n, k, a, b, c
185ce27a4eeSVijay Mahadevan    * x will vary from  m*A*q*block + a*q*block to m*A*q*block+(a+1)*q*block etc.
186ce27a4eeSVijay Mahadevan    */
187ce27a4eeSVijay Mahadevan   nnodes = genCtx.blockSizeVertexXYZ[0]*(genCtx.dim>1? genCtx.blockSizeVertexXYZ[1]*(genCtx.dim>2 ? genCtx.blockSizeVertexXYZ[2]:1) :1);
188ce27a4eeSVijay Mahadevan 
189ce27a4eeSVijay Mahadevan   merr = iface->get_node_coords(3, nnodes, 0, startv, arrays);MBERR("Can't get node coords.", merr);
190ce27a4eeSVijay Mahadevan 
191ce27a4eeSVijay Mahadevan   /* will start with the lower corner: */
192*49d66b22SVijay Mahadevan   x = ( m * genCtx.A + a ) * genCtx.q * genCtx.blockSizeElementXYZ[0];
193*49d66b22SVijay Mahadevan   y = ( n * genCtx.B + b ) * genCtx.q * genCtx.blockSizeElementXYZ[1];
194*49d66b22SVijay Mahadevan   z = ( k * genCtx.C + c ) * genCtx.q * genCtx.blockSizeElementXYZ[2];
195ce27a4eeSVijay Mahadevan   ix = 0;
196ce27a4eeSVijay Mahadevan   gids.resize(nnodes);
197ce27a4eeSVijay Mahadevan   moab::Range verts(startv, startv + nnodes - 1);
198ce27a4eeSVijay Mahadevan   for (kk = 0; kk < (genCtx.dim>2?genCtx.blockSizeVertexXYZ[2]:1); kk++) {
199ce27a4eeSVijay Mahadevan     for (jj = 0; jj < (genCtx.dim>1?genCtx.blockSizeVertexXYZ[1]:1); jj++) {
200ce27a4eeSVijay Mahadevan       for (ii = 0; ii < genCtx.blockSizeVertexXYZ[0]; ii++,ix++) {
201ce27a4eeSVijay Mahadevan         /* set coordinates for the vertices */
202ce27a4eeSVijay Mahadevan         arrays[0][ix] = (x + ii)*genCtx.dx + genCtx.xyzbounds[0];
203ce27a4eeSVijay Mahadevan         arrays[1][ix] = (y + jj)*genCtx.dy + genCtx.xyzbounds[2];
204ce27a4eeSVijay Mahadevan         arrays[2][ix] = (z + kk)*genCtx.dz + genCtx.xyzbounds[4];
205ce27a4eeSVijay Mahadevan 
206ce27a4eeSVijay Mahadevan         /* If we want to set some tags on the vertices -> use the following entity handle definition:
207ce27a4eeSVijay Mahadevan            moab::EntityHandle v = startv + ix;
208ce27a4eeSVijay Mahadevan         */
209ce27a4eeSVijay Mahadevan         /* compute the global ID for vertex */
210ce27a4eeSVijay Mahadevan         gids[ix] = 1 + (x + ii) + (y + jj) * genCtx.NX + (z + kk) * (genCtx.NX * genCtx.NY);
211ce27a4eeSVijay Mahadevan       }
212ce27a4eeSVijay Mahadevan     }
213ce27a4eeSVijay Mahadevan   }
214ce27a4eeSVijay Mahadevan   /* set global ID data on vertices */
215ce27a4eeSVijay Mahadevan   mbImpl->tag_set_data(global_id_tag, verts, &gids[0]);
216ce27a4eeSVijay Mahadevan   verts.swap(uverts);
217ce27a4eeSVijay Mahadevan   PetscFunctionReturn(0);
218ce27a4eeSVijay Mahadevan }
219ce27a4eeSVijay Mahadevan 
220ce27a4eeSVijay Mahadevan #undef __FUNCT__
221ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMoab_GenerateElements_Private"
222*49d66b22SVijay Mahadevan PetscErrorCode DMMoab_GenerateElements_Private(moab::Interface* mbImpl, moab::ReadUtilIface* iface, DMMoabMeshGeneratorCtx& genCtx, PetscInt m, PetscInt n, PetscInt k,
223*49d66b22SVijay Mahadevan     PetscInt a, PetscInt b, PetscInt c, moab::Tag& global_id_tag, moab::EntityHandle startv, moab::Range& cells)
224ce27a4eeSVijay Mahadevan {
225ce27a4eeSVijay Mahadevan   moab::ErrorCode merr;
226ce27a4eeSVijay Mahadevan   PetscInt ix,ie,xe,ye,ze;
227ce27a4eeSVijay Mahadevan   PetscInt ii,jj,kk,nvperelem;
228*49d66b22SVijay Mahadevan   PetscInt simplices_per_tensor[4] = {0,1,2,6};
229ce27a4eeSVijay 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);*/
230ce27a4eeSVijay Mahadevan   PetscInt nelems = ntensorelems;
231ce27a4eeSVijay Mahadevan   moab::EntityHandle starte; // connectivity
232ce27a4eeSVijay Mahadevan   moab::EntityHandle* conn;
233ce27a4eeSVijay Mahadevan 
234ce27a4eeSVijay Mahadevan   PetscFunctionBegin;
235ce27a4eeSVijay Mahadevan   switch(genCtx.dim) {
236ce27a4eeSVijay Mahadevan     case 1:
237ce27a4eeSVijay Mahadevan       nvperelem = 2;
238ce27a4eeSVijay Mahadevan       merr = iface->get_element_connect(nelems, 2, moab::MBEDGE, 0, starte, conn);MBERR("Can't get EDGE2 element connectivity.", merr);
239ce27a4eeSVijay Mahadevan       break;
240ce27a4eeSVijay Mahadevan     case 2:
241ce27a4eeSVijay Mahadevan       if (genCtx.simplex) {
242ce27a4eeSVijay Mahadevan         nvperelem = 3;
243*49d66b22SVijay Mahadevan         nelems = ntensorelems*simplices_per_tensor[genCtx.dim];
244ce27a4eeSVijay Mahadevan         merr = iface->get_element_connect(nelems, 3, moab::MBTRI, 0, starte, conn);MBERR("Can't get TRI3 element connectivity.", merr);
245ce27a4eeSVijay Mahadevan       }
246ce27a4eeSVijay Mahadevan       else {
247ce27a4eeSVijay Mahadevan         nvperelem = 4;
248ce27a4eeSVijay Mahadevan         merr = iface->get_element_connect(nelems, 4, moab::MBQUAD, 0, starte, conn);MBERR("Can't get QUAD4 element connectivity.", merr);
249ce27a4eeSVijay Mahadevan       }
250ce27a4eeSVijay Mahadevan       break;
251ce27a4eeSVijay Mahadevan     case 3:
252ce27a4eeSVijay Mahadevan     default:
253ce27a4eeSVijay Mahadevan       if (genCtx.simplex) {
254ce27a4eeSVijay Mahadevan         nvperelem = 4;
255*49d66b22SVijay Mahadevan         nelems = ntensorelems*simplices_per_tensor[genCtx.dim];
256ce27a4eeSVijay Mahadevan         merr = iface->get_element_connect(nelems, 4, moab::MBTET, 0, starte, conn);MBERR("Can't get TET4 element connectivity.", merr);
257ce27a4eeSVijay Mahadevan       }
258ce27a4eeSVijay Mahadevan       else {
259ce27a4eeSVijay Mahadevan         nvperelem = 8;
260ce27a4eeSVijay Mahadevan         merr = iface->get_element_connect(nelems, 8, moab::MBHEX, 0, starte, conn);MBERR("Can't get HEX8 element connectivity.", merr);
261ce27a4eeSVijay Mahadevan       }
262ce27a4eeSVijay Mahadevan       break;
263ce27a4eeSVijay Mahadevan   }
264ce27a4eeSVijay Mahadevan 
265ce27a4eeSVijay Mahadevan   ix = ie = 0; // index now in the elements, for global ids
266ce27a4eeSVijay Mahadevan 
267ce27a4eeSVijay Mahadevan   /* create a temporary range to store local element handles */
268ce27a4eeSVijay Mahadevan   moab::Range tmp(starte, starte + nelems - 1);
269*49d66b22SVijay Mahadevan   std::vector<PetscInt> gids(nelems);
270ce27a4eeSVijay Mahadevan 
271ce27a4eeSVijay Mahadevan   /* identify the elements at the lower corner, for their global ids */
272ce27a4eeSVijay Mahadevan   xe = m * genCtx.A * genCtx.blockSizeElementXYZ[0] + a * genCtx.blockSizeElementXYZ[0];
273ce27a4eeSVijay Mahadevan   ye = (genCtx.dim > 1 ? n * genCtx.B * genCtx.blockSizeElementXYZ[1] + b * genCtx.blockSizeElementXYZ[1] : 0);
274ce27a4eeSVijay Mahadevan   ze = (genCtx.dim > 2 ? k * genCtx.C * genCtx.blockSizeElementXYZ[2] + c * genCtx.blockSizeElementXYZ[2] : 0);
275ce27a4eeSVijay Mahadevan 
276ce27a4eeSVijay Mahadevan   /* create owned elements requested by genCtx */
277ce27a4eeSVijay Mahadevan   for (kk = 0; kk < (genCtx.dim>2?genCtx.blockSizeElementXYZ[2]:1); kk++) {
278ce27a4eeSVijay Mahadevan     for (jj = 0; jj < (genCtx.dim>1?genCtx.blockSizeElementXYZ[1]:1); jj++) {
279ce27a4eeSVijay Mahadevan       for (ii = 0; ii < genCtx.blockSizeElementXYZ[0]; ii++) {
280ce27a4eeSVijay Mahadevan 
281ce27a4eeSVijay Mahadevan         moab::EntityHandle corner = startv + genCtx.q * ii + genCtx.q * jj * genCtx.ystride + genCtx.q * kk * genCtx.zstride;
282ce27a4eeSVijay Mahadevan 
283*49d66b22SVijay Mahadevan         std::pair<PetscInt,PetscInt> entoffset = DMMoab_SetElementConnectivity_Private(genCtx, ix, corner, conn);
284ce27a4eeSVijay Mahadevan 
285*49d66b22SVijay Mahadevan         for (PetscInt j = 0; j < entoffset.second; j++) {
286ce27a4eeSVijay Mahadevan           /* The entity handle for the particular element -> if we want to set some tags is
287ce27a4eeSVijay Mahadevan              moab::EntityHandle eh = starte + ie + j;
288ce27a4eeSVijay Mahadevan           */
289ce27a4eeSVijay Mahadevan           gids[ie+j] = 1 + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney));
290ce27a4eeSVijay Mahadevan           //gids[ie+j] = ie + j + ((xe + ii) + (ye + jj) * genCtx.nex + (ze + kk) * (genCtx.nex * genCtx.ney)) ;
291ce27a4eeSVijay Mahadevan           //gids[ie+j] = 1 + ie;
292*49d66b22SVijay Mahadevan           //ie++;
293ce27a4eeSVijay Mahadevan         }
294ce27a4eeSVijay Mahadevan 
295ce27a4eeSVijay Mahadevan         ix += entoffset.first;
296*49d66b22SVijay Mahadevan         ie += entoffset.second;
297ce27a4eeSVijay Mahadevan       }
298ce27a4eeSVijay Mahadevan     }
299ce27a4eeSVijay Mahadevan   }
300ce27a4eeSVijay Mahadevan   if (genCtx.adjEnts) { /* we need to update adjacencies now, because some elements are new */
301ce27a4eeSVijay Mahadevan     merr = iface->update_adjacencies(starte, nelems, nvperelem, conn);MBERR("Can't update adjacencies", merr);
302ce27a4eeSVijay Mahadevan   }
303ce27a4eeSVijay Mahadevan   tmp.swap(cells);
304ce27a4eeSVijay Mahadevan   merr = mbImpl->tag_set_data(global_id_tag, cells, &gids[0]);MBERR("Can't set global ids to elements.", merr);
305ce27a4eeSVijay Mahadevan   PetscFunctionReturn(0);
306ce27a4eeSVijay Mahadevan }
307ce27a4eeSVijay Mahadevan 
308ce27a4eeSVijay Mahadevan #undef __FUNCT__
309ce27a4eeSVijay Mahadevan #define __FUNCT__ "DMMBUtil_InitializeOptions"
310ce27a4eeSVijay Mahadevan PetscErrorCode DMMBUtil_InitializeOptions(DMMoabMeshGeneratorCtx& genCtx, PetscInt dim, PetscBool simplex, PetscInt rank, PetscInt nprocs, const PetscReal* bounds, PetscInt nele)
311ce27a4eeSVijay Mahadevan {
312*49d66b22SVijay Mahadevan   PetscInt fraction=0,remainder=0;
313ce27a4eeSVijay Mahadevan   PetscFunctionBegin;
314ce27a4eeSVijay Mahadevan   if(nele<nprocs) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The dimensional discretization size should be greater or equal to number of processors: %D < %D",nele,nprocs);
315ce27a4eeSVijay Mahadevan   /* Initialize all genCtx data */
316ce27a4eeSVijay Mahadevan   genCtx.dim=dim;
317ce27a4eeSVijay Mahadevan   genCtx.simplex=simplex;
318*49d66b22SVijay Mahadevan   genCtx.newMergeMethod=genCtx.keep_skins=genCtx.adjEnts=true;
319ce27a4eeSVijay Mahadevan   /* determine other global quantities for the mesh used for nodes increments */
320ce27a4eeSVijay Mahadevan   genCtx.q = 1;
321ce27a4eeSVijay Mahadevan 
322ce27a4eeSVijay Mahadevan   if (!genCtx.usrxyzgrid) { /* not overridden by genCtx - assume nele equally and that genCtx wants a uniform cube mesh */
323ce27a4eeSVijay Mahadevan     //genCtx.blockSizeElementXYZ[0]=genCtx.blockSizeElementXYZ[1]=genCtx.blockSizeElementXYZ[2]=std::pow(nlocalele,1.0/dim);
324*49d66b22SVijay Mahadevan     if (dim==1) {
325*49d66b22SVijay Mahadevan       fraction=nele/nprocs;    /* partition only by the largest dimension */
326*49d66b22SVijay Mahadevan       remainder=nele%nprocs;   /* remainder after partition which gets evenly distributed by round-robin */
327*49d66b22SVijay Mahadevan       if(rank < remainder)     /* This process gets "fraction+1" elements */
328*49d66b22SVijay Mahadevan         fraction++;
329*49d66b22SVijay Mahadevan     }
330*49d66b22SVijay Mahadevan     genCtx.blockSizeElementXYZ[0]=(dim>1?nele:fraction);
331ce27a4eeSVijay Mahadevan     genCtx.blockSizeElementXYZ[1]=(dim>2?nele:nele/nprocs);
332ce27a4eeSVijay Mahadevan     genCtx.blockSizeElementXYZ[2]=(dim>2?nele/nprocs:1);
333ce27a4eeSVijay Mahadevan   }
334ce27a4eeSVijay Mahadevan 
335ce27a4eeSVijay Mahadevan   /* partition only by the largest dimension */
336ce27a4eeSVijay Mahadevan   /* Total number of local elements := genCtx.blockSizeElementXYZ[0]*(genCtx.dim>1? genCtx.blockSizeElementXYZ[1]*(genCtx.dim>2 ? genCtx.blockSizeElementXYZ[2]:1) :1); */
337ce27a4eeSVijay Mahadevan   if (bounds) {
338*49d66b22SVijay Mahadevan     for (PetscInt i=0; i<6; i++)
339ce27a4eeSVijay Mahadevan       genCtx.xyzbounds[i]=bounds[i];
340ce27a4eeSVijay Mahadevan   }
341ce27a4eeSVijay Mahadevan   else {
342ce27a4eeSVijay Mahadevan     genCtx.xyzbounds[0]=genCtx.xyzbounds[2]=genCtx.xyzbounds[4]=0.0;
343ce27a4eeSVijay Mahadevan     genCtx.xyzbounds[1]=genCtx.xyzbounds[3]=genCtx.xyzbounds[5]=1.0;
344ce27a4eeSVijay Mahadevan   }
345ce27a4eeSVijay Mahadevan 
346ce27a4eeSVijay Mahadevan   if (!genCtx.usrprocgrid) {
347ce27a4eeSVijay Mahadevan     genCtx.M=genCtx.N=genCtx.K=std::pow(nprocs,1.0/dim);
348ce27a4eeSVijay Mahadevan     switch(genCtx.dim) {
349ce27a4eeSVijay Mahadevan      case 1:
350ce27a4eeSVijay Mahadevan        genCtx.M=nprocs;
351ce27a4eeSVijay Mahadevan        genCtx.N=genCtx.K=1;
352ce27a4eeSVijay Mahadevan        break;
353ce27a4eeSVijay Mahadevan      case 2:
354ce27a4eeSVijay Mahadevan        genCtx.K=1;
355ce27a4eeSVijay Mahadevan        genCtx.N=nprocs/(genCtx.M);
356ce27a4eeSVijay Mahadevan        break;
357ce27a4eeSVijay Mahadevan      default:
358ce27a4eeSVijay Mahadevan        genCtx.K=nprocs/(genCtx.M*genCtx.N);
359ce27a4eeSVijay Mahadevan        break;
360ce27a4eeSVijay Mahadevan     }
361ce27a4eeSVijay Mahadevan   }
362ce27a4eeSVijay Mahadevan 
363ce27a4eeSVijay Mahadevan   if (!genCtx.usrrefgrid) {
364ce27a4eeSVijay Mahadevan     genCtx.A=genCtx.B=genCtx.C=1;
365ce27a4eeSVijay Mahadevan   }
366ce27a4eeSVijay Mahadevan 
367ce27a4eeSVijay Mahadevan   /* more default values */
368ce27a4eeSVijay Mahadevan   genCtx.nex=genCtx.ney=genCtx.nez=0;
369ce27a4eeSVijay Mahadevan   genCtx.xstride=genCtx.ystride=genCtx.zstride=0;
370ce27a4eeSVijay Mahadevan   genCtx.NX=genCtx.NY=genCtx.NZ=0;
371ce27a4eeSVijay Mahadevan   genCtx.nex=genCtx.ney=genCtx.nez=0;
372ce27a4eeSVijay Mahadevan   genCtx.blockSizeVertexXYZ[0]=genCtx.blockSizeVertexXYZ[1]=genCtx.blockSizeVertexXYZ[2]=1;
373ce27a4eeSVijay Mahadevan 
374ce27a4eeSVijay Mahadevan   //genCtx.blockSizeElementXYZ[0]=genCtx.blockSizeElementXYZ[1]=genCtx.blockSizeElementXYZ[2]=nele;
375ce27a4eeSVijay Mahadevan   //genCtx.blockSizeVertexXYZ[0]=genCtx.blockSizeVertexXYZ[1]=genCtx.blockSizeVertexXYZ[2]=nele+1;
376ce27a4eeSVijay Mahadevan 
377ce27a4eeSVijay Mahadevan   switch(genCtx.dim) {
378ce27a4eeSVijay Mahadevan    case 3:
379ce27a4eeSVijay Mahadevan      //genCtx.blockSizeElementXYZ[2]=std::max(1,genCtx.blockSizeElementXYZ[2]/nprocs);
380ce27a4eeSVijay Mahadevan      //genCtx.blockSizeElementXYZ[2]=std::max(1,nele/nprocs);
381ce27a4eeSVijay Mahadevan      genCtx.blockSizeVertexXYZ[0] = genCtx.q*genCtx.blockSizeElementXYZ[0] + 1;
382ce27a4eeSVijay Mahadevan      genCtx.blockSizeVertexXYZ[1] = genCtx.q*genCtx.blockSizeElementXYZ[1] + 1;
383ce27a4eeSVijay Mahadevan      //genCtx.blockSizeElementXYZ[2]=fraction/(genCtx.blockSizeElementXYZ[0]*genCtx.blockSizeElementXYZ[1]);
384ce27a4eeSVijay Mahadevan      genCtx.blockSizeVertexXYZ[2] = genCtx.q*genCtx.blockSizeElementXYZ[2] + 1;
385ce27a4eeSVijay Mahadevan 
386ce27a4eeSVijay Mahadevan      genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0];   /* number of elements in x direction, used for global id on element */
387ce27a4eeSVijay Mahadevan      genCtx.dx = (genCtx.xyzbounds[1]-genCtx.xyzbounds[0])/(genCtx.nex*genCtx.q);  /* distance between 2 nodes in x direction */
388ce27a4eeSVijay Mahadevan      genCtx.NX = (genCtx.q * genCtx.nex + 1);
389ce27a4eeSVijay Mahadevan      genCtx.xstride = 1;
390ce27a4eeSVijay Mahadevan      genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1];  /* number of elements in y direction  .... */
391ce27a4eeSVijay Mahadevan      genCtx.dy = (genCtx.xyzbounds[3]-genCtx.xyzbounds[2])/(genCtx.ney*genCtx.q);   /* distance between 2 nodes in y direction */
392ce27a4eeSVijay Mahadevan      genCtx.NY = (genCtx.q * genCtx.ney + 1);
393ce27a4eeSVijay Mahadevan      genCtx.ystride = genCtx.blockSizeVertexXYZ[0];
394ce27a4eeSVijay Mahadevan      genCtx.nez = genCtx.K * genCtx.C * genCtx.blockSizeElementXYZ[2];  /* number of elements in z direction  .... */
395ce27a4eeSVijay Mahadevan      genCtx.dz = (genCtx.xyzbounds[5]-genCtx.xyzbounds[4])/(genCtx.nez*genCtx.q);   /* distance between 2 nodes in z direction */
396ce27a4eeSVijay Mahadevan      genCtx.NZ = (genCtx.q * genCtx.nez + 1);
397ce27a4eeSVijay Mahadevan      genCtx.zstride = genCtx.blockSizeVertexXYZ[0] * genCtx.blockSizeVertexXYZ[1];
398ce27a4eeSVijay Mahadevan      break;
399ce27a4eeSVijay Mahadevan    case 2:
400ce27a4eeSVijay Mahadevan      //genCtx.blockSizeElementXYZ[1]=genCtx.blockSizeElementXYZ[2]=std::max(1,nele/nprocs);
401ce27a4eeSVijay Mahadevan      //genCtx.blockSizeElementXYZ[1]=fraction/genCtx.blockSizeElementXYZ[0];
402ce27a4eeSVijay Mahadevan      genCtx.blockSizeVertexXYZ[0] = genCtx.q*genCtx.blockSizeElementXYZ[0] + 1;
403ce27a4eeSVijay Mahadevan      genCtx.blockSizeVertexXYZ[1] = genCtx.q*genCtx.blockSizeElementXYZ[1] + 1;
404ce27a4eeSVijay Mahadevan      genCtx.blockSizeElementXYZ[2]=0;
405ce27a4eeSVijay Mahadevan 
406ce27a4eeSVijay Mahadevan      genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0];   /* number of elements in x direction, used for global id on element */
407ce27a4eeSVijay Mahadevan      genCtx.dx = (genCtx.xyzbounds[1]-genCtx.xyzbounds[0])/(genCtx.nex*genCtx.q);  /* distance between 2 nodes in x direction */
408ce27a4eeSVijay Mahadevan      genCtx.NX = (genCtx.q * genCtx.nex + 1);
409ce27a4eeSVijay Mahadevan      genCtx.xstride = 1;
410ce27a4eeSVijay Mahadevan      genCtx.ney = genCtx.N * genCtx.B * genCtx.blockSizeElementXYZ[1];  /* number of elements in y direction  .... */
411ce27a4eeSVijay Mahadevan      genCtx.dy = (genCtx.xyzbounds[3]-genCtx.xyzbounds[2])/(genCtx.ney*genCtx.q);   /* distance between 2 nodes in y direction */
412ce27a4eeSVijay Mahadevan      genCtx.NY = (genCtx.q * genCtx.ney + 1);
413ce27a4eeSVijay Mahadevan      genCtx.ystride = genCtx.blockSizeVertexXYZ[0];
414ce27a4eeSVijay Mahadevan      break;
415ce27a4eeSVijay Mahadevan    case 1:
416ce27a4eeSVijay Mahadevan      genCtx.blockSizeElementXYZ[1]=genCtx.blockSizeElementXYZ[2]=0;
417*49d66b22SVijay Mahadevan      //genCtx.blockSizeElementXYZ[0]=nele;
418ce27a4eeSVijay Mahadevan      genCtx.blockSizeVertexXYZ[0] = genCtx.q*genCtx.blockSizeElementXYZ[0] + 1;
419ce27a4eeSVijay Mahadevan 
420ce27a4eeSVijay Mahadevan      genCtx.nex = genCtx.M * genCtx.A * genCtx.blockSizeElementXYZ[0];   /* number of elements in x direction, used for global id on element */
421ce27a4eeSVijay Mahadevan      genCtx.dx = (genCtx.xyzbounds[1]-genCtx.xyzbounds[0])/(genCtx.nex*genCtx.q);  /* distance between 2 nodes in x direction */
422ce27a4eeSVijay Mahadevan      genCtx.NX = (genCtx.q * genCtx.nex + 1);
423ce27a4eeSVijay Mahadevan      genCtx.xstride = 1;
424ce27a4eeSVijay Mahadevan      break;
425ce27a4eeSVijay Mahadevan   }
426ce27a4eeSVijay Mahadevan 
427ce27a4eeSVijay Mahadevan   /*
428ce27a4eeSVijay Mahadevan   PetscInfo3(NULL, "Local elements:= %d, %d, %d\n",genCtx.blockSizeElementXYZ[0],genCtx.blockSizeElementXYZ[1],genCtx.blockSizeElementXYZ[2]);
429ce27a4eeSVijay Mahadevan   PetscInfo3(NULL, "Local vertices:= %d, %d, %d\n",genCtx.blockSizeVertexXYZ[0],genCtx.blockSizeVertexXYZ[1],genCtx.blockSizeVertexXYZ[2]);
430ce27a4eeSVijay Mahadevan   PetscInfo3(NULL, "Local nexyz:= %d, %d, %d\n",genCtx.nex,genCtx.ney,genCtx.nez);
431ce27a4eeSVijay Mahadevan   PetscInfo3(NULL, "Local delxyz:= %g, %g, %g\n",genCtx.dx,genCtx.dy,genCtx.dz);
432ce27a4eeSVijay Mahadevan   PetscInfo3(NULL, "Local strides:= %d, %d, %d\n",genCtx.xstride,genCtx.ystride,genCtx.zstride);
433ce27a4eeSVijay Mahadevan   */
434ce27a4eeSVijay Mahadevan   PetscFunctionReturn(0);
435ce27a4eeSVijay Mahadevan }
436ce27a4eeSVijay Mahadevan 
437aa859218SVijay Mahadevan /*@
438ce27a4eeSVijay Mahadevan   DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with genCtx specified bounds.
439aa859218SVijay Mahadevan 
440aa859218SVijay Mahadevan   Collective on MPI_Comm
441aa859218SVijay Mahadevan 
442aa859218SVijay Mahadevan   Input Parameters:
443aa859218SVijay Mahadevan + comm - The communicator for the DM object
444aa859218SVijay Mahadevan . dim - The spatial dimension
445b8ecf6d3SVijay 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
446b8ecf6d3SVijay Mahadevan . nele - The number of discrete elements in each direction
447b8ecf6d3SVijay Mahadevan . user_nghost - The number of ghosted layers needed in the partitioned mesh
448aa859218SVijay Mahadevan 
449aa859218SVijay Mahadevan   Output Parameter:
450aa859218SVijay Mahadevan . dm  - The DM object
451aa859218SVijay Mahadevan 
452aa859218SVijay Mahadevan   Level: beginner
453aa859218SVijay Mahadevan 
454aa859218SVijay Mahadevan .keywords: DM, create
455aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile()
456aa859218SVijay Mahadevan @*/
457ce27a4eeSVijay Mahadevan PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal* bounds, PetscInt nele, PetscInt nghost, DM *dm)
45851d15aeeSVijay Mahadevan {
45951d15aeeSVijay Mahadevan   PetscErrorCode         ierr;
46051d15aeeSVijay Mahadevan   moab::ErrorCode        merr;
461ce27a4eeSVijay Mahadevan   PetscInt               a,b,c,n,global_size,global_rank=-1;
46251d15aeeSVijay Mahadevan   DM_Moab               *dmmoab;
463ce27a4eeSVijay Mahadevan   moab::Interface       *mbImpl;
46451d15aeeSVijay Mahadevan   moab::ParallelComm    *pcomm;
465a4717116SVijay Mahadevan   moab::ReadUtilIface   *readMeshIface;
466ce27a4eeSVijay Mahadevan   moab::Range            verts,cells,edges,faces,adj,dim3,dim2;
467ce27a4eeSVijay Mahadevan   DMMoabMeshGeneratorCtx genCtx;
468a4717116SVijay Mahadevan   const PetscInt         npts=nele+1;        /* Number of points in every dimension */
469ce27a4eeSVijay Mahadevan 
470ce27a4eeSVijay Mahadevan   moab::Tag              global_id_tag,part_tag,geom_tag;
471ce27a4eeSVijay Mahadevan   moab::Range            ownedvtx,ownedelms,localvtxs,localelms;
472ce27a4eeSVijay Mahadevan   moab::EntityHandle     regionset;
473ce27a4eeSVijay Mahadevan   PetscInt               ml=0,nl=0,kl=0,leftover;
47451d15aeeSVijay Mahadevan 
47551d15aeeSVijay Mahadevan   PetscFunctionBegin;
476e5136372SVijay Mahadevan   if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n");
477e5136372SVijay Mahadevan 
478ce27a4eeSVijay Mahadevan   ierr = MPI_Comm_size(comm, &global_size);CHKERRQ(ierr);
479e5136372SVijay Mahadevan   /* total number of vertices in all dimensions */
480e5136372SVijay Mahadevan   n=pow(npts,dim);
481e5136372SVijay Mahadevan 
482e5136372SVijay Mahadevan   /* do some error checking */
483e5136372SVijay Mahadevan   if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n");
484ce27a4eeSVijay Mahadevan   if(global_size > n) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of processors must be less than or equal to number of elements.\n");
485e5136372SVijay Mahadevan   if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n");
486e5136372SVijay Mahadevan 
487a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
488c528d872SBarry Smith   ierr = DMMoabCreateMoab(comm, NULL, NULL, NULL, NULL, dm);CHKERRQ(ierr);
48951d15aeeSVijay Mahadevan 
490a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
49151d15aeeSVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
492ce27a4eeSVijay Mahadevan   mbImpl = dmmoab->mbiface;
49351d15aeeSVijay Mahadevan   pcomm = dmmoab->pcomm;
494ce27a4eeSVijay Mahadevan   global_id_tag = dmmoab->ltog_tag;
495ce27a4eeSVijay Mahadevan   global_rank = pcomm->rank();
496c8d5365dSVijay Mahadevan   dmmoab->dim = dim;
497*49d66b22SVijay Mahadevan   dmmoab->nghostrings=nghost;
49851d15aeeSVijay Mahadevan 
499b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
500ce27a4eeSVijay Mahadevan   merr = mbImpl->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
501b5410836SVijay Mahadevan 
502a4717116SVijay Mahadevan   /* No errors yet; proceed with building the mesh */
503ce27a4eeSVijay Mahadevan   merr = mbImpl->query_interface(readMeshIface);MBERRNM(merr);
50451d15aeeSVijay Mahadevan 
505ce27a4eeSVijay Mahadevan   //PetscObjectOptionsBegin(dm);
506ce27a4eeSVijay Mahadevan   PetscOptionsBegin(comm,"","DMMoab Creation Options","DMMOAB");
507ce27a4eeSVijay Mahadevan   /* Handle DMMoab spatial resolution */
508ce27a4eeSVijay Mahadevan   PetscOptionsInt("-dmb_grid_x","Number of grid points in x direction","DMMoabSetSizes",genCtx.blockSizeElementXYZ[0],&genCtx.blockSizeElementXYZ[0],&genCtx.usrxyzgrid);
509ce27a4eeSVijay Mahadevan   if (dim > 1) {PetscOptionsInt("-dmb_grid_y","Number of grid points in y direction","DMMoabSetSizes",genCtx.blockSizeElementXYZ[1],&genCtx.blockSizeElementXYZ[1],&genCtx.usrxyzgrid);}
510ce27a4eeSVijay Mahadevan   if (dim > 2) {PetscOptionsInt("-dmb_grid_z","Number of grid points in z direction","DMMoabSetSizes",genCtx.blockSizeElementXYZ[2],&genCtx.blockSizeElementXYZ[2],&genCtx.usrxyzgrid);}
511a4717116SVijay Mahadevan 
512ce27a4eeSVijay Mahadevan   /* Handle DMMoab parallel distibution */
513ce27a4eeSVijay Mahadevan   PetscOptionsInt("-dmb_processors_x","Number of processors in x direction","DMMoabSetNumProcs",genCtx.M,&genCtx.M,&genCtx.usrprocgrid);
514ce27a4eeSVijay Mahadevan   if (dim > 1) {PetscOptionsInt("-dmb_processors_y","Number of processors in y direction","DMMoabSetNumProcs",genCtx.N,&genCtx.N,&genCtx.usrprocgrid);}
515ce27a4eeSVijay Mahadevan   if (dim > 2) {PetscOptionsInt("-dmb_processors_z","Number of processors in z direction","DMMoabSetNumProcs",genCtx.K,&genCtx.K,&genCtx.usrprocgrid);}
51651d15aeeSVijay Mahadevan 
517ce27a4eeSVijay Mahadevan   /* Handle DMMoab block refinement */
518ce27a4eeSVijay Mahadevan   PetscOptionsInt("-dmb_refine_x","Number of refinement blocks in x direction","DMMoabSetRefinement",genCtx.A,&genCtx.A,&genCtx.usrrefgrid);
519ce27a4eeSVijay Mahadevan   if (dim > 1) {PetscOptionsInt("-dmb_refine_y","Number of refinement blocks in y direction","DMMoabSetRefinement",genCtx.B,&genCtx.B,&genCtx.usrrefgrid);}
520ce27a4eeSVijay Mahadevan   if (dim > 2) {PetscOptionsInt("-dmb_refine_z","Number of refinement blocks in z direction","DMMoabSetRefinement",genCtx.C,&genCtx.C,&genCtx.usrrefgrid);}
521ce27a4eeSVijay Mahadevan   PetscOptionsEnd();
52251d15aeeSVijay Mahadevan 
523ce27a4eeSVijay Mahadevan   ierr = DMMBUtil_InitializeOptions(genCtx, dim, useSimplex, global_rank, global_size, bounds, nele);CHKERRQ(ierr);
52451d15aeeSVijay Mahadevan 
525ce27a4eeSVijay Mahadevan   /* Lets check for some valid input */
526ce27a4eeSVijay Mahadevan   if (genCtx.dim < 1 || genCtx.dim > 3) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid topological dimension specified: %d.\n",genCtx.dim);
527ce27a4eeSVijay Mahadevan   if (genCtx.M * genCtx.N * genCtx.K != global_size) SETERRQ4(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid [m, n, k] data: %d, %d, %d. Product must be equal to global size = %d.\n",genCtx.M,genCtx.N,genCtx.K,global_size);
528ce27a4eeSVijay Mahadevan   if (genCtx.xyzbounds) {
52966f88a78SVijay Mahadevan     /* validate the bounds data */
530ce27a4eeSVijay Mahadevan     if(genCtx.xyzbounds[0] >= genCtx.xyzbounds[1]) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"X-dim: Left boundary cannot be greater than right. [%G >= %G]\n",genCtx.xyzbounds[0],genCtx.xyzbounds[1]);
531ce27a4eeSVijay Mahadevan     if(genCtx.dim > 1 && (genCtx.xyzbounds[2] >= genCtx.xyzbounds[3])) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Y-dim: Left boundary cannot be greater than right. [%G >= %G]\n",genCtx.xyzbounds[2],genCtx.xyzbounds[3]);
532ce27a4eeSVijay Mahadevan     if(genCtx.dim > 2 && (genCtx.xyzbounds[4] >= genCtx.xyzbounds[5])) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Z-dim: Left boundary cannot be greater than right. [%G >= %G]\n",genCtx.xyzbounds[4],genCtx.xyzbounds[5]);
53366f88a78SVijay Mahadevan   }
534ce27a4eeSVijay Mahadevan   if (genCtx.adjEnts) genCtx.keep_skins = true; /* do not delete anything - consumes more memory */
53566f88a78SVijay Mahadevan 
536ce27a4eeSVijay Mahadevan   /* determine m, n, k for processor rank */
537ce27a4eeSVijay Mahadevan   kl = (genCtx.dim>2?global_rank/(genCtx.M*genCtx.N):0);
538ce27a4eeSVijay Mahadevan   leftover = global_rank%(genCtx.M*genCtx.N);
539ce27a4eeSVijay Mahadevan   nl = leftover/genCtx.M;
540ce27a4eeSVijay Mahadevan   ml = leftover%genCtx.M;
541e5136372SVijay Mahadevan 
542ce27a4eeSVijay Mahadevan   /*
543ce27a4eeSVijay Mahadevan    * so there are a total of M * A * blockSizeElement elements in x direction (so M * A * blockSizeElement + 1 verts in x direction)
544ce27a4eeSVijay Mahadevan    * so there are a total of N * B * blockSizeElement elements in y direction (so N * B * blockSizeElement + 1 verts in y direction)
545ce27a4eeSVijay Mahadevan    * so there are a total of K * C * blockSizeElement elements in z direction (so K * C * blockSizeElement + 1 verts in z direction)
546c8d5365dSVijay Mahadevan 
547ce27a4eeSVijay Mahadevan    * there are ( M * A blockSizeElement )      *  ( N * B * blockSizeElement)      * (K * C * blockSizeElement )    hexas
548ce27a4eeSVijay Mahadevan    * there are ( M * A * blockSizeElement + 1) *  ( N * B * blockSizeElement + 1 ) * (K * C * blockSizeElement + 1) vertices
549ce27a4eeSVijay Mahadevan    * x is the first dimension that varies
550ce27a4eeSVijay Mahadevan    */
551e5136372SVijay Mahadevan 
552ce27a4eeSVijay Mahadevan   /* generate the block at (a, b, c); it will represent a partition , it will get a partition tag */
553*49d66b22SVijay Mahadevan   PetscInt dum_id=-1;
554ce27a4eeSVijay Mahadevan   merr = mbImpl->tag_get_handle("GLOBAL_ID", 1, moab::MB_TYPE_INTEGER, global_id_tag);MBERR("Getting Tag handle failed", merr);
55551d15aeeSVijay Mahadevan 
556ce27a4eeSVijay Mahadevan   merr = mbImpl->tag_get_handle("PARALLEL_PARTITION", 1, moab::MB_TYPE_INTEGER, part_tag, moab::MB_TAG_CREAT|moab::MB_TAG_SPARSE, &dum_id);MBERR("Getting Tag handle failed", merr);
55751d15aeeSVijay Mahadevan 
5583a4aeca1SVijay Mahadevan   /* lets create some sets */
559ce27a4eeSVijay Mahadevan   merr = mbImpl->tag_get_handle(GEOM_DIMENSION_TAG_NAME, 1, moab::MB_TYPE_INTEGER, geom_tag, moab::MB_TAG_CREAT|moab::MB_TAG_SPARSE, &dum_id);MBERRNM(merr);
560*49d66b22SVijay Mahadevan   merr = mbImpl->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr);
5613a4aeca1SVijay Mahadevan 
562ce27a4eeSVijay Mahadevan   for (a=0; a<(genCtx.dim>0?genCtx.A:genCtx.A); a++) {
563ce27a4eeSVijay Mahadevan     for (b=0; b<(genCtx.dim>1?genCtx.B:1); b++) {
564ce27a4eeSVijay Mahadevan       for (c=0; c<(genCtx.dim>2?genCtx.C:1); c++) {
565*49d66b22SVijay Mahadevan 
566*49d66b22SVijay Mahadevan         moab::EntityHandle startv;
567*49d66b22SVijay Mahadevan 
568ce27a4eeSVijay Mahadevan         ierr = DMMoab_GenerateVertices_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, verts);CHKERRQ(ierr);
5693a4aeca1SVijay Mahadevan 
570ce27a4eeSVijay Mahadevan         ierr = DMMoab_GenerateElements_Private(mbImpl, readMeshIface, genCtx, ml, nl, kl, a, b, c, global_id_tag, startv, cells);CHKERRQ(ierr);
5713a4aeca1SVijay Mahadevan 
572*49d66b22SVijay Mahadevan         PetscInt part_num=0;
573*49d66b22SVijay Mahadevan         switch(genCtx.dim) {
574*49d66b22SVijay Mahadevan           case 3:
575*49d66b22SVijay Mahadevan             part_num += (c+kl*genCtx.C)*(genCtx.M*genCtx.A * genCtx.N*genCtx.B);
576*49d66b22SVijay Mahadevan           case 2:
577*49d66b22SVijay Mahadevan             part_num += (b+nl*genCtx.B)*(genCtx.M*genCtx.A);
578*49d66b22SVijay Mahadevan           case 1:
579*49d66b22SVijay Mahadevan             part_num += (a+ml*genCtx.A);
580*49d66b22SVijay Mahadevan            break;
581*49d66b22SVijay Mahadevan         }
582*49d66b22SVijay Mahadevan 
583*49d66b22SVijay Mahadevan         moab::EntityHandle part_set;
584ce27a4eeSVijay Mahadevan         merr = mbImpl->create_meshset(moab::MESHSET_SET, part_set);MBERR("Can't create mesh set.", merr);
5853a4aeca1SVijay Mahadevan 
586ce27a4eeSVijay Mahadevan         merr = mbImpl->add_entities(part_set, verts);MBERR("Can't add vertices to set.", merr);
587*49d66b22SVijay Mahadevan         merr = mbImpl->add_entities(part_set, cells);MBERR("Can't add entities to set.", merr);
588*49d66b22SVijay Mahadevan         // PetscInfo2(NULL, "Generated local mesh: %D vertices and %D elements.\n", verts.size(), cells.size());
5893a4aeca1SVijay Mahadevan 
590ce27a4eeSVijay Mahadevan         merr = mbImpl->add_entities(regionset, cells);MBERR("Can't add entities to set.", merr);
591ce27a4eeSVijay Mahadevan 
592ce27a4eeSVijay Mahadevan         /* if needed, add all edges and faces */
593ce27a4eeSVijay Mahadevan         if (genCtx.adjEnts)
594ce27a4eeSVijay Mahadevan         {
595ce27a4eeSVijay Mahadevan           if (genCtx.dim > 1) {
596ce27a4eeSVijay Mahadevan             merr = mbImpl->get_adjacencies(cells, 1, true, edges, moab::Interface::UNION);MBERR("Can't get edges", merr);
597ce27a4eeSVijay Mahadevan             merr = mbImpl->add_entities(part_set, edges);MBERR("Can't add edges to partition set.", merr);
598ce27a4eeSVijay Mahadevan           }
599ce27a4eeSVijay Mahadevan           if (genCtx.dim > 2) {
600ce27a4eeSVijay Mahadevan             merr = mbImpl->get_adjacencies(cells, 2, true, faces, moab::Interface::UNION);MBERR("Can't get faces", merr);
601ce27a4eeSVijay Mahadevan             merr = mbImpl->add_entities(part_set, faces);MBERR("Can't add faces to partition set.", merr);
602ce27a4eeSVijay Mahadevan           }
603ce27a4eeSVijay Mahadevan           edges.clear();
604ce27a4eeSVijay Mahadevan           faces.clear();
605ce27a4eeSVijay Mahadevan         }
606*49d66b22SVijay Mahadevan         verts.clear(); cells.clear();
607ce27a4eeSVijay Mahadevan 
608ce27a4eeSVijay Mahadevan         merr = mbImpl->tag_set_data(part_tag, &part_set, 1, &part_num);MBERR("Can't set part tag on set", merr);
609*49d66b22SVijay Mahadevan         if (dmmoab->fileset) {
610ce27a4eeSVijay Mahadevan           merr = mbImpl->add_parent_child(dmmoab->fileset, part_set);MBERR("Can't add part set to file set.", merr);
611ce27a4eeSVijay Mahadevan           merr = mbImpl->unite_meshset(dmmoab->fileset, part_set);MBERRNM(merr);
612*49d66b22SVijay Mahadevan         }
613*49d66b22SVijay Mahadevan         merr = mbImpl->add_entities(dmmoab->fileset, &part_set, 1);MBERRNM(merr);
614ce27a4eeSVijay Mahadevan       }
615ce27a4eeSVijay Mahadevan     }
616ce27a4eeSVijay Mahadevan   }
617ce27a4eeSVijay Mahadevan 
618*49d66b22SVijay Mahadevan   /* set geometric dimension tag for regions */
619*49d66b22SVijay Mahadevan   merr = mbImpl->tag_set_data(geom_tag, &regionset, 1, &dmmoab->dim);MBERRNM(merr);
620*49d66b22SVijay Mahadevan   merr = mbImpl->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr);
621ce27a4eeSVijay Mahadevan 
622*49d66b22SVijay Mahadevan   /* Only in parallel: resolve shared entities between processors and exchange ghost layers */
623ce27a4eeSVijay Mahadevan   if (global_size>1) {
624ce27a4eeSVijay Mahadevan 
625ce27a4eeSVijay Mahadevan     merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, genCtx.dim, cells);MBERR("Can't get all d-dimensional elements.", merr);
626ce27a4eeSVijay Mahadevan     merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 0, verts);MBERR("Can't get all vertices.", merr);
627ce27a4eeSVijay Mahadevan 
628ce27a4eeSVijay Mahadevan     if (genCtx.A*genCtx.B*genCtx.C!=1) { //  merge needed
629*49d66b22SVijay Mahadevan       moab::MergeMesh mm(mbImpl);
630ce27a4eeSVijay Mahadevan       if (genCtx.newMergeMethod) {
631ce27a4eeSVijay Mahadevan         merr = mm.merge_using_integer_tag(verts, global_id_tag);MBERR("Can't merge with GLOBAL_ID tag", merr);
632ce27a4eeSVijay Mahadevan       }
633ce27a4eeSVijay Mahadevan       else {
634ce27a4eeSVijay Mahadevan         merr = mm.merge_entities(cells, 0.0001);MBERR("Can't merge with coordinates", merr);
6353a4aeca1SVijay Mahadevan       }
6363a4aeca1SVijay Mahadevan     }
6373a4aeca1SVijay Mahadevan 
638a4717116SVijay Mahadevan     /* check the handles */
639ce27a4eeSVijay Mahadevan     merr = pcomm->check_all_shared_handles();MBERRV(mbImpl,merr);
64051d15aeeSVijay Mahadevan 
641a4717116SVijay Mahadevan     /* resolve the shared entities by exchanging information to adjacent processors */
642ce27a4eeSVijay Mahadevan     merr = pcomm->resolve_shared_ents(dmmoab->fileset,cells,dim,dim-1,NULL,&global_id_tag);MBERRV(mbImpl,merr);
643*49d66b22SVijay Mahadevan     if (dmmoab->fileset) {
644ce27a4eeSVijay Mahadevan       merr = pcomm->exchange_ghost_cells(dim,0,nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbImpl,merr);
645*49d66b22SVijay Mahadevan     }
646*49d66b22SVijay Mahadevan     else {
647*49d66b22SVijay Mahadevan       merr = pcomm->exchange_ghost_cells(dim,0,nghost,dim,true,false);MBERRV(mbImpl,merr);
648*49d66b22SVijay Mahadevan     }
649c8d5365dSVijay Mahadevan 
650e427d9c9SVijay Mahadevan     /* Reassign global IDs on all entities. */
651e427d9c9SVijay Mahadevan     merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true,false);MBERRNM(merr);
652*49d66b22SVijay Mahadevan   }
6533f1c6e43SVijay Mahadevan 
654ce27a4eeSVijay Mahadevan   if (!genCtx.keep_skins) { // default is to delete the 1- and 2-dimensional entities
655ce27a4eeSVijay Mahadevan     // delete all quads and edges
656ce27a4eeSVijay Mahadevan     moab::Range toDelete;
657ce27a4eeSVijay Mahadevan     if (genCtx.dim > 1) {
658ce27a4eeSVijay Mahadevan       merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 1, toDelete);MBERR("Can't get edges", merr);
659ce27a4eeSVijay Mahadevan     }
6603f1c6e43SVijay Mahadevan 
661ce27a4eeSVijay Mahadevan     if (genCtx.dim > 2) {
662ce27a4eeSVijay Mahadevan       merr = mbImpl->get_entities_by_dimension(dmmoab->fileset, 2, toDelete);MBERR("Can't get faces", merr);
663ce27a4eeSVijay Mahadevan     }
664c8d5365dSVijay Mahadevan 
665ce27a4eeSVijay Mahadevan     merr = dmmoab->pcomm->delete_entities(toDelete) ;MBERR("Can't delete entities", merr);
666ce27a4eeSVijay Mahadevan   }
66751d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
66851d15aeeSVijay Mahadevan }
66951d15aeeSVijay Mahadevan 
67051d15aeeSVijay Mahadevan 
671*49d66b22SVijay Mahadevan #undef __FUNCT__
672*49d66b22SVijay Mahadevan #define __FUNCT__ "DMMoab_GetReadOptions_Private"
673*49d66b22SVijay Mahadevan PetscErrorCode DMMoab_GetReadOptions_Private(PetscBool by_rank, PetscInt numproc, PetscInt dim, PetscInt nghost, MoabReadMode mode, PetscInt dbglevel, const char* dm_opts, const char* extra_opts, const char** read_opts)
67451d15aeeSVijay Mahadevan {
6756acfe860SVijay Mahadevan   char           *ropts;
676*49d66b22SVijay Mahadevan   char           ropts_par[PETSC_MAX_PATH_LEN],ropts_pargh[PETSC_MAX_PATH_LEN];
67761eb6e27SVijay Mahadevan   char           ropts_dbg[PETSC_MAX_PATH_LEN];
678f90c3b0eSVijay Mahadevan   PetscErrorCode ierr;
67951d15aeeSVijay Mahadevan 
680a4717116SVijay Mahadevan   PetscFunctionBegin;
68195dccacaSBarry Smith   ierr = PetscMalloc1(PETSC_MAX_PATH_LEN,&ropts);CHKERRQ(ierr);
68261eb6e27SVijay Mahadevan   ierr = PetscMemzero(&ropts_par,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
683*49d66b22SVijay Mahadevan   ierr = PetscMemzero(&ropts_pargh,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
68461eb6e27SVijay Mahadevan   ierr = PetscMemzero(&ropts_dbg,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
68561eb6e27SVijay Mahadevan 
686e23c60ebSVijay Mahadevan   /* do parallel read unless using only one processor */
687a4717116SVijay Mahadevan   if (numproc > 1) {
688*49d66b22SVijay Mahadevan     // ierr = PetscSNPrintf(ropts_par, PETSC_MAX_PATH_LEN, "PARALLEL=%s;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=%d.0.1%s;",MoabReadModes[mode],dim,(by_rank ? ";PARTITION_BY_RANK":""));CHKERRQ(ierr);
689*49d66b22SVijay Mahadevan     ierr = PetscSNPrintf(ropts_par, PETSC_MAX_PATH_LEN, "PARALLEL=%s;PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;PARALLEL_RESOLVE_SHARED_ENTS;%s",MoabReadModes[mode],(by_rank ? "PARTITION_BY_RANK;":""));CHKERRQ(ierr);
690*49d66b22SVijay Mahadevan     if (nghost) {
691*49d66b22SVijay Mahadevan       ierr = PetscSNPrintf(ropts_pargh, PETSC_MAX_PATH_LEN, "PARALLEL_GHOSTS=%d.0.%d;",dim,nghost);CHKERRQ(ierr);
692*49d66b22SVijay Mahadevan     }
6932e4e7c01SVijay Mahadevan   }
6942e4e7c01SVijay Mahadevan 
695c8d5365dSVijay Mahadevan   if (dbglevel) {
696f90c3b0eSVijay Mahadevan     if (numproc>1) {
697*49d66b22SVijay Mahadevan       ierr = PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%d;DEBUG_PIO=%d;",dbglevel,dbglevel);CHKERRQ(ierr);
69861eb6e27SVijay Mahadevan     }
69961eb6e27SVijay Mahadevan     else {
700*49d66b22SVijay Mahadevan       ierr = PetscSNPrintf(ropts_dbg, PETSC_MAX_PATH_LEN, "CPUTIME;DEBUG_IO=%d;",dbglevel);CHKERRQ(ierr);
701f90c3b0eSVijay Mahadevan     }
702c8d5365dSVijay Mahadevan   }
70351d15aeeSVijay Mahadevan 
704*49d66b22SVijay Mahadevan   ierr = PetscSNPrintf(ropts, PETSC_MAX_PATH_LEN, "%s%s%s%s%s",ropts_par,(nghost?ropts_pargh:""),ropts_dbg,(extra_opts?extra_opts:""),(dm_opts?dm_opts:""));CHKERRQ(ierr);
705f90c3b0eSVijay Mahadevan   *read_opts = ropts;
70651d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
70751d15aeeSVijay Mahadevan }
70851d15aeeSVijay Mahadevan 
70951d15aeeSVijay Mahadevan 
710aa859218SVijay Mahadevan /*@
711aa859218SVijay Mahadevan   DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file.
712aa859218SVijay Mahadevan 
713aa859218SVijay Mahadevan   Collective on MPI_Comm
714aa859218SVijay Mahadevan 
715aa859218SVijay Mahadevan   Input Parameters:
716aa859218SVijay Mahadevan + comm - The communicator for the DM object
717aa859218SVijay Mahadevan . dim - The spatial dimension
718b8ecf6d3SVijay Mahadevan . filename - The name of the mesh file to be loaded
719b8ecf6d3SVijay Mahadevan . usrreadopts - The options string to read a MOAB mesh.
720aa859218SVijay Mahadevan   Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo)
721aa859218SVijay Mahadevan 
722aa859218SVijay Mahadevan   Output Parameter:
723aa859218SVijay Mahadevan . dm  - The DM object
724aa859218SVijay Mahadevan 
725aa859218SVijay Mahadevan   Level: beginner
726aa859218SVijay Mahadevan 
727aa859218SVijay Mahadevan .keywords: DM, create
728b8ecf6d3SVijay Mahadevan 
729aa859218SVijay Mahadevan .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh()
730aa859218SVijay Mahadevan @*/
731*49d66b22SVijay Mahadevan PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,PetscInt nghost,const char* filename,const char* usrreadopts,DM *dm)
73251d15aeeSVijay Mahadevan {
733a4717116SVijay Mahadevan   moab::ErrorCode     merr;
7342e4e7c01SVijay Mahadevan   PetscInt            nprocs;
735a4717116SVijay Mahadevan   DM_Moab            *dmmoab;
736a4717116SVijay Mahadevan   moab::Interface    *mbiface;
737a4717116SVijay Mahadevan   moab::ParallelComm *pcomm;
738a4717116SVijay Mahadevan   moab::Range         verts,elems;
7397023aa44SVijay Mahadevan   const char         *readopts;
740a4717116SVijay Mahadevan   PetscErrorCode      ierr;
74151d15aeeSVijay Mahadevan 
74251d15aeeSVijay Mahadevan   PetscFunctionBegin;
743*49d66b22SVijay Mahadevan   PetscValidPointer(dm,6);
74451d15aeeSVijay Mahadevan 
745a4717116SVijay Mahadevan   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
746c528d872SBarry Smith   ierr = DMMoabCreateMoab(comm, NULL, NULL, NULL, NULL, dm);CHKERRQ(ierr);
74751d15aeeSVijay Mahadevan 
748a4717116SVijay Mahadevan   /* get all the necessary handles from the private DM object */
749a4717116SVijay Mahadevan   dmmoab = (DM_Moab*)(*dm)->data;
750a4717116SVijay Mahadevan   mbiface = dmmoab->mbiface;
751a4717116SVijay Mahadevan   pcomm = dmmoab->pcomm;
752a4717116SVijay Mahadevan   nprocs = pcomm->size();
753aa859218SVijay Mahadevan   /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */
754c8d5365dSVijay Mahadevan   dmmoab->dim = dim;
755*49d66b22SVijay Mahadevan   dmmoab->nghostrings=nghost;
756a4717116SVijay Mahadevan 
757b5410836SVijay Mahadevan   /* create a file set to associate all entities in current mesh */
758b5410836SVijay Mahadevan   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
759b5410836SVijay Mahadevan 
760a4717116SVijay Mahadevan   /* add mesh loading options specific to the DM */
761*49d66b22SVijay Mahadevan   ierr = DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, nghost, dmmoab->read_mode,
7622e4e7c01SVijay Mahadevan                                         dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts);CHKERRQ(ierr);
7637023aa44SVijay Mahadevan 
764e5136372SVijay Mahadevan   PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts);
765a4717116SVijay Mahadevan 
766a4717116SVijay Mahadevan   /* Load the mesh from a file. */
767*49d66b22SVijay Mahadevan   if (dmmoab->fileset) {
7687023aa44SVijay Mahadevan     merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr);
769*49d66b22SVijay Mahadevan   }
770*49d66b22SVijay Mahadevan   else {
771*49d66b22SVijay Mahadevan     merr = mbiface->load_file(filename, 0, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr);
772*49d66b22SVijay Mahadevan   }
773a4717116SVijay Mahadevan 
7746e40195eSVijay Mahadevan   /* Reassign global IDs on all entities. */
775e5136372SVijay Mahadevan   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr);
776e5136372SVijay Mahadevan 
777e5136372SVijay Mahadevan   /* load the local vertices */
778e5136372SVijay Mahadevan   merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr);
779e5136372SVijay Mahadevan   /* load the local elements */
780e5136372SVijay Mahadevan   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr);
781e5136372SVijay Mahadevan 
782e5136372SVijay Mahadevan   /* Everything is set up, now just do a tag exchange to update tags
783e5136372SVijay Mahadevan      on all of the ghost vertexes */
784e5136372SVijay Mahadevan   merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr);
785e5136372SVijay Mahadevan   merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr);
786a4717116SVijay Mahadevan   merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr);
787a4717116SVijay Mahadevan 
788e5136372SVijay Mahadevan   PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size());
7896acfe860SVijay Mahadevan   ierr = PetscFree(readopts);CHKERRQ(ierr);
79051d15aeeSVijay Mahadevan   PetscFunctionReturn(0);
79151d15aeeSVijay Mahadevan }
79251d15aeeSVijay Mahadevan 
793