xref: /petsc/src/dm/impls/moab/dmmbutil.cxx (revision 061e9cde3ff1b999f985dc08fa89b32a32e267f2)
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/ScdInterface.hpp>
8 #include <moab/CN.hpp>
9 
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "DMMoabComputeDomainBounds_Private"
13 static PetscErrorCode DMMoabComputeDomainBounds_Private(moab::ParallelComm* pcomm, PetscInt dim, PetscInt neleglob, PetscInt *ise)
14 {
15   PetscInt size,rank;
16   PetscInt fraction,remainder;
17   PetscInt starts[3],sizes[3];
18 
19   PetscFunctionBegin;
20   if(dim<1 && dim>3) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The problem dimension is invalid: %D",dim);
21 
22   size=pcomm->size();
23   rank=pcomm->rank();
24   fraction=neleglob/size;    /* partition only by the largest dimension */
25   remainder=neleglob%size;   /* remainder after partition which gets evenly distributed by round-robin */
26 
27   if(fraction==0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"The leading dimension size should be greater than number of processors: %D < %D",neleglob,size);
28 
29   starts[0]=starts[1]=starts[2]=0;       /* default dimensional element = 1 */
30   sizes[0]=sizes[1]=sizes[2]=neleglob;   /* default dimensional element = 1 */
31 
32   if(rank < remainder) { /* This process gets "fraction+1" elements */
33     sizes[dim-1] = (fraction + 1);
34     starts[dim-1] = rank * (fraction+1);
35   } else { /* This process gets "fraction" elements */
36     sizes[dim-1] = fraction;
37     starts[dim-1] = (remainder*(fraction+1) + fraction*(rank-remainder));
38   }
39 
40   for(int i=dim-1; i>=0; --i) {
41     ise[2*i]=starts[i];ise[2*i+1]=starts[i]+sizes[i];
42   }
43   PetscFunctionReturn(0);
44 }
45 
46 #undef __FUNCT__
47 #define __FUNCT__ "DMMoab_SetStructuredCoords_Private"
48 static void DMMoab_SetStructuredCoords_Private(PetscInt i, PetscInt j, PetscInt k, PetscReal hx, PetscReal hy, PetscReal hz, PetscInt vcount, std::vector<double*>& vcoords)
49 {
50   vcoords[0][vcount] = i*hx;
51   vcoords[1][vcount] = j*hy;
52   vcoords[2][vcount] = k*hz;
53 }
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "DMMoab_SetTensorElementConnectivity_Private"
57 static void DMMoab_SetTensorElementConnectivity_Private(PetscInt dim, moab::EntityType etype, PetscInt offset, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity)
58 {
59   std::vector<int>    subent_conn(pow(2,dim));  /* only linear edge, quad, hex supported now */
60 
61   moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data());
62 
63   switch(dim) {
64     case 1:
65       connectivity[offset+subent_conn[0]] = vfirst+i;
66       connectivity[offset+subent_conn[1]] = vfirst+(i+1);
67       break;
68     case 2:
69       connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1);
70       connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1);
71       connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1);
72       connectivity[offset+subent_conn[3]] = vfirst+i+(j+1)*(nele+1);
73       break;
74     case 3:
75     default:
76       connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
77       connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k);
78       connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
79       connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k);
80       connectivity[offset+subent_conn[4]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1));
81       connectivity[offset+subent_conn[5]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
82       connectivity[offset+subent_conn[6]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1));
83       connectivity[offset+subent_conn[7]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
84       break;
85   }
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "DMMoab_SetSimplexElementConnectivity_Private"
90 static void DMMoab_SetSimplexElementConnectivity_Private(PetscInt dim, PetscInt subelem, moab::EntityType etype, PetscInt offset, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity)
91 {
92   std::vector<int>    subent_conn(pow(2,dim));  /* only linear edge, quad, hex supported now */
93 
94   moab::CN::SubEntityVertexIndices(etype, dim, 0, subent_conn.data());
95 
96   switch(dim) {
97     case 1:
98       connectivity[offset+subent_conn[0]] = vfirst+i;
99       connectivity[offset+subent_conn[1]] = vfirst+(i+1);
100       break;
101     case 2:
102       if (subelem) { /* 1 2 3 of a QUAD */
103         connectivity[offset+subent_conn[0]] = vfirst+i+j*(nele+1);
104         connectivity[offset+subent_conn[1]] = vfirst+(i+1)+j*(nele+1);
105         connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(j+1)*(nele+1);
106       }
107       else {        /* 3 4 1 of a QUAD */
108         connectivity[offset+subent_conn[0]] = vfirst+(i+1)+(j+1)*(nele+1);
109         connectivity[offset+subent_conn[1]] = vfirst+i+(j+1)*(nele+1);
110         connectivity[offset+subent_conn[2]] = vfirst+i+j*(nele+1);
111       }
112       break;
113     case 3:
114     default:
115       switch(subelem) {
116         case 0: /* 0 1 2 5 of a HEX */
117           connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
118           connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*k);
119           connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
120           connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
121           break;
122         case 1: /* 0 2 7 5 of a HEX */
123           connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
124           connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
125           connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
126           connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
127           break;
128         case 2: /* 0 2 3 7 of a HEX */
129           connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
130           connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
131           connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*k);
132           connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
133           break;
134         case 3: /* 0 5 7 4 of a HEX */
135           connectivity[offset+subent_conn[0]] = vfirst+i+(nele+1)*(j+(nele+1)*k);
136           connectivity[offset+subent_conn[1]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
137           connectivity[offset+subent_conn[2]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
138           connectivity[offset+subent_conn[3]] = vfirst+i+(nele+1)*(j+(nele+1)*(k+1));
139           break;
140         case 4: /* 2 7 5 6 of a HEX */
141           connectivity[offset+subent_conn[0]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*k);
142           connectivity[offset+subent_conn[1]] = vfirst+i+(nele+1)*((j+1)+(nele+1)*(k+1));
143           connectivity[offset+subent_conn[2]] = vfirst+(i+1)+(nele+1)*(j+(nele+1)*(k+1));
144           connectivity[offset+subent_conn[3]] = vfirst+(i+1)+(nele+1)*((j+1)+(nele+1)*(k+1));
145           break;
146       }
147       break;
148   }
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "DMMoab_SetElementConnectivity_Private"
153 static void DMMoab_SetElementConnectivity_Private(PetscBool useSimplex, PetscInt dim, moab::EntityType etype, PetscInt *ecount, PetscInt vpere, PetscInt nele, PetscInt i, PetscInt j, PetscInt k, PetscInt vfirst, moab::EntityHandle *connectivity)
154 {
155   PetscInt m,subelem;
156   if (useSimplex) {
157     switch (dim) {
158       case 1:
159         subelem=1;
160         DMMoab_SetSimplexElementConnectivity_Private(dim, 0, etype, (*ecount)*vpere, nele, i, j, k, vfirst, connectivity);
161         break;
162       case 2:
163         subelem=2;
164         for (m=0; m<subelem; m++)
165           DMMoab_SetSimplexElementConnectivity_Private(dim, m, etype, (*ecount+m)*vpere, nele, i, j, k, vfirst, connectivity);
166         break;
167       case 3:
168         subelem=5;
169         for (m=0; m<subelem; m++)
170           DMMoab_SetSimplexElementConnectivity_Private(dim, m, etype, (*ecount+m)*vpere, nele, i, j, k, vfirst, connectivity);
171         break;
172     }
173   }
174   else {
175     subelem=1;
176     DMMoab_SetTensorElementConnectivity_Private(dim, etype, (*ecount)*vpere, nele, i, j, k, vfirst, connectivity);
177   }
178   *ecount+=subelem;
179 }
180 
181 
182 #undef __FUNCT__
183 #define __FUNCT__ "DMMoabCreateBoxMesh"
184 /*@
185   DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with user specified bounds.
186 
187   Collective on MPI_Comm
188 
189   Input Parameters:
190 + comm - The communicator for the DM object
191 . dim - The spatial dimension
192 . 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
193 . nele - The number of discrete elements in each direction
194 . user_nghost - The number of ghosted layers needed in the partitioned mesh
195 
196   Output Parameter:
197 . dm  - The DM object
198 
199   Level: beginner
200 
201 .keywords: DM, create
202 .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile()
203 @*/
204 PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool useSimplex, const PetscReal* bounds, PetscInt nele, PetscInt user_nghost, DM *dm)
205 {
206   PetscErrorCode  ierr;
207   moab::ErrorCode merr;
208   PetscInt        i,j,k,n,nprocs;
209   DM_Moab        *dmmoab;
210   moab::Interface *mbiface;
211   moab::ParallelComm *pcomm;
212   moab::ReadUtilIface* readMeshIface;
213 
214   moab::Tag  id_tag=PETSC_NULL;
215   moab::Range         ownedvtx,ownedelms;
216   moab::EntityHandle  vfirst,efirst,regionset,faceset,edgeset,vtxset;
217   std::vector<double*> vcoords;
218   moab::EntityHandle  *connectivity = 0;
219   moab::EntityType etype;
220   PetscInt    ise[6];
221   PetscReal   xse[6],defbounds[6];
222   /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */
223   const PetscInt nghost=0;
224 
225   moab::Tag geom_tag;
226 
227   moab::Range adj,dim3,dim2;
228   bool build_adjacencies=false;
229 
230   const PetscInt npts=nele+1;        /* Number of points in every dimension */
231   PetscInt vpere,locnele,locnpts,ghnele,ghnpts;    /* Number of verts/element, vertices, elements owned by this process */
232 
233   PetscFunctionBegin;
234   if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n");
235 
236   ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
237   /* total number of vertices in all dimensions */
238   n=pow(npts,dim);
239 
240   /* do some error checking */
241   if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n");
242   if(nprocs > n) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of processors must be less than or equal to number of elements.\n");
243   if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n");
244 
245   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
246   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
247 
248   /* get all the necessary handles from the private DM object */
249   dmmoab = (DM_Moab*)(*dm)->data;
250   mbiface = dmmoab->mbiface;
251   pcomm = dmmoab->pcomm;
252   id_tag = dmmoab->ltog_tag;
253   nprocs = pcomm->size();
254   dmmoab->dim = dim;
255 
256   /* create a file set to associate all entities in current mesh */
257   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
258 
259   /* No errors yet; proceed with building the mesh */
260   merr = mbiface->query_interface(readMeshIface);MBERRNM(merr);
261 
262   ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr);
263 
264   /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */
265   ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr);
266 
267   /* set some variables based on dimension */
268   switch(dim) {
269    case 1:
270     vpere = 2;
271     locnele = (ise[1]-ise[0]);
272     locnpts = (ise[1]-ise[0]+1);
273     ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 );
274     ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0);
275     etype = moab::MBEDGE;
276     break;
277    case 2:
278     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1);
279     ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0);
280     if (useSimplex) {
281       vpere = 3;
282       locnele = 2*(ise[1]-ise[0])*(ise[3]-ise[2]);
283       ghnele = 2*(nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0);
284       etype = moab::MBTRI;
285     }
286     else {
287       vpere = 4;
288       locnele = (ise[1]-ise[0])*(ise[3]-ise[2]);
289       ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0);
290       etype = moab::MBQUAD;
291     }
292     break;
293    case 3:
294     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1);
295     ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0);
296     if (useSimplex) {
297       vpere = 4;
298       locnele = 5*(ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]);
299       ghnele = 5*(nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0);
300       etype = moab::MBTET;
301     }
302     else {
303       vpere = 8;
304       locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]);
305       ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0);
306       etype = moab::MBHEX;
307     }
308     break;
309   }
310 
311   /* we have a domain of size [1,1,1] - now compute local co-ordinate box */
312   ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr);
313   for(i=0; i<6; ++i) {
314     xse[i]=(PetscReal)ise[i]/nele;
315   }
316 
317   /* Create vertexes and set the coodinate of each vertex */
318   merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr);
319 
320   /* Compute the co-ordinates of vertices and global IDs */
321   std::vector<int>    vgid(locnpts+ghnpts);
322   int vcount=0;
323 
324   if (!bounds) { /* default box mesh is defined on a unit-cube */
325     defbounds[0]=0.0; defbounds[1]=1.0;
326     defbounds[2]=0.0; defbounds[3]=1.0;
327     defbounds[4]=0.0; defbounds[5]=1.0;
328     bounds=defbounds;
329   }
330   else {
331     /* validate the bounds data */
332     if(bounds[0] >= bounds[1]) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"X-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[0],bounds[1]);
333     if(dim > 1 && (bounds[2] >= bounds[3])) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Y-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[2],bounds[3]);
334     if(dim > 2 && (bounds[4] >= bounds[5])) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Z-dim: Left boundary cannot be greater than right. [%G >= %G]\n",bounds[4],bounds[5]);
335   }
336 
337   const double hx=(bounds[1]-bounds[0])/nele;
338   const double hy=(dim > 1 ? (bounds[3]-bounds[2])/nele : 0.0);
339   const double hz=(dim > 2 ? (bounds[5]-bounds[4])/nele : 0.0);
340 
341   /* create all the owned vertices */
342   for (k = ise[4]; k <= ise[5]; k++) {
343     for (j = ise[2]; j <= ise[3]; j++) {
344       for (i = ise[0]; i <= ise[1]; i++, vcount++) {
345         DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
346         vgid[vcount] = (k*npts+j)*npts+i+1;
347       }
348     }
349   }
350 
351   /* create ghosted vertices requested by user - below the current plane */
352   if (ise[2*dim-2] > 0) {
353     for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) {
354       for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) {
355         for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) {
356           DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
357           vgid[vcount] = (k*npts+j)*npts+i+1;
358         }
359       }
360     }
361   }
362 
363   /* create ghosted vertices requested by user - above the current plane */
364   if (ise[2*dim-1] < nele) {
365     for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) {
366       for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) {
367         for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) {
368           DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
369           vgid[vcount] = (k*npts+j)*npts+i+1;
370         }
371       }
372     }
373   }
374 
375   if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount);
376 
377   merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
378 
379   /* The global ID tag is applied to each owned
380      vertex. It acts as an global identifier which MOAB uses to
381      assemble the individual pieces of the mesh:
382      Set the global ID indices */
383   merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr);
384 
385   /* Create elements between mesh points using the ReadUtilInterface
386      get the reference to element connectivities for all local elements from the ReadUtilInterface */
387   merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr);
388 
389   /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */
390   vfirst-=vgid[0]-1;
391 
392    /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k):
393          and then set the connectivity for each element appropriately */
394   int ecount=0;
395 
396   /* create ghosted elements requested by user - below the current plane */
397   if (ise[2*dim-2] >= nghost) {
398     for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) {
399       for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) {
400         for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++) {
401           DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity);
402         }
403       }
404     }
405   }
406 
407   /* create owned elements requested by user */
408   for (k = ise[4]; k < std::max(ise[5],1); k++) {
409     for (j = ise[2]; j < std::max(ise[3],1); j++) {
410       for (i = ise[0]; i < std::max(ise[1],1); i++) {
411         DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity);
412       }
413     }
414   }
415 
416   /* create ghosted elements requested by user - above the current plane */
417   if (ise[2*dim-1] <= nele-nghost) {
418     for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) {
419       for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) {
420         for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++) {
421           DMMoab_SetElementConnectivity_Private(useSimplex, dim, etype, &ecount, vpere, nele, i, j, k, vfirst, connectivity);
422         }
423       }
424     }
425   }
426 
427   merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr);
428 
429   /* 2. Get the vertices and hexes from moab and check their numbers against I*J*K and (I-1)*(J-1)*(K-1), resp.
430         first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */
431   merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr);
432 
433   if (locnele+ghnele != (int) ownedelms.size())
434     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size());
435   else if(locnpts+ghnpts != (int) ownedvtx.size())
436     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size());
437   else
438     PetscInfo2(NULL, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size());
439 
440   /* lets create some sets */
441   merr = mbiface->tag_get_handle(GEOM_DIMENSION_TAG_NAME, 1, moab::MB_TYPE_INTEGER, geom_tag, moab::MB_TAG_SPARSE|moab::MB_TAG_CREAT);MBERRNM(merr);
442 
443   merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr);
444   merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr);
445   merr = mbiface->tag_set_data(geom_tag, &regionset, 1, &dmmoab->dim);MBERRNM(merr);
446   merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr);
447   merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr);
448 
449   merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr);
450   merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr);
451   merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr);
452   merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr);
453 
454   if (build_adjacencies) {
455     // generate all lower dimensional adjacencies
456     merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr);
457     merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr);
458     adj.merge(dim2);
459 
460     /* create face sets */
461     merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr);
462     merr = mbiface->add_entities(faceset, adj);MBERRNM(merr);
463     merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr);
464     i=dim-1;
465     merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr);
466     merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr);
467     PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-1);
468 
469     if (dim > 2) {
470       dim2.clear();
471       /* create edge sets, if appropriate i.e., if dim=3 */
472       merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr);
473       merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr);
474       merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr);
475       merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr);
476       i=dim-2;
477       merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr);
478       merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr);
479       PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-2);
480     }
481   }
482 
483   /* check the handles */
484   merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr);
485 
486   /* resolve the shared entities by exchanging information to adjacent processors */
487   merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr);
488   merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,NULL,&id_tag);MBERRV(mbiface,merr);
489 
490   merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr);
491 
492   /* Reassign global IDs on all entities. */
493   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true,false);MBERRNM(merr);
494 
495   /* Everything is set up, now just do a tag exchange to update tags
496      on all of the ghost vertexes */
497   merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
498   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr);
499 
500   merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr);
501   merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr);
502   PetscFunctionReturn(0);
503 }
504 
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "DMMoab_GetReadOptions_Private"
508 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)
509 {
510   char           ropts[PETSC_MAX_PATH_LEN];
511   char           ropts_par[PETSC_MAX_PATH_LEN];
512   char           ropts_dbg[PETSC_MAX_PATH_LEN];
513   PetscErrorCode ierr;
514 
515   PetscFunctionBegin;
516   ierr = PetscMemzero(&ropts,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
517   ierr = PetscMemzero(&ropts_par,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
518   ierr = PetscMemzero(&ropts_dbg,PETSC_MAX_PATH_LEN);CHKERRQ(ierr);
519 
520   /* do parallel read unless using only one processor */
521   if (numproc > 1) {
522     ierr = PetscSNPrintf(ropts_par, sizeof(ropts_par), "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);
523   }
524 
525   if (dbglevel) {
526     if (numproc>1) {
527       ierr = PetscSNPrintf(ropts_dbg, sizeof(ropts_dbg), "%sCPUTIME;DEBUG_IO=%d;DEBUG_PIO=%d;",dbglevel,dbglevel);CHKERRQ(ierr);
528     }
529     else {
530       ierr = PetscSNPrintf(ropts_dbg, sizeof(ropts_dbg), "%sCPUTIME;DEBUG_IO=%d;",dbglevel);CHKERRQ(ierr);
531     }
532   }
533 
534   ierr = PetscSNPrintf(ropts, sizeof(ropts), "%s%s%s%s",ropts_par,ropts_dbg,(extra_opts?extra_opts:""),(dm_opts?dm_opts:""));CHKERRQ(ierr);
535   *read_opts = ropts;
536   PetscFunctionReturn(0);
537 }
538 
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "DMMoabLoadFromFile"
542 /*@
543   DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file.
544 
545   Collective on MPI_Comm
546 
547   Input Parameters:
548 + comm - The communicator for the DM object
549 . dim - The spatial dimension
550 . filename - The name of the mesh file to be loaded
551 . usrreadopts - The options string to read a MOAB mesh.
552   Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo)
553 
554   Output Parameter:
555 . dm  - The DM object
556 
557   Level: beginner
558 
559 .keywords: DM, create
560 
561 .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh()
562 @*/
563 PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm)
564 {
565   moab::ErrorCode merr;
566   PetscInt        nprocs;
567   DM_Moab        *dmmoab;
568   moab::Interface *mbiface;
569   moab::ParallelComm *pcomm;
570   moab::Range verts,elems;
571   const char *readopts;
572   PetscErrorCode ierr;
573 
574   PetscFunctionBegin;
575   PetscValidPointer(dm,5);
576 
577   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
578   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
579 
580   /* get all the necessary handles from the private DM object */
581   dmmoab = (DM_Moab*)(*dm)->data;
582   mbiface = dmmoab->mbiface;
583   pcomm = dmmoab->pcomm;
584   nprocs = pcomm->size();
585   /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */
586   dmmoab->dim = dim;
587 
588   /* create a file set to associate all entities in current mesh */
589   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
590 
591   /* add mesh loading options specific to the DM */
592   ierr = DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, dmmoab->read_mode,
593                                         dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts);CHKERRQ(ierr);
594 
595   PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts);
596 
597   /* Load the mesh from a file. */
598   merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr);
599 
600   /* Reassign global IDs on all entities. */
601   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr);
602 
603   /* load the local vertices */
604   merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr);
605   /* load the local elements */
606   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr);
607 
608   /* Everything is set up, now just do a tag exchange to update tags
609      on all of the ghost vertexes */
610   merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr);
611   merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr);
612 
613   merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr);
614 
615   merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr);
616 
617   PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size());
618   PetscFunctionReturn(0);
619 }
620 
621