xref: /petsc/src/dm/impls/moab/dmmoab.cxx (revision 6935a3f3eabe640fe4e66f153e1295282b78a4cc)
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 
31   ierr = PetscFree(dmmoab->isbndyvtx);CHKERRQ(ierr);
32   ierr = PetscFree(dmmoab->isbndyfaces);CHKERRQ(ierr);
33   ierr = PetscFree(dmmoab->isbndyelems);CHKERRQ(ierr);
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     moab::Range bndyfaces, bndyvtx, bndyelems;
187     /* get the skin vertices of boundary faces for the current partition and then filter
188        the local, boundary faces, vertices and elements alone via PSTATUS flags;
189        this should not give us any ghosted boundary, but if user needs such a functionality
190        it would be easy to add it based on the find_skin query below */
191     moab::Skinner skinner(dmmoab->mbiface);
192 
193     /* get the entities on the skin - only the faces */
194     merr = skinner.find_skin(dmmoab->fileset, *dmmoab->elocal, false, bndyfaces, NULL, false, true, false, false);MBERRNM(merr); // 'false' param indicates we want faces back, not vertices
195 
196     /* filter all the non-owned and shared entities out of the list */
197     merr = dmmoab->pcomm->filter_pstatus(bndyfaces,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
198     merr = dmmoab->pcomm->filter_pstatus(bndyfaces,PSTATUS_SHARED,PSTATUS_NOT);MBERRNM(merr);
199 
200     /* get all the nodes via connectivity and the parent elements via adjacency information */
201     merr = dmmoab->mbiface->get_connectivity(bndyfaces, bndyvtx, false);MBERRNM(ierr);
202     merr = dmmoab->mbiface->get_adjacencies(bndyfaces, dmmoab->dim, false, bndyelems, moab::Interface::UNION);MBERRNM(ierr);
203     PetscInfo3(dm, "Found %D boundary vertices, %D boundary faces and %D boundary elements.\n", bndyvtx.size(), bndyvtx.size(), bndyelems.size());
204 
205     /* cache a bit-vector for easy query */
206     ierr = PetscMalloc(sizeof(PetscBool)*((PetscInt)(*bndyvtx.rbegin())+1),&dmmoab->isbndyvtx);CHKERRQ(ierr);
207     ierr = PetscMemzero(dmmoab->isbndyvtx,sizeof(PetscBool)*((PetscInt)(*bndyvtx.rbegin())+1));CHKERRQ(ierr);
208     for(moab::Range::iterator iter = bndyvtx.begin(); iter != bndyvtx.end(); iter++) {
209       dmmoab->isbndyvtx[(PetscInt)*iter]=PETSC_TRUE;
210     }
211 
212     ierr = PetscMalloc(sizeof(PetscBool)*((PetscInt)(*bndyelems.rbegin())+1),&dmmoab->isbndyelems);CHKERRQ(ierr);
213     ierr = PetscMemzero(dmmoab->isbndyelems,sizeof(PetscBool)*((PetscInt)(*bndyelems.rbegin())+1));CHKERRQ(ierr);
214     for(moab::Range::iterator iter = bndyelems.begin(); iter != bndyelems.end(); iter++) {
215       dmmoab->isbndyelems[(PetscInt)*iter]=PETSC_TRUE;
216     }
217 
218     ierr = PetscMalloc(sizeof(PetscBool)*((PetscInt)(*bndyfaces.rbegin())+1),&dmmoab->isbndyfaces);CHKERRQ(ierr);
219     ierr = PetscMemzero(dmmoab->isbndyfaces,sizeof(PetscBool)*((PetscInt)(*bndyfaces.rbegin())+1));CHKERRQ(ierr);
220     for(moab::Range::iterator iter = bndyfaces.begin(); iter != bndyfaces.end(); iter++) {
221       dmmoab->isbndyfaces[(PetscInt)*iter]=PETSC_TRUE;
222     }
223   }
224   PetscFunctionReturn(0);
225 }
226 
227 
228 #undef __FUNCT__
229 #define __FUNCT__ "DMCreate_Moab"
230 PETSC_EXTERN PetscErrorCode DMCreate_Moab(DM dm)
231 {
232   PetscErrorCode ierr;
233 
234   PetscFunctionBegin;
235   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
236   ierr = PetscNewLog(dm,&dm->data);CHKERRQ(ierr);
237 
238   ((DM_Moab*)dm->data)->bs = 1;
239   ((DM_Moab*)dm->data)->nfields = 1;
240   ((DM_Moab*)dm->data)->n = 0;
241   ((DM_Moab*)dm->data)->nloc = 0;
242   ((DM_Moab*)dm->data)->nele = 0;
243   ((DM_Moab*)dm->data)->neleloc = 0;
244   ((DM_Moab*)dm->data)->nghost = 0;
245   ((DM_Moab*)dm->data)->ltog_map = PETSC_NULL;
246   ((DM_Moab*)dm->data)->ltog_sendrecv = PETSC_NULL;
247 
248   ((DM_Moab*)dm->data)->vlocal = new moab::Range();
249   ((DM_Moab*)dm->data)->vowned = new moab::Range();
250   ((DM_Moab*)dm->data)->vghost = new moab::Range();
251   ((DM_Moab*)dm->data)->elocal = new moab::Range();
252   ((DM_Moab*)dm->data)->eghost = new moab::Range();
253 
254   dm->ops->createglobalvector              = DMCreateGlobalVector_Moab;
255   dm->ops->createlocalvector               = DMCreateLocalVector_Moab;
256   dm->ops->creatematrix                    = DMCreateMatrix_Moab;
257   dm->ops->setup                           = DMSetUp_Moab;
258   dm->ops->destroy                         = DMDestroy_Moab;
259   dm->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Moab;
260   dm->ops->globaltolocalend                = DMGlobalToLocalEnd_Moab;
261   dm->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Moab;
262   dm->ops->localtoglobalend                = DMLocalToGlobalEnd_Moab;
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "DMMoabCreate"
268 /*@
269   DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance
270 
271   Collective on MPI_Comm
272 
273   Input Parameter:
274 . comm - The communicator for the DMMoab object
275 
276   Output Parameter:
277 . dmb  - The DMMoab object
278 
279   Level: beginner
280 
281 .keywords: DMMoab, create
282 @*/
283 PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *dmb)
284 {
285   PetscErrorCode ierr;
286 
287   PetscFunctionBegin;
288   PetscValidPointer(dmb,2);
289   ierr = DMCreate(comm, dmb);CHKERRQ(ierr);
290   ierr = DMSetType(*dmb, DMMOAB);CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "DMMoabCreateMoab"
296 /*@
297   DMMoabCreate - Creates a DMMoab object, optionally from an instance and other data
298 
299   Collective on MPI_Comm
300 
301   Input Parameter:
302 . comm - The communicator for the DMMoab object
303 . mbiface - (ptr to) the MOAB Instance; if passed in NULL, MOAB instance is created inside PETSc, and destroyed
304          along with the DMMoab
305 . pcomm - (ptr to) a ParallelComm; if NULL, creates one internally for the whole communicator
306 . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag
307 . range - If non-NULL, contains range of entities to which DOFs will be assigned
308 
309   Output Parameter:
310 . dmb  - The DMMoab object
311 
312   Level: intermediate
313 
314 .keywords: DMMoab, create
315 @*/
316 PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *dmb)
317 {
318   PetscErrorCode ierr;
319   moab::ErrorCode merr;
320   moab::EntityHandle partnset;
321   PetscInt rank, nprocs;
322   DM_Moab        *dmmoab;
323 
324   PetscFunctionBegin;
325   PetscValidPointer(dmb,6);
326   ierr = DMMoabCreate(comm, dmb);CHKERRQ(ierr);
327   dmmoab = (DM_Moab*)(*dmb)->data;
328 
329   if (!mbiface) {
330     dmmoab->mbiface = new moab::Core();
331     dmmoab->icreatedinstance = PETSC_TRUE;
332   }
333   else {
334     dmmoab->mbiface = mbiface;
335     dmmoab->icreatedinstance = PETSC_FALSE;
336   }
337 
338   /* create a fileset to store the hierarchy of entities belonging to current DM */
339   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_ORDERED, dmmoab->fileset);MBERR("Creating file set failed", merr);
340 
341   if (!pcomm) {
342     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
343     ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
344 
345     /* Create root sets for each mesh.  Then pass these
346        to the load_file functions to be populated. */
347     merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, partnset);MBERR("Creating partition set failed", merr);
348 
349     /* Create the parallel communicator object with the partition handle associated with MOAB */
350     dmmoab->pcomm = moab::ParallelComm::get_pcomm(dmmoab->mbiface, partnset, &comm);
351   }
352   else {
353     ierr = DMMoabSetParallelComm(*dmb, pcomm);CHKERRQ(ierr);
354   }
355 
356   /* do the remaining initializations for DMMoab */
357   dmmoab->bs = 1;
358   dmmoab->nfields = 1;
359 
360   /* set global ID tag handle */
361   if (!ltog_tag) {
362     merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr);
363   }
364   else {
365     ierr = DMMoabSetLocalToGlobalTag(*dmb, *ltog_tag);CHKERRQ(ierr);
366   }
367 
368   /* set the local range of entities (vertices) of interest */
369   if (range) {
370     ierr = DMMoabSetLocalVertices(*dmb, range);CHKERRQ(ierr);
371   }
372   PetscFunctionReturn(0);
373 }
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "DMMoabSetParallelComm"
377 /*@
378   DMMoabSetParallelComm - Set the ParallelComm used with this DMMoab
379 
380   Collective on MPI_Comm
381 
382   Input Parameter:
383 . dm    - The DMMoab object being set
384 . pcomm - The ParallelComm being set on the DMMoab
385 
386   Level: beginner
387 
388 .keywords: DMMoab, create
389 @*/
390 PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm)
391 {
392   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
393 
394   PetscFunctionBegin;
395   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
396   PetscValidPointer(pcomm,2);
397   dmmoab->pcomm = pcomm;
398   dmmoab->mbiface = pcomm->get_moab();
399   dmmoab->icreatedinstance = PETSC_FALSE;
400   PetscFunctionReturn(0);
401 }
402 
403 
404 #undef __FUNCT__
405 #define __FUNCT__ "DMMoabGetParallelComm"
406 /*@
407   DMMoabGetParallelComm - Get the ParallelComm used with this DMMoab
408 
409   Collective on MPI_Comm
410 
411   Input Parameter:
412 . dm    - The DMMoab object being set
413 
414   Output Parameter:
415 . pcomm - The ParallelComm for the DMMoab
416 
417   Level: beginner
418 
419 .keywords: DMMoab, create
420 @*/
421 PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm)
422 {
423   PetscFunctionBegin;
424   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
425   *pcomm = ((DM_Moab*)(dm)->data)->pcomm;
426   PetscFunctionReturn(0);
427 }
428 
429 
430 #undef __FUNCT__
431 #define __FUNCT__ "DMMoabSetInterface"
432 /*@
433   DMMoabSetInterface - Set the MOAB instance used with this DMMoab
434 
435   Collective on MPI_Comm
436 
437   Input Parameter:
438 . dm      - The DMMoab object being set
439 . mbiface - The MOAB instance being set on this DMMoab
440 
441   Level: beginner
442 
443 .keywords: DMMoab, create
444 @*/
445 PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *mbiface)
446 {
447   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
448 
449   PetscFunctionBegin;
450   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
451   PetscValidPointer(mbiface,2);
452   dmmoab->pcomm = NULL;
453   dmmoab->mbiface = mbiface;
454   dmmoab->icreatedinstance = PETSC_FALSE;
455   PetscFunctionReturn(0);
456 }
457 
458 
459 #undef __FUNCT__
460 #define __FUNCT__ "DMMoabGetInterface"
461 /*@
462   DMMoabGetInterface - Get the MOAB instance used with this DMMoab
463 
464   Collective on MPI_Comm
465 
466   Input Parameter:
467 . dm      - The DMMoab object being set
468 
469   Output Parameter:
470 . mbiface - The MOAB instance set on this DMMoab
471 
472   Level: beginner
473 
474 .keywords: DMMoab, create
475 @*/
476 PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbiface)
477 {
478   PetscErrorCode   ierr;
479   static PetscBool cite = PETSC_FALSE;
480 
481   PetscFunctionBegin;
482   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
483   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);
484   *mbiface = ((DM_Moab*)dm->data)->mbiface;
485   PetscFunctionReturn(0);
486 }
487 
488 
489 #undef __FUNCT__
490 #define __FUNCT__ "DMMoabSetLocalVertices"
491 /*@
492   DMMoabSetLocalVertices - Set the entities having DOFs on this DMMoab
493 
494   Collective on MPI_Comm
495 
496   Input Parameter:
497 . dm    - The DMMoab object being set
498 . range - The entities treated by this DMMoab
499 
500   Level: beginner
501 
502 .keywords: DMMoab, create
503 @*/
504 PetscErrorCode DMMoabSetLocalVertices(DM dm,moab::Range *range)
505 {
506   moab::ErrorCode merr;
507   PetscErrorCode  ierr;
508   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
509 
510   PetscFunctionBegin;
511   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
512   dmmoab->vlocal->clear();
513   dmmoab->vowned->clear();
514   dmmoab->vlocal->insert(range->begin(), range->end());
515   *dmmoab->vowned = *dmmoab->vlocal;
516   merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
517   *dmmoab->vghost = moab::subtract(*range, *dmmoab->vowned);
518   dmmoab->nloc=dmmoab->vowned->size();
519   dmmoab->nghost=dmmoab->vghost->size();
520   ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
521   PetscFunctionReturn(0);
522 }
523 
524 
525 #undef __FUNCT__
526 #define __FUNCT__ "DMMoabGetAllVertices"
527 /*@
528   DMMoabGetAllVertices - Get the entities having DOFs on this DMMoab
529 
530   Collective on MPI_Comm
531 
532   Input Parameter:
533 . dm    - The DMMoab object being set
534 
535   Output Parameter:
536 . owned - The local vertex entities in this DMMoab = (owned+ghosted)
537 
538   Level: beginner
539 
540 .keywords: DMMoab, create
541 @*/
542 PetscErrorCode DMMoabGetAllVertices(DM dm,moab::Range *local)
543 {
544   PetscFunctionBegin;
545   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
546   if (local) *local = *((DM_Moab*)dm->data)->vlocal;
547   PetscFunctionReturn(0);
548 }
549 
550 
551 
552 #undef __FUNCT__
553 #define __FUNCT__ "DMMoabGetLocalVertices"
554 /*@
555   DMMoabGetLocalVertices - Get the entities having DOFs on this DMMoab
556 
557   Collective on MPI_Comm
558 
559   Input Parameter:
560 . dm    - The DMMoab object being set
561 
562   Output Parameter:
563 . owned - The owned vertex entities in this DMMoab
564 . ghost - The ghosted entities (non-owned) stored locally in this partition
565 
566   Level: beginner
567 
568 .keywords: DMMoab, create
569 @*/
570 PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost)
571 {
572   PetscFunctionBegin;
573   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
574   if (owned) *owned = *((DM_Moab*)dm->data)->vowned;
575   if (ghost) *ghost = *((DM_Moab*)dm->data)->vghost;
576   PetscFunctionReturn(0);
577 }
578 
579 #undef __FUNCT__
580 #define __FUNCT__ "DMMoabGetLocalElements"
581 /*@
582   DMMoabGetLocalElements - Get the higher-dimensional entities that are locally owned
583 
584   Collective on MPI_Comm
585 
586   Input Parameter:
587 . dm    - The DMMoab object being set
588 
589   Output Parameter:
590 . range - The entities owned locally
591 
592   Level: beginner
593 
594 .keywords: DMMoab, create
595 @*/
596 PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range)
597 {
598   PetscFunctionBegin;
599   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
600   if (range) *range = *((DM_Moab*)dm->data)->elocal;
601   PetscFunctionReturn(0);
602 }
603 
604 
605 #undef __FUNCT__
606 #define __FUNCT__ "DMMoabSetLocalElements"
607 /*@
608   DMMoabSetLocalElements - Set the entities having DOFs on this DMMoab
609 
610   Collective on MPI_Comm
611 
612   Input Parameter:
613 . dm    - The DMMoab object being set
614 . range - The entities treated by this DMMoab
615 
616   Level: beginner
617 
618 .keywords: DMMoab, create
619 @*/
620 PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range)
621 {
622   moab::ErrorCode merr;
623   PetscErrorCode  ierr;
624   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
625 
626   PetscFunctionBegin;
627   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
628   dmmoab->elocal->clear();
629   dmmoab->eghost->clear();
630   dmmoab->elocal->insert(range->begin(), range->end());
631   PetscInfo2(dm, "Range size = %D; elocal size = %D.\n", range->size(), dmmoab->elocal->size());
632   merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
633   *dmmoab->eghost = moab::subtract(*range, *dmmoab->elocal);
634   dmmoab->neleloc=dmmoab->elocal->size();
635   ierr = MPI_Allreduce(&dmmoab->nele, &dmmoab->neleloc, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
636   PetscInfo2(dm, "Created %D local and %D glocal elements.\n", dmmoab->neleloc, dmmoab->nele);
637   PetscFunctionReturn(0);
638 }
639 
640 
641 #undef __FUNCT__
642 #define __FUNCT__ "DMMoabSetLocalToGlobalTag"
643 /*@
644   DMMoabSetLocalToGlobalTag - Set the tag used for local to global numbering
645 
646   Collective on MPI_Comm
647 
648   Input Parameter:
649 . dm      - The DMMoab object being set
650 . ltogtag - The MOAB tag used for local to global ids
651 
652   Level: beginner
653 
654 .keywords: DMMoab, create
655 @*/
656 PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag)
657 {
658   PetscFunctionBegin;
659   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
660   ((DM_Moab*)dm->data)->ltog_tag = ltogtag;
661   PetscFunctionReturn(0);
662 }
663 
664 
665 #undef __FUNCT__
666 #define __FUNCT__ "DMMoabGetLocalToGlobalTag"
667 /*@
668   DMMoabGetLocalToGlobalTag - Get the tag used for local to global numbering
669 
670   Collective on MPI_Comm
671 
672   Input Parameter:
673 . dm      - The DMMoab object being set
674 
675   Output Parameter:
676 . ltogtag - The MOAB tag used for local to global ids
677 
678   Level: beginner
679 
680 .keywords: DMMoab, create
681 @*/
682 PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag)
683 {
684   PetscFunctionBegin;
685   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
686   *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag;
687   PetscFunctionReturn(0);
688 }
689 
690 
691 #undef __FUNCT__
692 #define __FUNCT__ "DMMoabSetBlockSize"
693 /*@
694   DMMoabSetBlockSize - Set the block size used with this DMMoab
695 
696   Collective on MPI_Comm
697 
698   Input Parameter:
699 . dm - The DMMoab object being set
700 . bs - The block size used with this DMMoab
701 
702   Level: beginner
703 
704 .keywords: DMMoab, create
705 @*/
706 PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs)
707 {
708   PetscFunctionBegin;
709   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
710   ((DM_Moab*)dm->data)->bs = bs;
711   PetscFunctionReturn(0);
712 }
713 
714 
715 #undef __FUNCT__
716 #define __FUNCT__ "DMMoabGetBlockSize"
717 /*@
718   DMMoabGetBlockSize - Get the block size used with this DMMoab
719 
720   Collective on MPI_Comm
721 
722   Input Parameter:
723 . dm - The DMMoab object being set
724 
725   Output Parameter:
726 . bs - The block size used with this DMMoab
727 
728   Level: beginner
729 
730 .keywords: DMMoab, create
731 @*/
732 PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs)
733 {
734   PetscFunctionBegin;
735   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
736   *bs = ((DM_Moab*)dm->data)->bs;
737   PetscFunctionReturn(0);
738 }
739 
740 
741 #undef __FUNCT__
742 #define __FUNCT__ "DMMoabGetSize"
743 /*@
744   DMMoabGetSize - Get the global vertex size used with this DMMoab
745 
746   Collective on MPI_Comm
747 
748   Input Parameter:
749 . dm - The DMMoab object being set
750 
751   Output Parameter:
752 . ng - The global size of the DMMoab instance
753 
754   Level: beginner
755 
756 .keywords: DMMoab, create
757 @*/
758 PetscErrorCode DMMoabGetSize(DM dm,PetscInt *ng)
759 {
760   PetscFunctionBegin;
761   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
762   if(ng) *ng = ((DM_Moab*)dm->data)->n;
763   PetscFunctionReturn(0);
764 }
765 
766 
767 #undef __FUNCT__
768 #define __FUNCT__ "DMMoabGetLocalSize"
769 /*@
770   DMMoabGetLocalSize - Get the local and ghosted vertex size used with this DMMoab
771 
772   Collective on MPI_Comm
773 
774   Input Parameter:
775 . dm - The DMMoab object being set
776 
777   Output Parameter:
778 . nl - The local size of the DMMoab instance
779 . ng - The ghosted size of the DMMoab instance
780 
781   Level: beginner
782 
783 .keywords: DMMoab, create
784 @*/
785 PetscErrorCode DMMoabGetLocalSize(DM dm,PetscInt *nl,PetscInt *ng)
786 {
787   PetscFunctionBegin;
788   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
789   if(nl) *nl = ((DM_Moab*)dm->data)->nloc;
790   if(ng) *ng = ((DM_Moab*)dm->data)->nghost;
791   PetscFunctionReturn(0);
792 }
793 
794 
795 #undef __FUNCT__
796 #define __FUNCT__ "DMMoabGetDimension"
797 /*@
798   DMMoabGetDimension - Get the dimension of the DM Mesh
799 
800   Collective on MPI_Comm
801 
802   Input Parameter:
803 . dm - The DMMoab object being set
804 
805   Output Parameter:
806 . dim - The dimension of DM
807 
808   Level: beginner
809 
810 .keywords: DMMoab, create
811 @*/
812 PetscErrorCode DMMoabGetDimension(DM dm,PetscInt *dim)
813 {
814   PetscFunctionBegin;
815   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
816   *dim = ((DM_Moab*)dm->data)->dim;
817   PetscFunctionReturn(0);
818 }
819 
820 
821 
822 #undef __FUNCT__
823 #define __FUNCT__ "DMMoabSetFieldVector"
824 PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec)
825 {
826   DM_Moab        *dmmoab;
827   moab::Tag     vtag,ntag;
828   const PetscScalar *varray;
829   PetscScalar *farray;
830   moab::ErrorCode merr;
831   PetscErrorCode  ierr;
832   std::string tag_name;
833 
834   PetscFunctionBegin;
835   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
836   dmmoab = (DM_Moab*)(dm)->data;
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_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr);
841 
842   ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr);
843 
844   merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
845   if (!tag_name.length() && merr !=moab::MB_SUCCESS) {
846     ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr);
847     /* use the entity handle and the Dof index to set the right value */
848     merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)varray);MBERRNM(merr);
849     ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr);
850   }
851   else {
852     ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr);
853     /* we are using a MOAB Vec - directly copy the tag data to new one */
854     merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)farray);MBERRNM(merr);
855     merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
856     /* make sure the parallel exchange for ghosts are done appropriately */
857     ierr = PetscFree(farray);CHKERRQ(ierr);
858   }
859   merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vowned);MBERRNM(merr);
860   PetscFunctionReturn(0);
861 }
862 
863 
864 #undef __FUNCT__
865 #define __FUNCT__ "DMMoabSetGlobalFieldVector"
866 PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec)
867 {
868   DM_Moab        *dmmoab;
869   moab::Tag     vtag,ntag;
870   const PetscScalar   *varray;
871   PetscScalar   *farray;
872   moab::ErrorCode merr;
873   PetscErrorCode  ierr;
874   PetscSection section;
875   PetscInt i,ifield;
876   std::string tag_name;
877   moab::Range::iterator iter;
878 
879   PetscFunctionBegin;
880   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
881   dmmoab = (DM_Moab*)(dm)->data;
882 
883   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
884 
885   /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */
886   ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr);
887   merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
888   ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr);
889   if (!tag_name.length() && merr !=moab::MB_SUCCESS) {
890     /* not a MOAB vector - use VecGetSubVector to get the parts as needed */
891 
892     ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr);
893     for (ifield=0; ifield<dmmoab->nfields; ++ifield) {
894 
895       /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
896       merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
897                                             moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr);
898 
899       for(i=0;i<dmmoab->nloc;i++) {
900         if (dmmoab->bs == 1)
901           farray[i]=varray[ifield*dmmoab->nloc+i];
902         else
903           farray[i]=varray[i*dmmoab->nfields+ifield];
904       }
905 
906       /* use the entity handle and the Dof index to set the right value */
907       merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
908     }
909     ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr);
910   }
911   else {
912     ierr = PetscMalloc(dmmoab->nloc*dmmoab->bs*sizeof(PetscScalar),&varray);CHKERRQ(ierr);
913 
914     /* we are using a MOAB Vec - directly copy the tag data to new one */
915     merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr);
916     for (ifield=0; ifield<dmmoab->nfields; ++ifield) {
917 
918       /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
919       merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
920                                             moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr);
921 
922       /* we are using a MOAB Vec - directly copy the tag data to new one */
923       for(i=0; i < dmmoab->nloc; i++) {
924         farray[i] = varray[i*dmmoab->bs+ifield];
925       }
926 
927       merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
928       /* make sure the parallel exchange for ghosts are done appropriately */
929       merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr);
930     }
931     ierr = PetscFree(varray);CHKERRQ(ierr);
932   }
933   ierr = PetscFree(farray);CHKERRQ(ierr);
934   PetscFunctionReturn(0);
935 }
936 
937 
938 
939 #undef __FUNCT__
940 #define __FUNCT__ "DMMoabGetVertexCoordinates"
941 PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos)
942 {
943   DM_Moab         *dmmoab;
944   PetscErrorCode  ierr;
945   moab::ErrorCode merr;
946 
947   PetscFunctionBegin;
948   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
949   PetscValidPointer(conn,3);
950   dmmoab = (DM_Moab*)(dm)->data;
951 
952   if (!vpos) {
953     ierr = PetscMalloc(sizeof(PetscScalar)*nconn*3, &vpos);CHKERRQ(ierr);
954   }
955 
956   /* Get connectivity information in MOAB canonical ordering */
957   merr = dmmoab->mbiface->get_coords(conn, nconn, vpos);MBERRNM(merr);
958   PetscFunctionReturn(0);
959 }
960 
961 
962 #undef __FUNCT__
963 #define __FUNCT__ "DMMoabGetVertexConnectivity"
964 PetscErrorCode DMMoabGetVertexConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn, moab::EntityHandle **conn)
965 {
966   DM_Moab        *dmmoab;
967   std::vector<moab::EntityHandle> adj_entities,connect;
968   PetscErrorCode  ierr;
969   moab::ErrorCode merr;
970 
971   PetscFunctionBegin;
972   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
973   PetscValidPointer(conn,4);
974   dmmoab = (DM_Moab*)(dm)->data;
975 
976   /* Get connectivity information in MOAB canonical ordering */
977   merr = dmmoab->mbiface->get_adjacencies(&ehandle, 1, 1, true, adj_entities, moab::Interface::UNION);MBERRNM(merr);
978   merr = dmmoab->mbiface->get_connectivity(&adj_entities[0],adj_entities.size(),connect);MBERRNM(merr);
979 
980 #if 0
981   for(unsigned int jter = 0; jter < connect.size(); jter++) {
982     PetscPrintf(PETSC_COMM_SELF,"Handle=%D\tAdj_Size=%D\tAdj_Entity=%D\n",ehandle,connect.size(),connect[jter]);
983   }
984 #endif
985 
986   if (conn) {
987     ierr = PetscMalloc(sizeof(moab::EntityHandle)*connect.size(), conn);CHKERRQ(ierr);
988     ierr = PetscMemcpy(*conn, &connect[0], sizeof(moab::EntityHandle)*connect.size());CHKERRQ(ierr);
989   }
990   if (nconn) *nconn=connect.size();
991   PetscFunctionReturn(0);
992 }
993 
994 
995 #undef __FUNCT__
996 #define __FUNCT__ "DMMoabRestoreVertexConnectivity"
997 PetscErrorCode DMMoabRestoreVertexConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn, moab::EntityHandle **conn)
998 {
999   PetscErrorCode  ierr;
1000 
1001   PetscFunctionBegin;
1002   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1003   PetscValidPointer(conn,4);
1004 
1005   if (conn) {
1006     ierr = PetscFree(*conn);CHKERRQ(ierr);
1007   }
1008   if (nconn) *nconn=0;
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 
1013 
1014 #undef __FUNCT__
1015 #define __FUNCT__ "DMMoabGetElementConnectivity"
1016 PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn)
1017 {
1018   DM_Moab        *dmmoab;
1019   const moab::EntityHandle *connect;
1020   moab::ErrorCode merr;
1021   PetscInt nnodes;
1022 
1023   PetscFunctionBegin;
1024   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1025   PetscValidPointer(conn,4);
1026   dmmoab = (DM_Moab*)(dm)->data;
1027 
1028   /* Get connectivity information in MOAB canonical ordering */
1029   merr = dmmoab->mbiface->get_connectivity(ehandle, connect, nnodes);MBERRNM(merr);
1030   if (conn) *conn=connect;
1031   if (nconn) *nconn=nnodes;
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 
1036 #undef __FUNCT__
1037 #define __FUNCT__ "DMMoabIsEntityOnBoundary"
1038 PetscErrorCode DMMoabIsEntityOnBoundary(DM dm,const moab::EntityHandle ent,PetscBool* ent_on_boundary)
1039 {
1040   moab::EntityType etype;
1041   DM_Moab         *dmmoab;
1042   PetscInt         edim;
1043 
1044   PetscFunctionBegin;
1045   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1046   PetscValidPointer(ent_on_boundary,3);
1047   dmmoab = (DM_Moab*)(dm)->data;
1048 
1049   /* get the entity type and handle accordingly */
1050   etype=dmmoab->mbiface->type_from_handle(ent);
1051   if(etype >= moab::MBPOLYHEDRON) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Entity type on the boundary skin is invalid. EntityType = %D\n",etype);
1052 
1053   /* get the entity dimension */
1054   edim=dmmoab->mbiface->dimension_from_handle(ent);
1055 
1056   *ent_on_boundary=PETSC_FALSE;
1057   if(etype == moab::MBVERTEX && edim == 0) {
1058     if (ent < (*dmmoab->vlocal)[0] || ent > (*dmmoab->vlocal)[dmmoab->nloc-1]) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid boundary vertex entity handle: %D\n",ent);
1059     *ent_on_boundary=dmmoab->isbndyvtx[(PetscInt)ent];
1060   }
1061   else {
1062     if (edim == dmmoab->dim) {  /* check the higher-dimensional elements first */
1063       if (ent < (*dmmoab->elocal)[0] || ent > (*dmmoab->elocal)[dmmoab->neleloc-1]) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid boundary element entity handle: %D\n",ent);
1064       *ent_on_boundary=dmmoab->isbndyelems[(PetscInt)ent];
1065     }
1066     else {                      /* next check the lower-dimensional faces */
1067       /* how do we check the bounds before accessing ? will segfault for non-boundary faces */
1068       *ent_on_boundary=dmmoab->isbndyfaces[(PetscInt)ent];
1069     }
1070   }
1071   PetscFunctionReturn(0);
1072 }
1073 
1074 
1075 #undef __FUNCT__
1076 #define __FUNCT__ "DMMoabCheckBoundaryVertices"
1077 PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx)
1078 {
1079   DM_Moab        *dmmoab;
1080   PetscInt       i;
1081 
1082   PetscFunctionBegin;
1083   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1084   PetscValidPointer(cnt,3);
1085   PetscValidPointer(isbdvtx,4);
1086   dmmoab = (DM_Moab*)(dm)->data;
1087 
1088   for (i=0; i < nconn; ++i) {
1089     isbdvtx[i]=dmmoab->isbndyvtx[(PetscInt)cnt[i]];
1090   }
1091   PetscFunctionReturn(0);
1092 }
1093 
1094 
1095 #undef __FUNCT__
1096 #define __FUNCT__ "DMMoabGetBoundaryMarkers"
1097 PetscErrorCode DMMoabGetBoundaryMarkers(DM dm,PetscBool **bdvtx,PetscBool** bdelems,PetscBool** bdfaces)
1098 {
1099   DM_Moab        *dmmoab;
1100 
1101   PetscFunctionBegin;
1102   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1103   dmmoab = (DM_Moab*)(dm)->data;
1104 
1105   if (bdvtx)  *bdvtx = dmmoab->isbndyvtx;
1106   if (bdfaces)  *bdfaces = dmmoab->isbndyfaces;
1107   if (bdelems)  *bdfaces = dmmoab->isbndyelems;
1108   PetscFunctionReturn(0);
1109 }
1110 
1111 
1112 #undef __FUNCT__
1113 #define __FUNCT__ "DMMoabSetFields"
1114 PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields)
1115 {
1116   DM_Moab        *dmmoab;
1117 
1118   PetscFunctionBegin;
1119   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1120   dmmoab = (DM_Moab*)(dm)->data;
1121 
1122   dmmoab->fields = fields;
1123   dmmoab->nfields = nfields;
1124   PetscFunctionReturn(0);
1125 }
1126 
1127 
1128 #undef __FUNCT__
1129 #define __FUNCT__ "DMMoabGetFieldDof"
1130 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof)
1131 {
1132   PetscSection section;
1133   PetscInt gid;
1134   PetscErrorCode ierr;
1135   moab::ErrorCode merr;
1136   DM_Moab        *dmmoab;
1137 
1138   PetscFunctionBegin;
1139   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1140   dmmoab = (DM_Moab*)(dm)->data;
1141 
1142   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1143 
1144   /* first get the global ID for the point */
1145   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&point,1,&gid);MBERRNM(merr);
1146 
1147   /* get the dof value for the field */
1148   ierr = PetscSectionGetFieldDof(section, gid, field, dof);CHKERRQ(ierr);
1149   PetscFunctionReturn(0);
1150 }
1151 
1152 
1153 #undef __FUNCT__
1154 #define __FUNCT__ "DMMoabGetFieldDofs"
1155 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
1156 {
1157   PetscInt i,gid;
1158   PetscSection section;
1159   PetscErrorCode  ierr;
1160   moab::ErrorCode merr;
1161   DM_Moab        *dmmoab;
1162 
1163   PetscFunctionBegin;
1164   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1165   PetscValidPointer(points,2);
1166   dmmoab = (DM_Moab*)(dm)->data;
1167 
1168   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1169   if (!dof) {
1170     ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr);
1171   }
1172 
1173   /* first get the local indices */
1174   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr);
1175 
1176   for (i=0; i<npoints; ++i) {
1177     gid=dof[i];
1178     ierr = PetscSectionGetFieldDof(section, gid, field, &dof[i]);CHKERRQ(ierr);
1179   }
1180   PetscFunctionReturn(0);
1181 }
1182 
1183 
1184 #undef __FUNCT__
1185 #define __FUNCT__ "DMMoabGetFieldDofsLocal"
1186 PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
1187 {
1188   PetscInt i,offset;
1189   PetscErrorCode  ierr;
1190   DM_Moab        *dmmoab;
1191 
1192   PetscFunctionBegin;
1193   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1194   PetscValidPointer(points,2);
1195   dmmoab = (DM_Moab*)(dm)->data;
1196 
1197   if (!dof) {
1198     ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr);
1199   }
1200 
1201   if (dmmoab->bs > 1) {
1202     for (i=0; i<npoints; ++i)
1203       dof[i] = (points[i]-1)*dmmoab->bs+field;
1204   }
1205   else {
1206     offset = field*dmmoab->n; /* assume all fields have equal distribution */
1207     for (i=0; i<npoints; ++i)
1208       dof[i] = offset+points[i]-1;
1209   }
1210   PetscFunctionReturn(0);
1211 }
1212 
1213 
1214 #undef __FUNCT__
1215 #define __FUNCT__ "DMMoabGetDofs"
1216 PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
1217 {
1218   PetscInt i,f,gid;
1219   PetscSection section;
1220   PetscErrorCode  ierr;
1221   moab::ErrorCode merr;
1222   DM_Moab        *dmmoab;
1223 
1224   PetscFunctionBegin;
1225   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1226   PetscValidPointer(points,2);
1227   dmmoab = (DM_Moab*)(dm)->data;
1228 
1229   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1230   if (!dof) {
1231     ierr = PetscMalloc(sizeof(PetscInt)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr);
1232   }
1233 
1234   /* first get the local indices */
1235   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr);
1236 
1237   for (i=0; i<npoints; ++i) {
1238     gid=dof[i];
1239     for (f=0; f<dmmoab->nfields; ++f) {
1240       ierr = PetscSectionGetFieldDof(section, gid, f, &dof[i*dmmoab->nfields+f]);CHKERRQ(ierr);
1241     }
1242   }
1243   PetscFunctionReturn(0);
1244 }
1245 
1246 
1247 #undef __FUNCT__
1248 #define __FUNCT__ "DMMoabGetDofsLocal"
1249 PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
1250 {
1251   PetscInt        i,f,offset;
1252   PetscErrorCode  ierr;
1253   DM_Moab        *dmmoab;
1254 
1255   PetscFunctionBegin;
1256   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1257   PetscValidPointer(points,2);
1258   dmmoab = (DM_Moab*)(dm)->data;
1259 
1260   if (!dof) {
1261     ierr = PetscMalloc(sizeof(PetscInt)*dmmoab->nfields*npoints, &dof);CHKERRQ(ierr);
1262   }
1263 
1264   if (dmmoab->bs > 1) {
1265     for (f=0; f<dmmoab->nfields; ++f)
1266       for (i=0; i<npoints; ++i)
1267         dof[i*dmmoab->nfields+f] = (points[i]-1)*dmmoab->bs+f;
1268   }
1269   else {
1270     for (f=0; f<dmmoab->nfields; ++f) {
1271       offset = f*dmmoab->n;     /* assume all fields have equal distribution - say all vertex based */
1272       for (i=0; i<npoints; ++i)
1273         dof[i*dmmoab->nfields+f] = offset+points[i]-1;
1274     }
1275   }
1276   PetscFunctionReturn(0);
1277 }
1278 
1279 
1280 #undef __FUNCT__
1281 #define __FUNCT__ "DMMoabGetDofsBlocked"
1282 PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
1283 {
1284   PetscInt i,gid,dofindx;
1285   PetscSection section;
1286   PetscErrorCode  ierr;
1287   moab::ErrorCode merr;
1288   DM_Moab        *dmmoab;
1289 
1290   PetscFunctionBegin;
1291   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1292   PetscValidPointer(points,2);
1293   dmmoab = (DM_Moab*)(dm)->data;
1294 
1295   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1296   if (!dof) {
1297     ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr);
1298   }
1299 
1300   /* first get the local indices */
1301   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr);
1302 
1303   for (i=0; i<npoints; ++i) {
1304     gid=dof[i];
1305     ierr = PetscSectionGetFieldDof(section, gid, 0, &dofindx);CHKERRQ(ierr);
1306     if (dmmoab->bs > 1)  dof[i]=dofindx/dmmoab->bs;
1307     else dof[i]=dofindx;
1308   }
1309   PetscFunctionReturn(0);
1310 }
1311 
1312 
1313 #undef __FUNCT__
1314 #define __FUNCT__ "DMMoabGetDofsBlockedLocal"
1315 PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
1316 {
1317   PetscInt        i;
1318   PetscErrorCode  ierr;
1319 
1320   PetscFunctionBegin;
1321   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1322   PetscValidPointer(points,2);
1323 
1324   if (!dof) {
1325     ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr);
1326   }
1327 
1328   for (i=0; i<npoints; ++i)
1329     dof[i] = points[i]-1;
1330   PetscFunctionReturn(0);
1331 }
1332 
1333 
1334 #undef __FUNCT__
1335 #define __FUNCT__ "DMMoabGetVertexDofsBlocked"
1336 PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm,PetscInt** dof)
1337 {
1338   PetscInt        i,gid;
1339   DM_Moab        *dmmoab;
1340   PetscSection section;
1341   PetscErrorCode  ierr;
1342   moab::ErrorCode merr;
1343 
1344   PetscFunctionBegin;
1345   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1346   dmmoab = (DM_Moab*)(dm)->data;
1347 
1348   *dof = dmmoab->gsindices;
1349 
1350   if (false) {
1351     if (!dof) {
1352       ierr = PetscMalloc(sizeof(PetscInt)*(dmmoab->nloc+dmmoab->nghost), dof);CHKERRQ(ierr);
1353     }
1354 
1355     /* first get the local indices */
1356     merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vlocal,*dof);MBERRNM(merr);
1357 
1358     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1359 
1360     /* Compute function over the locally owned part of the grid */
1361     for(i=0; i<dmmoab->nloc+dmmoab->nghost; i++) {
1362       gid=(*dof)[i];
1363       ierr = PetscSectionGetFieldDof(section, gid, 0, &(*dof)[i]);CHKERRQ(ierr);
1364     }
1365   }
1366   PetscFunctionReturn(0);
1367 }
1368 
1369 
1370 #undef __FUNCT__
1371 #define __FUNCT__ "DMMoabGetVertexDofsBlockedLocal"
1372 PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm,PetscInt** dof)
1373 {
1374   PetscInt        i;
1375   DM_Moab        *dmmoab;
1376   PetscErrorCode  ierr;
1377 
1378   PetscFunctionBegin;
1379   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1380   PetscValidPointer(dof,2);
1381   dmmoab = (DM_Moab*)(dm)->data;
1382 
1383   if (!(*dof)) {
1384     ierr = PetscMalloc(sizeof(PetscInt)*(dmmoab->nloc+dmmoab->nghost), dof);CHKERRQ(ierr);
1385   }
1386 
1387   i=0;
1388   /* Compute function over the locally owned part of the grid */
1389   for(moab::Range::iterator iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++,i++) {
1390     (*dof)[i] = (*iter)-1;
1391   }
1392   for(moab::Range::iterator iter = dmmoab->vghost->begin(); iter != dmmoab->vghost->end(); iter++,i++) {
1393     (*dof)[i] = (*iter)-1;
1394   }
1395   PetscFunctionReturn(0);
1396 }
1397 
1398 
1399 #undef __FUNCT__
1400 #define __FUNCT__ "DMMoab_GetWriteOptions_Private"
1401 PetscErrorCode DMMoab_GetWriteOptions_Private(PetscInt fsetid, PetscInt numproc, PetscInt dim, MoabWriteMode mode, PetscInt dbglevel, const char* extra_opts, const char** write_opts)
1402 {
1403   std::ostringstream str;
1404 
1405   PetscFunctionBegin;
1406 
1407   // do parallel read unless only one processor
1408   if (numproc > 1) {
1409     str << "PARALLEL=" << mode << ";";
1410     if (fsetid>=0) str << "PARALLEL_COMM=" << fsetid << ";";
1411   }
1412 
1413   if (dbglevel)
1414     str << "CPUTIME;DEBUG_IO=" << dbglevel << ";";
1415 
1416   if (extra_opts)
1417     str << extra_opts ;
1418 
1419   *write_opts = str.str().c_str();
1420   PetscFunctionReturn(0);
1421 }
1422 
1423 
1424 #undef __FUNCT__
1425 #define __FUNCT__ "DMMoabOutput"
1426 PetscErrorCode DMMoabOutput(DM dm,const char* filename,const char* usrwriteopts)
1427 {
1428   DM_Moab        *dmmoab;
1429   PetscInt       dbglevel=0;
1430   const char *writeopts;
1431 
1432   PetscErrorCode ierr;
1433   moab::ErrorCode merr;
1434 
1435   PetscFunctionBegin;
1436   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1437   dmmoab = (DM_Moab*)(dm)->data;
1438 
1439   PetscBarrier((PetscObject)dm);
1440 
1441   /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */
1442   ierr  = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab");
1443   ierr  = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr);
1444   ierr  = PetscOptionsEnd();
1445 
1446   /* add mesh loading options specific to the DM */
1447   ierr = DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, MOAB_PARWOPTS_WRITE_PART, dbglevel, usrwriteopts, &writeopts);CHKERRQ(ierr);
1448   PetscInfo2(dm, "Writing file %s with options: %s\n",filename,writeopts);
1449 
1450   /* output file, using parallel write */
1451   merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1);MBERRVM(dmmoab->mbiface,"Writing output of DMMoab failed.",merr);
1452   PetscFunctionReturn(0);
1453 }
1454 
1455