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