xref: /petsc/src/dm/impls/moab/dmmbutil.cxx (revision 2e4e7c01f2fc776e49aa0e3bac11897d06784bfb)
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_SetElementConnectivity_Private"
57 static void DMMoab_SetElementConnectivity_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__ "DMMoabCreateBoxMesh"
90 /*@
91   DMMoabCreateBoxMesh - Creates a mesh on the tensor product (box) of intervals with user specified bounds.
92 
93   Collective on MPI_Comm
94 
95   Input Parameters:
96 + comm - The communicator for the DM object
97 . dim - The spatial dimension
98 . 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
99 . nele - The number of discrete elements in each direction
100 . user_nghost - The number of ghosted layers needed in the partitioned mesh
101 
102   Output Parameter:
103 . dm  - The DM object
104 
105   Level: beginner
106 
107 .keywords: DM, create
108 .seealso: DMSetType(), DMCreate(), DMMoabLoadFromFile()
109 @*/
110 PetscErrorCode DMMoabCreateBoxMesh(MPI_Comm comm, PetscInt dim, const PetscReal* bounds, PetscInt nele, PetscInt user_nghost, DM *dm)
111 {
112   PetscErrorCode  ierr;
113   moab::ErrorCode merr;
114   PetscInt        i,j,k,n,nprocs;
115   DM_Moab        *dmmoab;
116   moab::Interface *mbiface;
117   moab::ParallelComm *pcomm;
118   moab::ReadUtilIface* readMeshIface;
119 
120   moab::Tag  id_tag=PETSC_NULL;
121   moab::Range         ownedvtx,ownedelms;
122   moab::EntityHandle  vfirst,efirst,regionset,faceset,edgeset,vtxset;
123   std::vector<double*> vcoords;
124   moab::EntityHandle  *connectivity = 0;
125   moab::EntityType etype;
126   PetscInt    ise[6];
127   PetscReal   xse[6],defbounds[6];
128   /* TODO: Fix nghost > 0 - now relying on exchange_ghost_cells */
129   const PetscInt nghost=0;
130 
131   moab::Tag geom_tag;
132 
133   moab::Range adj,dim3,dim2;
134   bool build_adjacencies=false;
135 
136   const PetscInt npts=nele+1;        /* Number of points in every dimension */
137   PetscInt vpere,locnele,locnpts,ghnele,ghnpts;    /* Number of verts/element, vertices, elements owned by this process */
138 
139   PetscFunctionBegin;
140   if(dim < 1 || dim > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid dimension argument for mesh: dim=[1,3].\n");
141 
142   ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
143   /* total number of vertices in all dimensions */
144   n=pow(npts,dim);
145 
146   /* do some error checking */
147   if(n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be >= 2.\n");
148   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");
149   if(nghost < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Number of ghost layers cannot be negative.\n");
150 
151   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
152   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
153 
154   /* get all the necessary handles from the private DM object */
155   dmmoab = (DM_Moab*)(*dm)->data;
156   mbiface = dmmoab->mbiface;
157   pcomm = dmmoab->pcomm;
158   id_tag = dmmoab->ltog_tag;
159   nprocs = pcomm->size();
160   dmmoab->dim = dim;
161 
162   /* create a file set to associate all entities in current mesh */
163   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
164 
165   /* No errors yet; proceed with building the mesh */
166   merr = mbiface->query_interface(readMeshIface);MBERRNM(merr);
167 
168   ierr = PetscMemzero(ise,sizeof(PetscInt)*6);CHKERRQ(ierr);
169 
170   /* call the collective routine that computes the domain bounds for a structured mesh using MOAB */
171   ierr = DMMoabComputeDomainBounds_Private(pcomm, dim, nele, ise);CHKERRQ(ierr);
172 
173   /* set some variables based on dimension */
174   switch(dim) {
175    case 1:
176     vpere = 2;
177     locnele = (ise[1]-ise[0]);
178     locnpts = (ise[1]-ise[0]+1);
179     ghnele = (nghost > 0 ? (ise[0] > nghost ? 1 : 0) + (ise[1] < nele - nghost ? 1 : 0) : 0 );
180     ghnpts = (nghost > 0 ? (ise[0] > 0 ? 1 : 0) + (ise[1] < nele ? 1 : 0) : 0);
181     etype = moab::MBEDGE;
182     break;
183    case 2:
184     vpere = 4;
185     locnele = (ise[1]-ise[0])*(ise[3]-ise[2]);
186     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1);
187     ghnele = (nghost > 0 ? (ise[2] > 0 ? nele : 0) + (ise[3] < nele ? nele : 0) : 0);
188     ghnpts = (nghost > 0 ? (ise[2] > 0 ? npts : 0) + (ise[3] < nele ? npts : 0) : 0);
189     etype = moab::MBQUAD;
190     break;
191    case 3:
192     vpere = 8;
193     locnele = (ise[1]-ise[0])*(ise[3]-ise[2])*(ise[5]-ise[4]);
194     locnpts = (ise[1]-ise[0]+1)*(ise[3]-ise[2]+1)*(ise[5]-ise[4]+1);
195     ghnele = (nghost > 0 ? (ise[4] > 0 ? nele*nele : 0) + (ise[5] < nele ? nele*nele : 0) : 0);
196     ghnpts = (nghost > 0 ? (ise[4] > 0 ? npts*npts : 0) + (ise[5] < nele ? npts*npts : 0) : 0);
197     etype = moab::MBHEX;
198     break;
199   }
200 
201   /* we have a domain of size [1,1,1] - now compute local co-ordinate box */
202   ierr = PetscMemzero(xse,sizeof(PetscReal)*6);CHKERRQ(ierr);
203   for(i=0; i<6; ++i) {
204     xse[i]=(PetscReal)ise[i]/nele;
205   }
206 
207   /* Create vertexes and set the coodinate of each vertex */
208   merr = readMeshIface->get_node_coords(3,locnpts+ghnpts,0,vfirst,vcoords,n);MBERRNM(merr);
209 
210   /* Compute the co-ordinates of vertices and global IDs */
211   std::vector<int>    vgid(locnpts+ghnpts);
212   int vcount=0;
213 
214   if (!bounds) { /* default box mesh is defined on a unit-cube */
215     defbounds[0]=0.0; defbounds[1]=1.0;
216     defbounds[2]=0.0; defbounds[3]=1.0;
217     defbounds[4]=0.0; defbounds[5]=1.0;
218     bounds=defbounds;
219   }
220   else {
221     /* validate the bounds data */
222     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]);
223     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]);
224     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]);
225   }
226 
227   const double hx=(bounds[1]-bounds[0])/nele;
228   const double hy=(dim > 1 ? (bounds[3]-bounds[2])/nele : 0.0);
229   const double hz=(dim > 2 ? (bounds[5]-bounds[4])/nele : 0.0);
230 
231   /* create all the owned vertices */
232   for (k = ise[4]; k <= ise[5]; k++) {
233     for (j = ise[2]; j <= ise[3]; j++) {
234       for (i = ise[0]; i <= ise[1]; i++, vcount++) {
235         DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
236         vgid[vcount] = (k*npts+j)*npts+i+1;
237       }
238     }
239   }
240 
241   /* create ghosted vertices requested by user - below the current plane */
242   if (ise[2*dim-2] > 0) {
243     for (k = (dim==3?ise[4]-nghost:ise[4]); k <= (dim==3?ise[4]-1:ise[5]); k++) {
244       for (j = (dim==2?ise[2]-nghost:ise[2]); j <= (dim==2?ise[2]-1:ise[3]); j++) {
245         for (i = (dim>1?ise[0]:ise[0]-nghost); i <= (dim>1?ise[1]:ise[0]-1); i++, vcount++) {
246           DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
247           vgid[vcount] = (k*npts+j)*npts+i+1;
248         }
249       }
250     }
251   }
252 
253   /* create ghosted vertices requested by user - above the current plane */
254   if (ise[2*dim-1] < nele) {
255     for (k = (dim==3?ise[5]+1:ise[4]); k <= (dim==3?ise[5]+nghost:ise[5]); k++) {
256       for (j = (dim==2?ise[3]+1:ise[2]); j <= (dim==2?ise[3]+nghost:ise[3]); j++) {
257         for (i = (dim>1?ise[0]:ise[1]+1); i <= (dim>1?ise[1]:ise[1]+nghost); i++, vcount++) {
258           DMMoab_SetStructuredCoords_Private(i,j,k,hx,hy,hz,vcount,vcoords);
259           vgid[vcount] = (k*npts+j)*npts+i+1;
260         }
261       }
262     }
263   }
264 
265   if (locnpts+ghnpts != vcount) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,vcount);
266 
267   merr = mbiface->get_entities_by_type(0,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
268 
269   /* The global ID tag is applied to each owned
270      vertex. It acts as an global identifier which MOAB uses to
271      assemble the individual pieces of the mesh:
272      Set the global ID indices */
273   merr = mbiface->tag_set_data(id_tag,ownedvtx,vgid.data());MBERRNM(merr);
274 
275   /* Create elements between mesh points using the ReadUtilInterface
276      get the reference to element connectivities for all local elements from the ReadUtilInterface */
277   merr = readMeshIface->get_element_connect (locnele+ghnele,vpere,etype,1,efirst,connectivity);MBERRNM(merr);
278 
279   /* offset appropriately so that only local ID and not global ID numbers are set for connectivity array */
280 //  PetscPrintf(PETSC_COMM_SELF, "[%D] first local handle %D\n", rank, vfirst);
281   vfirst-=vgid[0]-1;
282 
283    /* 3. Loop over elements in 3 nested loops over i, j, k; for each (i,j,k):
284          and then set the connectivity for each element appropriately */
285   int ecount=0;
286 
287   /* create ghosted elements requested by user - below the current plane */
288   if (ise[2*dim-2] >= nghost) {
289     for (k = (dim==3?ise[4]-nghost:ise[4]); k < (dim==3?ise[4]:std::max(ise[5],1)); k++) {
290       for (j = (dim==2?ise[2]-nghost:ise[2]); j < (dim==2?ise[2]:std::max(ise[3],1)); j++) {
291         for (i = (dim>1?ise[0]:ise[0]-nghost); i < (dim>1?std::max(ise[1],1):ise[0]); i++, ecount++) {
292           DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
293         }
294       }
295     }
296   }
297 
298   /* create owned elements requested by user */
299   for (k = ise[4]; k < std::max(ise[5],1); k++) {
300     for (j = ise[2]; j < std::max(ise[3],1); j++) {
301       for (i = ise[0]; i < std::max(ise[1],1); i++,ecount++) {
302         DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
303       }
304     }
305   }
306 
307   /* create ghosted elements requested by user - above the current plane */
308   if (ise[2*dim-1] <= nele-nghost) {
309     for (k = (dim==3?ise[5]:ise[4]); k < (dim==3?ise[5]+nghost:std::max(ise[5],1)); k++) {
310       for (j = (dim==2?ise[3]:ise[2]); j < (dim==2?ise[3]+nghost:std::max(ise[3],1)); j++) {
311         for (i = (dim>1?ise[0]:ise[1]); i < (dim>1?std::max(ise[1],1):ise[1]+nghost); i++, ecount++) {
312           DMMoab_SetElementConnectivity_Private(dim, etype, ecount*vpere, nele, i, j, k, vfirst, connectivity);
313         }
314       }
315     }
316   }
317 
318   merr = readMeshIface->update_adjacencies(efirst,locnele+ghnele,vpere,connectivity);MBERRNM(merr);
319 
320   /* 2. Get the vertices and hexes from moab and check their numbers against I*J*K and (I-1)*(J-1)*(K-1), resp.
321         first '0' specifies "root set", or entire MOAB instance, second the entity dimension being requested */
322   merr = mbiface->get_entities_by_dimension(0, dim, ownedelms);MBERRNM(merr);
323 
324   if (locnele+ghnele != (int) ownedelms.size())
325     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of elements! (%D!=%D)",locnele+ghnele,ownedelms.size());
326   else if(locnpts+ghnpts != (int) ownedvtx.size())
327     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Created the wrong number of vertices! (%D!=%D)",locnpts+ghnpts,ownedvtx.size());
328   else
329     PetscInfo2(NULL, "Created %D elements and %D vertices.\n", ownedelms.size(), ownedvtx.size());
330 
331   /* lets create some sets */
332   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);
333 
334   merr = mbiface->create_meshset(moab::MESHSET_SET, regionset);MBERRNM(merr);
335   merr = mbiface->add_entities(regionset, ownedelms);MBERRNM(merr);
336   merr = mbiface->tag_set_data(geom_tag, &regionset, 1, &dmmoab->dim);MBERRNM(merr);
337   merr = mbiface->add_parent_child(dmmoab->fileset,regionset);MBERRNM(merr);
338   merr = mbiface->unite_meshset(dmmoab->fileset, regionset);MBERRNM(merr);
339 
340   merr = mbiface->create_meshset(moab::MESHSET_SET, vtxset);MBERRNM(merr);
341   merr = mbiface->add_entities(vtxset, ownedvtx);MBERRNM(merr);
342   merr = mbiface->add_parent_child(dmmoab->fileset,vtxset);MBERRNM(merr);
343   merr = mbiface->unite_meshset(dmmoab->fileset, vtxset);MBERRNM(merr);
344 
345   if (build_adjacencies) {
346     // generate all lower dimensional adjacencies
347     merr = mbiface->get_adjacencies( ownedelms, dim-1, true, adj, moab::Interface::UNION );MBERRNM(merr);
348     merr = dmmoab->pcomm->get_part_entities(dim2, dim-1);MBERRNM(merr);
349     adj.merge(dim2);
350 
351     /* create face sets */
352     merr = mbiface->create_meshset(moab::MESHSET_SET, faceset);MBERRNM(merr);
353     merr = mbiface->add_entities(faceset, adj);MBERRNM(merr);
354     merr = mbiface->add_parent_child(dmmoab->fileset,faceset);MBERRNM(merr);
355     i=dim-1;
356     merr = mbiface->tag_set_data(geom_tag, &faceset, 1, &i);MBERRNM(merr);
357     merr = mbiface->unite_meshset(dmmoab->fileset, faceset);MBERRNM(merr);
358     PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-1);
359 
360     if (dim > 2) {
361       dim2.clear();
362       /* create edge sets, if appropriate i.e., if dim=3 */
363       merr = mbiface->create_meshset(moab::MESHSET_SET, edgeset);MBERRNM(merr);
364       merr = mbiface->get_adjacencies(adj, dim-1, true, dim2, moab::Interface::UNION );MBERRNM(merr);
365       merr = mbiface->add_entities(edgeset, dim2);MBERRNM(merr);
366       merr = mbiface->add_parent_child(dmmoab->fileset,edgeset);MBERRNM(merr);
367       i=dim-2;
368       merr = mbiface->tag_set_data(geom_tag, &edgeset, 1, &i);MBERRNM(merr);
369       merr = mbiface->unite_meshset(dmmoab->fileset, edgeset);MBERRNM(merr);
370       PetscInfo2(NULL, "Found %d %d-Dim quantities.\n", adj.size(), dim-2);
371     }
372   }
373 
374   /* check the handles */
375   merr = pcomm->check_all_shared_handles();MBERRV(mbiface,merr);
376 
377   /* resolve the shared entities by exchanging information to adjacent processors */
378   merr = mbiface->get_entities_by_type(dmmoab->fileset,etype,ownedelms,true);MBERRNM(merr);
379   merr = pcomm->resolve_shared_ents(dmmoab->fileset,ownedelms,dim,dim-1,NULL,&id_tag);MBERRV(mbiface,merr);
380 
381   merr = pcomm->exchange_ghost_cells(dim,0,user_nghost,dim,true,false,&dmmoab->fileset);MBERRV(mbiface,merr);
382 
383   /* Reassign global IDs on all entities. */
384   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,false,true,false);MBERRNM(merr);
385 
386   /* Everything is set up, now just do a tag exchange to update tags
387      on all of the ghost vertexes */
388   merr = mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,ownedvtx,true);MBERRNM(merr);
389   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, ownedelms);MBERRNM(merr);
390 
391   merr = pcomm->exchange_tags(id_tag,ownedvtx);MBERRV(mbiface,merr);
392   merr = pcomm->exchange_tags(id_tag,ownedelms);MBERRV(mbiface,merr);
393   PetscFunctionReturn(0);
394 }
395 
396 
397 #undef __FUNCT__
398 #define __FUNCT__ "DMMoab_GetReadOptions_Private"
399 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)
400 {
401   std::ostringstream str;
402 
403   PetscFunctionBegin;
404   /* do parallel read unless using only one processor */
405   if (numproc > 1) {
406     str << "PARALLEL=" << MoabReadModes[mode] << ";PARTITION=PARALLEL_PARTITION;PARTITION_DISTRIBUTE;";
407     str << "PARALLEL_RESOLVE_SHARED_ENTS;PARALLEL_GHOSTS=" << dim << ".0.1;";
408     if (by_rank)
409       str << "PARTITION_BY_RANK;";
410   }
411 
412   if (strlen(dm_opts)) {
413     str << dm_opts << ";";
414   }
415 
416   if (dbglevel) {
417     str << "CPUTIME;DEBUG_IO=" << dbglevel << ";";
418     if (numproc>1)
419       str << "DEBUG_PIO=" << dbglevel << ";";
420   }
421 
422   if (extra_opts)
423     str << extra_opts;
424 
425   *read_opts = str.str().c_str();
426   PetscFunctionReturn(0);
427 }
428 
429 
430 #undef __FUNCT__
431 #define __FUNCT__ "DMMoabLoadFromFile"
432 /*@
433   DMMoabLoadFromFile - Creates a DM object by loading the mesh from a user specified file.
434 
435   Collective on MPI_Comm
436 
437   Input Parameters:
438 + comm - The communicator for the DM object
439 . dim - The spatial dimension
440 . filename - The name of the mesh file to be loaded
441 . usrreadopts - The options string to read a MOAB mesh.
442   Reference (Parallel Mesh Initialization: http://www.mcs.anl.gov/~fathom/moab-docs/html/contents.html#fivetwo)
443 
444   Output Parameter:
445 . dm  - The DM object
446 
447   Level: beginner
448 
449 .keywords: DM, create
450 
451 .seealso: DMSetType(), DMCreate(), DMMoabCreateBoxMesh()
452 @*/
453 PetscErrorCode DMMoabLoadFromFile(MPI_Comm comm,PetscInt dim,const char* filename, const char* usrreadopts, DM *dm)
454 {
455   moab::ErrorCode merr;
456   PetscInt        nprocs;
457   DM_Moab        *dmmoab;
458   moab::Interface *mbiface;
459   moab::ParallelComm *pcomm;
460   moab::Range verts,elems;
461   const char *readopts;
462   PetscErrorCode ierr;
463 
464   PetscFunctionBegin;
465   PetscValidPointer(dm,5);
466 
467   /* Create the basic DMMoab object and keep the default parameters created by DM impls */
468   ierr = DMMoabCreateMoab(comm, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, dm);CHKERRQ(ierr);
469 
470   /* get all the necessary handles from the private DM object */
471   dmmoab = (DM_Moab*)(*dm)->data;
472   mbiface = dmmoab->mbiface;
473   pcomm = dmmoab->pcomm;
474   nprocs = pcomm->size();
475   /* TODO: Decipher dimension based on the loaded mesh instead of getting from user */
476   dmmoab->dim = dim;
477 
478   /* create a file set to associate all entities in current mesh */
479   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, dmmoab->fileset);MBERR("Creating file set failed", merr);
480 
481   /* add mesh loading options specific to the DM */
482   ierr = DMMoab_GetReadOptions_Private(dmmoab->partition_by_rank, nprocs, dim, dmmoab->read_mode,
483                                         dmmoab->rw_dbglevel, dmmoab->extra_read_options, usrreadopts, &readopts);CHKERRQ(ierr);
484 
485   PetscInfo2(*dm, "Reading file %s with options: %s\n",filename,readopts);
486 
487   /* Load the mesh from a file. */
488   merr = mbiface->load_file(filename, &dmmoab->fileset, readopts);MBERRVM(mbiface,"Reading MOAB file failed.", merr);
489 
490   /* Reassign global IDs on all entities. */
491   merr = pcomm->assign_global_ids(dmmoab->fileset,dim,1,true,true,true);MBERRNM(merr);
492 
493   /* load the local vertices */
494   merr = mbiface->get_entities_by_type(dmmoab->fileset, moab::MBVERTEX, verts, true);MBERRNM(merr);
495   /* load the local elements */
496   merr = mbiface->get_entities_by_dimension(dmmoab->fileset, dim, elems, true);MBERRNM(merr);
497 
498   /* Everything is set up, now just do a tag exchange to update tags
499      on all of the ghost vertexes */
500   merr = pcomm->exchange_tags(dmmoab->ltog_tag,verts);MBERRV(mbiface,merr);
501   merr = pcomm->exchange_tags(dmmoab->ltog_tag,elems);MBERRV(mbiface,merr);
502 
503   merr = pcomm->exchange_ghost_cells(dim,0,1,0,true,true,&dmmoab->fileset);MBERRV(mbiface,merr);
504 
505   merr = pcomm->collective_sync_partition();MBERR("Collective sync failed", merr);
506 
507   PetscInfo3(*dm, "MOAB file '%s' was successfully loaded. Found %D vertices and %D elements.\n", filename, verts.size(), elems.size());
508   PetscFunctionReturn(0);
509 }
510 
511