xref: /petsc/src/dm/impls/moab/dmmoab.cxx (revision 324f1edfa1c27a5f313d88127d17fd9181f21be7)
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/Skinner.hpp>
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "DMDestroy_Moab"
10 PetscErrorCode DMDestroy_Moab(DM dm)
11 {
12   PetscErrorCode ierr;
13   DM_Moab        *dmmoab = (DM_Moab*)dm->data;
14   PetscSection   section;
15 
16   PetscFunctionBegin;
17   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
18   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
19   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
20   if (dmmoab->icreatedinstance) {
21     delete dmmoab->mbiface;
22   }
23   dmmoab->mbiface = NULL;
24   dmmoab->pcomm = NULL;
25   delete dmmoab->vlocal;
26   delete dmmoab->vowned;
27   delete dmmoab->vghost;
28   delete dmmoab->elocal;
29   delete dmmoab->eghost;
30   delete dmmoab->bndyvtx;
31   delete dmmoab->bndyfaces;
32   delete dmmoab->bndyelems;
33 
34   ierr = PetscFree(dmmoab->gsindices);CHKERRQ(ierr);
35   ierr = VecScatterDestroy(&dmmoab->ltog_sendrecv);CHKERRQ(ierr);
36   ierr = ISLocalToGlobalMappingDestroy(&dmmoab->ltog_map);CHKERRQ(ierr);
37   ierr = PetscFree(dm->data);CHKERRQ(ierr);
38   PetscFunctionReturn(0);
39 }
40 
41 #undef __FUNCT__
42 #define __FUNCT__ "DMSetUp_Moab"
43 PetscErrorCode DMSetUp_Moab(DM dm)
44 {
45   PetscErrorCode          ierr;
46   moab::ErrorCode         merr;
47   Vec                     local, global;
48   IS                      from;
49   moab::Range::iterator   iter;
50   PetscInt                i,j,bs,gsiz,lsiz;
51   DM_Moab                *dmmoab = (DM_Moab*)dm->data;
52   PetscInt                totsize;
53   PetscSection            section;
54   PetscInt                gmin,lmin,lmax;
55 
56   moab::Range adj;
57 
58   PetscFunctionBegin;
59   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
60   /* Get the local and shared vertices and cache it */
61   if (dmmoab->mbiface == PETSC_NULL || dmmoab->pcomm == PETSC_NULL) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ORDER, "Set the MOAB Interface and ParallelComm objects before calling SetUp.");
62 
63   /* Get the entities recursively in the current part of the mesh, if user did not set the local vertices explicitly */
64   if (dmmoab->vlocal->empty()) {
65     merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,*dmmoab->vlocal,true);MBERRNM(merr);
66 
67     /* filter based on parallel status */
68     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,dmmoab->vowned);MBERRNM(merr);
69     *dmmoab->vghost = moab::subtract(*dmmoab->vlocal, *dmmoab->vowned);
70 
71     dmmoab->nloc = dmmoab->vowned->size();
72     dmmoab->nghost = dmmoab->vghost->size();
73     ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
74 
75 #if 0
76     if(dmmoab->pcomm->rank() || dmmoab->pcomm->size()==1) {
77       PetscPrintf(PETSC_COMM_SELF, "Vertices: global: %D, local: %D", dmmoab->n, dmmoab->nloc+dmmoab->nghost);
78       dmmoab->vlocal->print(0);
79       PetscPrintf(PETSC_COMM_SELF, "Vertices: owned: %D", dmmoab->nloc);
80       dmmoab->vowned->print(0);
81       PetscPrintf(PETSC_COMM_SELF, "Vertices: ghost: %D", dmmoab->nghost);
82       dmmoab->vghost->print(0);
83     }
84 #endif
85   }
86 
87   /* get the information about the local elements in the mesh */
88   {
89     dmmoab->eghost->clear();
90 
91     /* first decipher the leading dimension */
92     for (i=3;i>0;i--) {
93       dmmoab->elocal->clear();
94       merr = dmmoab->mbiface->get_entities_by_dimension(dmmoab->fileset, i, *dmmoab->elocal, true);CHKERRQ(merr);
95 
96       /* store the current mesh dimension */
97       if (dmmoab->elocal->size()) {
98         dmmoab->dim=i;
99         break;
100       }
101     }
102 
103     *dmmoab->eghost = *dmmoab->elocal;
104     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
105     *dmmoab->eghost = moab::subtract(*dmmoab->eghost, *dmmoab->elocal);
106 
107     dmmoab->neleloc = dmmoab->elocal->size();
108     ierr = MPI_Allreduce(&dmmoab->neleloc, &dmmoab->nele, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
109   }
110 
111   bs = dmmoab->bs;
112   if (!dmmoab->ltog_tag) {
113     /* Get the global ID tag. The global ID tag is applied to each
114        vertex. It acts as an global identifier which MOAB uses to
115        assemble the individual pieces of the mesh */
116     merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr);
117   }
118 
119   totsize=dmmoab->vlocal->size();
120   ierr = PetscMalloc(totsize*sizeof(PetscInt), &dmmoab->gsindices);CHKERRQ(ierr);
121   {
122     /* first get the local indices */
123     merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vowned,&dmmoab->gsindices[0]);MBERRNM(merr);
124     /* next get the ghosted indices */
125     if (dmmoab->nghost) {
126       merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vghost,&dmmoab->gsindices[dmmoab->nloc]);MBERRNM(merr);
127     }
128 
129     /* find out the local and global minima of GLOBAL_ID */
130     lmin=lmax=dmmoab->gsindices[0];
131     for (i=0; i<totsize; ++i) {
132       if(lmin>dmmoab->gsindices[i]) lmin=dmmoab->gsindices[i];
133       if(lmax<dmmoab->gsindices[i]) lmax=dmmoab->gsindices[i];
134     }
135 
136     ierr = MPI_Allreduce(&lmin, &gmin, 1, MPI_INT, MPI_MIN, ((PetscObject)dm)->comm);CHKERRQ(ierr);
137     PetscInfo3(dm, "GLOBAL_ID: Local minima - %D, Local maxima - %D, Global minima - %D.\n", lmin, lmax, gmin);
138   }
139 
140   {
141     ierr = PetscSectionCreate(((PetscObject)dm)->comm, &section);CHKERRQ(ierr);
142     ierr = PetscSectionSetNumFields(section, dmmoab->nfields);CHKERRQ(ierr);
143     ierr = PetscSectionSetChart(section, lmin, lmax+1);CHKERRQ(ierr);
144     for (j=0; j<totsize; ++j) {
145       PetscInt locgid = dmmoab->gsindices[j];
146       for (i=0; i < dmmoab->nfields; ++i) {
147         ierr = PetscSectionSetFieldName(section, i, dmmoab->fields[i]);CHKERRQ(ierr);
148         if (bs>1) {
149           ierr = PetscSectionSetFieldDof(section, locgid, i, (locgid-gmin)*dmmoab->nfields+i);CHKERRQ(ierr);
150           ierr = PetscSectionSetFieldOffset(section, locgid, i, (locgid-gmin)*dmmoab->nfields);
151         }
152         else {
153           ierr = PetscSectionSetFieldDof(section, locgid, i, dmmoab->n*i+locgid-gmin);CHKERRQ(ierr);
154           ierr = PetscSectionSetFieldOffset(section, locgid, i, i*dmmoab->n);
155         }
156       }
157       ierr = PetscSectionSetDof(section, locgid, dmmoab->nfields);CHKERRQ(ierr);
158     }
159     ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
160     ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr);
161   }
162 
163   {
164     for (i=0; i<totsize; ++i) {
165       dmmoab->gsindices[i]-=gmin;   /* zero based index needed for IS */
166     }
167 
168     /* Create Global to Local Vector Scatter Context */
169     ierr = DMCreateGlobalVector_Moab(dm, &global);CHKERRQ(ierr);
170     ierr = DMCreateLocalVector_Moab(dm, &local);CHKERRQ(ierr);
171 
172     /* global to local must retrieve ghost points */
173     ierr = ISCreateBlock(((PetscObject)dm)->comm,bs,totsize,&dmmoab->gsindices[0],PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
174 
175     ierr = VecGetLocalSize(global,&gsiz);CHKERRQ(ierr);
176     ierr = VecGetLocalSize(local,&lsiz);CHKERRQ(ierr);
177 
178     ierr = VecScatterCreate(local,from,global,from,&dmmoab->ltog_sendrecv);CHKERRQ(ierr);
179     ierr = ISDestroy(&from);CHKERRQ(ierr);
180     ierr = VecDestroy(&local);CHKERRQ(ierr);
181     ierr = VecDestroy(&global);CHKERRQ(ierr);
182   }
183 
184   /* skin the boundary and store nodes */
185   {
186     /* get the skin vertices of boundary faces for the current partition and then filter
187        the local, boundary faces, vertices and elements alone via PSTATUS flags;
188        this should not give us any ghosted boundary, but if user needs such a functionality
189        it would be easy to add it based on the find_skin query below */
190     moab::Skinner skinner(dmmoab->mbiface);
191     dmmoab->bndyvtx = new moab::Range();
192     dmmoab->bndyfaces = new moab::Range();
193     dmmoab->bndyelems = new moab::Range();
194 
195     /* get the entities on the skin - only the faces */
196     merr = skinner.find_skin(dmmoab->fileset, *dmmoab->elocal, false, *dmmoab->bndyfaces, NULL, false, true, false, false);MBERRNM(merr); // 'false' param indicates we want faces back, not vertices
197 
198     /* filter all the non-owned and shared entities out of the list */
199     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->bndyfaces,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
200     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->bndyfaces,PSTATUS_SHARED,PSTATUS_NOT);MBERRNM(merr);
201 
202     if (dmmoab->dim == 3) {
203       // get the edges from faces and then do the same if needed
204     }
205 
206     /* get all the nodes via connectivity and the parent elements via adjacency information */
207     merr = dmmoab->mbiface->get_connectivity(*dmmoab->bndyfaces, *dmmoab->bndyvtx, false);MBERRNM(ierr);
208     merr = dmmoab->mbiface->get_adjacencies(*dmmoab->bndyfaces, dmmoab->dim, false, *dmmoab->bndyelems, moab::Interface::UNION);MBERRNM(ierr);
209     PetscInfo3(dm, "Found %D boundary vertices, %D boundary faces and %D boundary elements.\n", dmmoab->bndyvtx->size(), dmmoab->bndyvtx->size(), dmmoab->bndyelems->size());
210   }
211   PetscFunctionReturn(0);
212 }
213 
214 
215 #undef __FUNCT__
216 #define __FUNCT__ "DMCreate_Moab"
217 PETSC_EXTERN PetscErrorCode DMCreate_Moab(DM dm)
218 {
219   PetscErrorCode ierr;
220 
221   PetscFunctionBegin;
222   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
223   ierr = PetscNewLog(dm,&dm->data);CHKERRQ(ierr);
224 
225   ((DM_Moab*)dm->data)->bs = 1;
226   ((DM_Moab*)dm->data)->nfields = 1;
227   ((DM_Moab*)dm->data)->n = 0;
228   ((DM_Moab*)dm->data)->nloc = 0;
229   ((DM_Moab*)dm->data)->nele = 0;
230   ((DM_Moab*)dm->data)->neleloc = 0;
231   ((DM_Moab*)dm->data)->nghost = 0;
232   ((DM_Moab*)dm->data)->ltog_map = PETSC_NULL;
233   ((DM_Moab*)dm->data)->ltog_sendrecv = PETSC_NULL;
234 
235   ((DM_Moab*)dm->data)->vlocal = new moab::Range();
236   ((DM_Moab*)dm->data)->vowned = new moab::Range();
237   ((DM_Moab*)dm->data)->vghost = new moab::Range();
238   ((DM_Moab*)dm->data)->elocal = new moab::Range();
239   ((DM_Moab*)dm->data)->eghost = new moab::Range();
240 
241   dm->ops->createglobalvector              = DMCreateGlobalVector_Moab;
242   dm->ops->createlocalvector               = DMCreateLocalVector_Moab;
243   dm->ops->creatematrix                    = DMCreateMatrix_Moab;
244   dm->ops->setup                           = DMSetUp_Moab;
245   dm->ops->destroy                         = DMDestroy_Moab;
246   dm->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Moab;
247   dm->ops->globaltolocalend                = DMGlobalToLocalEnd_Moab;
248   dm->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Moab;
249   dm->ops->localtoglobalend                = DMLocalToGlobalEnd_Moab;
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "DMMoabCreate"
255 /*@
256   DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance
257 
258   Collective on MPI_Comm
259 
260   Input Parameter:
261 . comm - The communicator for the DMMoab object
262 
263   Output Parameter:
264 . dmb  - The DMMoab object
265 
266   Level: beginner
267 
268 .keywords: DMMoab, create
269 @*/
270 PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *dmb)
271 {
272   PetscErrorCode ierr;
273 
274   PetscFunctionBegin;
275   PetscValidPointer(dmb,2);
276   ierr = DMCreate(comm, dmb);CHKERRQ(ierr);
277   ierr = DMSetType(*dmb, DMMOAB);CHKERRQ(ierr);
278   PetscFunctionReturn(0);
279 }
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "DMMoabCreateMoab"
283 /*@
284   DMMoabCreate - Creates a DMMoab object, optionally from an instance and other data
285 
286   Collective on MPI_Comm
287 
288   Input Parameter:
289 . comm - The communicator for the DMMoab object
290 . mbiface - (ptr to) the MOAB Instance; if passed in NULL, MOAB instance is created inside PETSc, and destroyed
291          along with the DMMoab
292 . pcomm - (ptr to) a ParallelComm; if NULL, creates one internally for the whole communicator
293 . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag
294 . range - If non-NULL, contains range of entities to which DOFs will be assigned
295 
296   Output Parameter:
297 . dmb  - The DMMoab object
298 
299   Level: intermediate
300 
301 .keywords: DMMoab, create
302 @*/
303 PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *dmb)
304 {
305   PetscErrorCode ierr;
306   moab::ErrorCode merr;
307   moab::EntityHandle partnset;
308   PetscInt rank, nprocs;
309   DM_Moab        *dmmoab;
310 
311   PetscFunctionBegin;
312   PetscValidPointer(dmb,6);
313   ierr = DMMoabCreate(comm, dmb);CHKERRQ(ierr);
314   dmmoab = (DM_Moab*)(*dmb)->data;
315 
316   if (!mbiface) {
317     dmmoab->mbiface = new moab::Core();
318     dmmoab->icreatedinstance = PETSC_TRUE;
319   }
320   else {
321     dmmoab->mbiface = mbiface;
322     dmmoab->icreatedinstance = PETSC_FALSE;
323   }
324 
325   /* create a fileset to store the hierarchy of entities belonging to current DM */
326   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_ORDERED, dmmoab->fileset);MBERR("Creating file set failed", merr);
327 
328   if (!pcomm) {
329     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
330     ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
331 
332     /* Create root sets for each mesh.  Then pass these
333        to the load_file functions to be populated. */
334     merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, partnset);
335     MBERR("Creating partition set failed", merr);
336 
337     /* Create the parallel communicator object with the partition handle associated with MOAB */
338     dmmoab->pcomm = moab::ParallelComm::get_pcomm(dmmoab->mbiface, partnset, &comm);
339   }
340   else {
341     ierr = DMMoabSetParallelComm(*dmb, pcomm);CHKERRQ(ierr);
342   }
343 
344   /* do the remaining initializations for DMMoab */
345   dmmoab->bs = 1;
346   dmmoab->nfields = 1;
347 
348   /* set global ID tag handle */
349   if (!ltog_tag) {
350     merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr);
351   }
352   else {
353     ierr = DMMoabSetLocalToGlobalTag(*dmb, *ltog_tag);CHKERRQ(ierr);
354   }
355 
356   /* set the local range of entities (vertices) of interest */
357   if (range) {
358     ierr = DMMoabSetLocalVertices(*dmb, range);CHKERRQ(ierr);
359   }
360   PetscFunctionReturn(0);
361 }
362 
363 #undef __FUNCT__
364 #define __FUNCT__ "DMMoabSetParallelComm"
365 /*@
366   DMMoabSetParallelComm - Set the ParallelComm used with this DMMoab
367 
368   Collective on MPI_Comm
369 
370   Input Parameter:
371 . dm    - The DMMoab object being set
372 . pcomm - The ParallelComm being set on the DMMoab
373 
374   Level: beginner
375 
376 .keywords: DMMoab, create
377 @*/
378 PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm)
379 {
380   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
384   PetscValidPointer(pcomm,2);
385   dmmoab->pcomm = pcomm;
386   dmmoab->mbiface = pcomm->get_moab();
387   dmmoab->icreatedinstance = PETSC_FALSE;
388   PetscFunctionReturn(0);
389 }
390 
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "DMMoabGetParallelComm"
394 /*@
395   DMMoabGetParallelComm - Get the ParallelComm used with this DMMoab
396 
397   Collective on MPI_Comm
398 
399   Input Parameter:
400 . dm    - The DMMoab object being set
401 
402   Output Parameter:
403 . pcomm - The ParallelComm for the DMMoab
404 
405   Level: beginner
406 
407 .keywords: DMMoab, create
408 @*/
409 PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm)
410 {
411   PetscFunctionBegin;
412   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
413   *pcomm = ((DM_Moab*)(dm)->data)->pcomm;
414   PetscFunctionReturn(0);
415 }
416 
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "DMMoabSetInterface"
420 /*@
421   DMMoabSetInterface - Set the MOAB instance used with this DMMoab
422 
423   Collective on MPI_Comm
424 
425   Input Parameter:
426 . dm      - The DMMoab object being set
427 . mbiface - The MOAB instance being set on this DMMoab
428 
429   Level: beginner
430 
431 .keywords: DMMoab, create
432 @*/
433 PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *mbiface)
434 {
435   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
436 
437   PetscFunctionBegin;
438   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
439   PetscValidPointer(mbiface,2);
440   dmmoab->pcomm = NULL;
441   dmmoab->mbiface = mbiface;
442   dmmoab->icreatedinstance = PETSC_FALSE;
443   PetscFunctionReturn(0);
444 }
445 
446 
447 #undef __FUNCT__
448 #define __FUNCT__ "DMMoabGetInterface"
449 /*@
450   DMMoabGetInterface - Get the MOAB instance used with this DMMoab
451 
452   Collective on MPI_Comm
453 
454   Input Parameter:
455 . dm      - The DMMoab object being set
456 
457   Output Parameter:
458 . mbiface - The MOAB instance set on this DMMoab
459 
460   Level: beginner
461 
462 .keywords: DMMoab, create
463 @*/
464 PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbiface)
465 {
466   PetscErrorCode   ierr;
467   static PetscBool cite = PETSC_FALSE;
468 
469   PetscFunctionBegin;
470   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
471   ierr = PetscCitationsRegister("@techreport{tautges_moab:_2004,\n  type = {{SAND2004-1592}},\n  title = {{MOAB:} A Mesh-Oriented Database},  institution = {Sandia National Laboratories},\n  author = {Tautges, T. J. and Meyers, R. and Merkley, K. and Stimpson, C. and Ernst, C.},\n  year = {2004},  note = {Report}\n}\n",&cite);CHKERRQ(ierr);
472   *mbiface = ((DM_Moab*)dm->data)->mbiface;
473   PetscFunctionReturn(0);
474 }
475 
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "DMMoabSetLocalVertices"
479 /*@
480   DMMoabSetLocalVertices - Set the entities having DOFs on this DMMoab
481 
482   Collective on MPI_Comm
483 
484   Input Parameter:
485 . dm    - The DMMoab object being set
486 . range - The entities treated by this DMMoab
487 
488   Level: beginner
489 
490 .keywords: DMMoab, create
491 @*/
492 PetscErrorCode DMMoabSetLocalVertices(DM dm,moab::Range *range)
493 {
494   moab::ErrorCode merr;
495   PetscErrorCode  ierr;
496   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
497 
498   PetscFunctionBegin;
499   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
500   dmmoab->vlocal->clear();
501   dmmoab->vowned->clear();
502   dmmoab->vlocal->insert(range->begin(), range->end());
503   *dmmoab->vowned = *dmmoab->vlocal;
504   merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
505   *dmmoab->vghost = moab::subtract(*range, *dmmoab->vowned);
506   dmmoab->nloc=dmmoab->vowned->size();
507   dmmoab->nghost=dmmoab->vghost->size();
508   ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
509   PetscFunctionReturn(0);
510 }
511 
512 
513 #undef __FUNCT__
514 #define __FUNCT__ "DMMoabGetLocalVertices"
515 /*@
516   DMMoabGetLocalVertices - Get the entities having DOFs on this DMMoab
517 
518   Collective on MPI_Comm
519 
520   Input Parameter:
521 . dm    - The DMMoab object being set
522 
523   Output Parameter:
524 . owned - The owned vertex entities in this DMMoab
525 . ghost - The ghosted entities (non-owned) stored locally in this partition
526 
527   Level: beginner
528 
529 .keywords: DMMoab, create
530 @*/
531 PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost)
532 {
533   PetscFunctionBegin;
534   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
535   if (owned) *owned = *((DM_Moab*)dm->data)->vowned;
536   if (ghost) *ghost = *((DM_Moab*)dm->data)->vghost;
537   PetscFunctionReturn(0);
538 }
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "DMMoabGetLocalElements"
542 /*@
543   DMMoabGetLocalElements - Get the higher-dimensional entities that are locally owned
544 
545   Collective on MPI_Comm
546 
547   Input Parameter:
548 . dm    - The DMMoab object being set
549 
550   Output Parameter:
551 . range - The entities owned locally
552 
553   Level: beginner
554 
555 .keywords: DMMoab, create
556 @*/
557 PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range)
558 {
559   PetscFunctionBegin;
560   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
561   if (range) *range = *((DM_Moab*)dm->data)->elocal;
562   PetscFunctionReturn(0);
563 }
564 
565 
566 #undef __FUNCT__
567 #define __FUNCT__ "DMMoabSetLocalElements"
568 /*@
569   DMMoabSetLocalElements - Set the entities having DOFs on this DMMoab
570 
571   Collective on MPI_Comm
572 
573   Input Parameter:
574 . dm    - The DMMoab object being set
575 . range - The entities treated by this DMMoab
576 
577   Level: beginner
578 
579 .keywords: DMMoab, create
580 @*/
581 PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range)
582 {
583   moab::ErrorCode merr;
584   PetscErrorCode  ierr;
585   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
586 
587   PetscFunctionBegin;
588   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
589   dmmoab->elocal->clear();
590   dmmoab->eghost->clear();
591   dmmoab->elocal->insert(range->begin(), range->end());
592   PetscInfo2(dm, "Range size = %D; elocal size = %D.\n", range->size(), dmmoab->elocal->size());
593   merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
594   *dmmoab->eghost = moab::subtract(*range, *dmmoab->elocal);
595   dmmoab->neleloc=dmmoab->elocal->size();
596   ierr = MPI_Allreduce(&dmmoab->nele, &dmmoab->neleloc, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
597   PetscInfo2(dm, "Created %D local and %D glocal elements.\n", dmmoab->neleloc, dmmoab->nele);
598   PetscFunctionReturn(0);
599 }
600 
601 
602 #undef __FUNCT__
603 #define __FUNCT__ "DMMoabSetLocalToGlobalTag"
604 /*@
605   DMMoabSetLocalToGlobalTag - Set the tag used for local to global numbering
606 
607   Collective on MPI_Comm
608 
609   Input Parameter:
610 . dm      - The DMMoab object being set
611 . ltogtag - The MOAB tag used for local to global ids
612 
613   Level: beginner
614 
615 .keywords: DMMoab, create
616 @*/
617 PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag)
618 {
619   PetscFunctionBegin;
620   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
621   ((DM_Moab*)dm->data)->ltog_tag = ltogtag;
622   PetscFunctionReturn(0);
623 }
624 
625 
626 #undef __FUNCT__
627 #define __FUNCT__ "DMMoabGetLocalToGlobalTag"
628 /*@
629   DMMoabGetLocalToGlobalTag - Get the tag used for local to global numbering
630 
631   Collective on MPI_Comm
632 
633   Input Parameter:
634 . dm      - The DMMoab object being set
635 
636   Output Parameter:
637 . ltogtag - The MOAB tag used for local to global ids
638 
639   Level: beginner
640 
641 .keywords: DMMoab, create
642 @*/
643 PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag)
644 {
645   PetscFunctionBegin;
646   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
647   *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag;
648   PetscFunctionReturn(0);
649 }
650 
651 
652 #undef __FUNCT__
653 #define __FUNCT__ "DMMoabSetBlockSize"
654 /*@
655   DMMoabSetBlockSize - Set the block size used with this DMMoab
656 
657   Collective on MPI_Comm
658 
659   Input Parameter:
660 . dm - The DMMoab object being set
661 . bs - The block size used with this DMMoab
662 
663   Level: beginner
664 
665 .keywords: DMMoab, create
666 @*/
667 PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs)
668 {
669   PetscFunctionBegin;
670   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
671   ((DM_Moab*)dm->data)->bs = bs;
672   PetscFunctionReturn(0);
673 }
674 
675 
676 #undef __FUNCT__
677 #define __FUNCT__ "DMMoabGetBlockSize"
678 /*@
679   DMMoabGetBlockSize - Get the block size used with this DMMoab
680 
681   Collective on MPI_Comm
682 
683   Input Parameter:
684 . dm - The DMMoab object being set
685 
686   Output Parameter:
687 . bs - The block size used with this DMMoab
688 
689   Level: beginner
690 
691 .keywords: DMMoab, create
692 @*/
693 PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs)
694 {
695   PetscFunctionBegin;
696   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
697   *bs = ((DM_Moab*)dm->data)->bs;
698   PetscFunctionReturn(0);
699 }
700 
701 
702 #undef __FUNCT__
703 #define __FUNCT__ "DMMoabGetSize"
704 /*@
705   DMMoabGetSize - Get the global vertex size used with this DMMoab
706 
707   Collective on MPI_Comm
708 
709   Input Parameter:
710 . dm - The DMMoab object being set
711 
712   Output Parameter:
713 . ng - The global size of the DMMoab instance
714 
715   Level: beginner
716 
717 .keywords: DMMoab, create
718 @*/
719 PetscErrorCode DMMoabGetSize(DM dm,PetscInt *ng)
720 {
721   PetscFunctionBegin;
722   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
723   if(ng) *ng = ((DM_Moab*)dm->data)->n;
724   PetscFunctionReturn(0);
725 }
726 
727 
728 #undef __FUNCT__
729 #define __FUNCT__ "DMMoabGetLocalSize"
730 /*@
731   DMMoabGetLocalSize - Get the local and ghosted vertex size used with this DMMoab
732 
733   Collective on MPI_Comm
734 
735   Input Parameter:
736 . dm - The DMMoab object being set
737 
738   Output Parameter:
739 . nl - The local size of the DMMoab instance
740 . ng - The ghosted size of the DMMoab instance
741 
742   Level: beginner
743 
744 .keywords: DMMoab, create
745 @*/
746 PetscErrorCode DMMoabGetLocalSize(DM dm,PetscInt *nl,PetscInt *ng)
747 {
748   PetscFunctionBegin;
749   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
750   if(nl) *nl = ((DM_Moab*)dm->data)->nloc;
751   if(ng) *ng = ((DM_Moab*)dm->data)->nghost;
752   PetscFunctionReturn(0);
753 }
754 
755 
756 #undef __FUNCT__
757 #define __FUNCT__ "DMMoabSetFieldVector"
758 PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec)
759 {
760   DM_Moab        *dmmoab;
761   moab::Tag     vtag,ntag;
762   const PetscScalar *varray;
763   PetscScalar *farray;
764   moab::ErrorCode merr;
765   PetscErrorCode  ierr;
766   PetscInt doff;
767   std::string tag_name;
768   moab::Range::iterator iter;
769 
770   PetscFunctionBegin;
771   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
772   dmmoab = (DM_Moab*)(dm)->data;
773 
774   /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
775   merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
776                                           moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr);
777 
778   ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr);
779 
780   merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
781   if (!tag_name.length() && merr !=moab::MB_SUCCESS) {
782     ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr);
783 
784     for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) {
785       moab::EntityHandle vtx = (*iter);
786 
787       /* get field dof index */
788       ierr = DMMoabGetFieldDof(dm, vtx, ifield, &doff);
789 
790       /* use the entity handle and the Dof index to set the right value */
791       merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr);
792     }
793     ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr);
794   }
795   else {
796     ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr);
797     /* we are using a MOAB Vec - directly copy the tag data to new one */
798     merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)farray);MBERRNM(merr);
799     merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
800     /* make sure the parallel exchange for ghosts are done appropriately */
801     merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr);
802     ierr = PetscFree(farray);CHKERRQ(ierr);
803   }
804   PetscFunctionReturn(0);
805 }
806 
807 
808 #undef __FUNCT__
809 #define __FUNCT__ "DMMoabSetGlobalFieldVector"
810 PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec)
811 {
812   DM_Moab        *dmmoab;
813   moab::Tag     vtag,ntag;
814   const PetscScalar   *varray;
815   PetscScalar   *farray;
816   moab::ErrorCode merr;
817   PetscErrorCode  ierr;
818   PetscSection section;
819   PetscInt i,doff,ifield;
820   std::string tag_name;
821   moab::Range::iterator iter;
822 
823   PetscFunctionBegin;
824   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
825   dmmoab = (DM_Moab*)(dm)->data;
826 
827   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
828 
829   /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */
830   ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr);
831   merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
832   if (!tag_name.length() && merr !=moab::MB_SUCCESS) {
833     /* not a MOAB vector - use VecGetSubVector to get the parts as needed */
834 
835     ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr);
836     for (ifield=0; ifield<dmmoab->nfields; ++ifield) {
837 
838       /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
839       merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
840                                             moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr);
841 
842       for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) {
843         moab::EntityHandle vtx = (*iter);
844 
845         /* get field dof index */
846         ierr = DMMoabGetFieldDof(dm, vtx, ifield, &doff);
847 
848         /* use the entity handle and the Dof index to set the right value */
849         merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr);
850       }
851     }
852     ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr);
853   }
854   else {
855     ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr);
856     ierr = PetscMalloc(dmmoab->nloc*dmmoab->bs*sizeof(PetscScalar),&varray);CHKERRQ(ierr);
857 
858     /* we are using a MOAB Vec - directly copy the tag data to new one */
859     merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr);
860     for (ifield=0; ifield<dmmoab->nfields; ++ifield) {
861 
862       /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
863       merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
864                                             moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr);
865 
866       /* we are using a MOAB Vec - directly copy the tag data to new one */
867       for(i=0; i < dmmoab->nloc; i++) {
868         farray[i] = varray[i*dmmoab->bs+ifield];
869       }
870 
871       merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
872       /* make sure the parallel exchange for ghosts are done appropriately */
873       merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr);
874     }
875     ierr = PetscFree(farray);CHKERRQ(ierr);
876     ierr = PetscFree(varray);CHKERRQ(ierr);
877   }
878   PetscFunctionReturn(0);
879 }
880 
881 
882 
883 #undef __FUNCT__
884 #define __FUNCT__ "DMMoabGetVertexCoordinates"
885 PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos)
886 {
887   DM_Moab         *dmmoab;
888   PetscErrorCode  ierr;
889   moab::ErrorCode merr;
890 
891   PetscFunctionBegin;
892   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
893   PetscValidPointer(conn,3);
894   dmmoab = (DM_Moab*)(dm)->data;
895 
896   if (!vpos) {
897     ierr = PetscMalloc(sizeof(PetscScalar)*nconn*3, &vpos);CHKERRQ(ierr);
898   }
899 
900   /* Get connectivity information in MOAB canonical ordering */
901   merr = dmmoab->mbiface->get_coords(conn, nconn, vpos);MBERRNM(merr);
902   PetscFunctionReturn(0);
903 }
904 
905 
906 #undef __FUNCT__
907 #define __FUNCT__ "DMMoabGetElementConnectivity"
908 PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn)
909 {
910   DM_Moab        *dmmoab;
911   const moab::EntityHandle *connect;
912   moab::ErrorCode merr;
913   PetscInt nnodes;
914 
915   PetscFunctionBegin;
916   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
917   PetscValidPointer(conn,4);
918   dmmoab = (DM_Moab*)(dm)->data;
919 
920   /* Get connectivity information in MOAB canonical ordering */
921   merr = dmmoab->mbiface->get_connectivity(ehandle, connect, nnodes);MBERRNM(merr);
922   if (conn) *conn=connect;
923   if (nconn) *nconn=nnodes;
924   PetscFunctionReturn(0);
925 }
926 
927 
928 #undef __FUNCT__
929 #define __FUNCT__ "DMMoabIsEntityOnBoundary"
930 PetscErrorCode DMMoabIsEntityOnBoundary(DM dm,const moab::EntityHandle ent,PetscBool* ent_on_boundary)
931 {
932   moab::EntityType etype;
933   DM_Moab         *dmmoab;
934   PetscInt         edim;
935 
936   PetscFunctionBegin;
937   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
938   PetscValidPointer(ent_on_boundary,3);
939   dmmoab = (DM_Moab*)(dm)->data;
940 
941   /* get the entity type and handle accordingly */
942   etype=dmmoab->mbiface->type_from_handle(ent);
943   if(etype >= moab::MBPOLYHEDRON) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Entity type on the boundary skin is invalid. EntityType = %D\n",etype);
944 
945   /* get the entity dimension */
946   edim=dmmoab->mbiface->dimension_from_handle(ent);
947 
948   *ent_on_boundary=PETSC_FALSE;
949   if(etype == moab::MBVERTEX && edim == 0) {
950     moab::Range::const_iterator giter = dmmoab->bndyvtx->find(ent);
951     if (giter != dmmoab->bndyvtx->end()) *ent_on_boundary=PETSC_TRUE;
952   }
953   else {
954     if (edim == dmmoab->dim) {  /* check the higher-dimensional elements first */
955       moab::Range::const_iterator geiter = dmmoab->bndyelems->find(ent);
956       if (geiter != dmmoab->bndyelems->end()) *ent_on_boundary=PETSC_TRUE;
957     }
958     else {                      /* next check the lower-dimensional faces */
959       moab::Range::const_iterator gfiter = dmmoab->bndyfaces->find(ent);
960       if (gfiter != dmmoab->bndyfaces->end()) *ent_on_boundary=PETSC_TRUE;
961     }
962   }
963   PetscFunctionReturn(0);
964 }
965 
966 
967 #undef __FUNCT__
968 #define __FUNCT__ "DMMoabCheckBoundaryVertices"
969 PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx)
970 {
971   DM_Moab        *dmmoab;
972   PetscInt       i;
973 
974   PetscFunctionBegin;
975   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
976   PetscValidPointer(cnt,3);
977   PetscValidPointer(isbdvtx,4);
978   dmmoab = (DM_Moab*)(dm)->data;
979 
980   for (i=0; i < nconn; ++i) {
981     moab::Range::const_iterator giter = dmmoab->bndyvtx->find(cnt[i]);
982     if (giter != dmmoab->bndyvtx->end()) isbdvtx[i] = PETSC_TRUE;
983     else isbdvtx[i] = PETSC_FALSE;
984   }
985   PetscFunctionReturn(0);
986 }
987 
988 
989 #undef __FUNCT__
990 #define __FUNCT__ "DMMoabGetBoundaryEntities"
991 PetscErrorCode DMMoabGetBoundaryEntities(DM dm,moab::Range *bdvtx,moab::Range* bdfaces,moab::Range* bdelems)
992 {
993   DM_Moab        *dmmoab;
994 
995   PetscFunctionBegin;
996   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
997   dmmoab = (DM_Moab*)(dm)->data;
998 
999   if (bdvtx)  *bdvtx = *dmmoab->bndyvtx;
1000   if (bdfaces)  *bdfaces = *dmmoab->bndyfaces;
1001   if (bdelems)  *bdfaces = *dmmoab->bndyelems;
1002   PetscFunctionReturn(0);
1003 }
1004 
1005 
1006 #undef __FUNCT__
1007 #define __FUNCT__ "DMMoabSetFields"
1008 PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields)
1009 {
1010   DM_Moab        *dmmoab;
1011 
1012   PetscFunctionBegin;
1013   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1014   dmmoab = (DM_Moab*)(dm)->data;
1015 
1016   dmmoab->fields = fields;
1017   dmmoab->nfields = nfields;
1018   PetscFunctionReturn(0);
1019 }
1020 
1021 
1022 #undef __FUNCT__
1023 #define __FUNCT__ "DMMoabGetFieldDof"
1024 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof)
1025 {
1026   PetscSection section;
1027   PetscInt gid;
1028   PetscErrorCode ierr;
1029   moab::ErrorCode merr;
1030   DM_Moab        *dmmoab;
1031 
1032   PetscFunctionBegin;
1033   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1034   dmmoab = (DM_Moab*)(dm)->data;
1035 
1036   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1037 
1038   /* first get the global ID for the point */
1039   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&point,1,&gid);MBERRNM(merr);
1040 
1041   /* get the dof value for the field */
1042   ierr = PetscSectionGetFieldDof(section, gid, field, dof);CHKERRQ(ierr);
1043   PetscFunctionReturn(0);
1044 }
1045 
1046 
1047 #undef __FUNCT__
1048 #define __FUNCT__ "DMMoabGetFieldDofs"
1049 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
1050 {
1051   PetscInt i,gid;
1052   PetscSection section;
1053   PetscErrorCode  ierr;
1054   moab::ErrorCode merr;
1055   DM_Moab        *dmmoab;
1056 
1057   PetscFunctionBegin;
1058   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1059   PetscValidPointer(points,2);
1060   dmmoab = (DM_Moab*)(dm)->data;
1061 
1062   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1063   if (!dof) {
1064     ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr);
1065   }
1066 
1067   /* first get the local indices */
1068   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr);
1069 
1070   for (i=0; i<npoints; ++i) {
1071     gid=dof[i];
1072     ierr = PetscSectionGetFieldDof(section, gid, field, &dof[i]);CHKERRQ(ierr);
1073   }
1074   PetscFunctionReturn(0);
1075 }
1076 
1077 
1078 #undef __FUNCT__
1079 #define __FUNCT__ "DMMoabGetFieldDofsLocal"
1080 PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
1081 {
1082   PetscInt i,offset;
1083   PetscErrorCode  ierr;
1084   DM_Moab        *dmmoab;
1085 
1086   PetscFunctionBegin;
1087   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1088   PetscValidPointer(points,2);
1089   dmmoab = (DM_Moab*)(dm)->data;
1090 
1091   if (!dof) {
1092     ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr);
1093   }
1094 
1095   if (dmmoab->bs > 1) {
1096     for (i=0; i<npoints; ++i)
1097       dof[i] = (points[i]-1)*dmmoab->bs+field;
1098   }
1099   else {
1100     offset = field*dmmoab->n; /* assume all fields have equal distribution */
1101     for (i=0; i<npoints; ++i)
1102       dof[i] = offset+points[i]-1;
1103   }
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 
1108 #undef __FUNCT__
1109 #define __FUNCT__ "DMMoabGetDofs"
1110 PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
1111 {
1112   PetscInt i,f,gid;
1113   PetscSection section;
1114   PetscErrorCode  ierr;
1115   moab::ErrorCode merr;
1116   DM_Moab        *dmmoab;
1117 
1118   PetscFunctionBegin;
1119   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1120   PetscValidPointer(points,2);
1121   dmmoab = (DM_Moab*)(dm)->data;
1122 
1123   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1124   if (!dof) {
1125     ierr = PetscMalloc(sizeof(PetscScalar)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr);
1126   }
1127 
1128   /* first get the local indices */
1129   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr);
1130 
1131   for (i=0; i<npoints; ++i) {
1132     gid=dof[i];
1133     for (f=0; f<dmmoab->nfields; ++f) {
1134       ierr = PetscSectionGetFieldDof(section, gid, f, &dof[i*dmmoab->nfields+f]);CHKERRQ(ierr);
1135     }
1136   }
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 
1141 #undef __FUNCT__
1142 #define __FUNCT__ "DMMoabGetDofsLocal"
1143 PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
1144 {
1145   PetscInt        i,f,offset;
1146   PetscErrorCode  ierr;
1147   DM_Moab        *dmmoab;
1148 
1149   PetscFunctionBegin;
1150   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1151   PetscValidPointer(points,2);
1152   dmmoab = (DM_Moab*)(dm)->data;
1153 
1154   if (!dof) {
1155     ierr = PetscMalloc(sizeof(PetscScalar)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr);
1156   }
1157 
1158   if (dmmoab->bs > 1) {
1159     for (f=0; f<dmmoab->nfields; ++f)
1160       for (i=0; i<npoints; ++i)
1161         dof[i*dmmoab->nfields+f] = (points[i]-1)*dmmoab->bs+f;
1162   }
1163   else {
1164     for (f=0; f<dmmoab->nfields; ++f) {
1165       offset = f*dmmoab->n;     /* assume all fields have equal distribution - say all vertex based */
1166       for (i=0; i<npoints; ++i)
1167         dof[i*dmmoab->nfields+f] = offset+points[i]-1;
1168     }
1169   }
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 
1174 #undef __FUNCT__
1175 #define __FUNCT__ "DMMoabGetDofsBlocked"
1176 PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
1177 {
1178   PetscInt i,gid,dofindx;
1179   PetscSection section;
1180   PetscErrorCode  ierr;
1181   moab::ErrorCode merr;
1182   DM_Moab        *dmmoab;
1183 
1184   PetscFunctionBegin;
1185   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1186   PetscValidPointer(points,2);
1187   dmmoab = (DM_Moab*)(dm)->data;
1188 
1189   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1190   if (!dof) {
1191     ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr);
1192   }
1193 
1194   /* first get the local indices */
1195   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr);
1196 
1197   for (i=0; i<npoints; ++i) {
1198     gid=dof[i];
1199     ierr = PetscSectionGetFieldDof(section, gid, 0, &dofindx);CHKERRQ(ierr);
1200     if (dmmoab->bs > 1)  dof[i]=dofindx/dmmoab->bs;
1201     else dof[i]=dofindx;
1202   }
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 
1207 #undef __FUNCT__
1208 #define __FUNCT__ "DMMoabGetDofsBlockedLocal"
1209 PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
1210 {
1211   PetscInt        i;
1212   PetscErrorCode  ierr;
1213 
1214   PetscFunctionBegin;
1215   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1216   PetscValidPointer(points,2);
1217 
1218   if (!dof) {
1219     ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr);
1220   }
1221 
1222   for (i=0; i<npoints; ++i)
1223     dof[i] = points[i]-1;
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 #undef __FUNCT__
1228 #define __FUNCT__ "DMMoab_GetWriteOptions_Private"
1229 PetscErrorCode DMMoab_GetWriteOptions_Private(PetscInt fsetid, PetscInt numproc, PetscInt dim, MoabWriteMode mode, PetscInt dbglevel, const char* extra_opts, const char** write_opts)
1230 {
1231   std::ostringstream str;
1232 
1233   PetscFunctionBegin;
1234 
1235   // do parallel read unless only one processor
1236   if (numproc > 1) {
1237     str << "PARALLEL=" << mode << ";";
1238     if (fsetid>=0) str << "PARALLEL_COMM=" << fsetid << ";";
1239   }
1240 
1241   if (dbglevel)
1242     str << "CPUTIME;DEBUG_IO=" << dbglevel << ";";
1243 
1244   if (extra_opts)
1245     str << extra_opts ;
1246 
1247   *write_opts = str.str().c_str();
1248   PetscFunctionReturn(0);
1249 }
1250 
1251 
1252 #undef __FUNCT__
1253 #define __FUNCT__ "DMMoabOutput"
1254 PetscErrorCode DMMoabOutput(DM dm,const char* filename,const char* usrwriteopts)
1255 {
1256   DM_Moab        *dmmoab;
1257   PetscInt       dbglevel=0;
1258   const char *writeopts;
1259 
1260   PetscErrorCode ierr;
1261   moab::ErrorCode merr;
1262 
1263   PetscFunctionBegin;
1264   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1265   dmmoab = (DM_Moab*)(dm)->data;
1266 
1267   PetscBarrier((PetscObject)dm);
1268 
1269   /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */
1270   ierr  = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab");
1271   ierr  = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr);
1272   ierr  = PetscOptionsEnd();
1273 
1274   /* add mesh loading options specific to the DM */
1275   ierr = DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, MOAB_PARWOPTS_WRITE_PART, dbglevel, usrwriteopts, &writeopts);CHKERRQ(ierr);
1276   PetscInfo2(dm, "Writing file %s with options: %s\n",filename,writeopts);
1277 
1278   /* output file, using parallel write */
1279   merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1);MBERRVM(dmmoab->mbiface,"Writing output of DMMoab failed.",merr);
1280   PetscFunctionReturn(0);
1281 }
1282 
1283