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