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