xref: /petsc/src/dm/impls/moab/dmmoab.cxx (revision 07da4cae19445fca7142c05c299201ea5e718035)
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 
15   PetscFunctionBegin;
16   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
17   if (dmmoab->icreatedinstance) {
18     delete dmmoab->mbiface;
19   }
20   dmmoab->mbiface = NULL;
21   dmmoab->pcomm = NULL;
22   delete dmmoab->vlocal;
23   delete dmmoab->vowned;
24   delete dmmoab->vghost;
25   delete dmmoab->elocal;
26   delete dmmoab->eghost;
27   delete dmmoab->bndyvtx;
28   delete dmmoab->bndyfaces;
29   delete dmmoab->bndyelems;
30 
31   ierr = PetscFree(dmmoab->gsindices);CHKERRQ(ierr);
32   ierr = VecScatterDestroy(&dmmoab->ltog_sendrecv);CHKERRQ(ierr);
33   ierr = ISLocalToGlobalMappingDestroy(&dmmoab->ltog_map);CHKERRQ(ierr);
34   ierr = PetscFree(dm->data);CHKERRQ(ierr);
35   PetscFunctionReturn(0);
36 }
37 
38 #undef __FUNCT__
39 #define __FUNCT__ "DMSetUp_Moab"
40 PetscErrorCode DMSetUp_Moab(DM dm)
41 {
42   PetscErrorCode          ierr;
43   moab::ErrorCode         merr;
44   Vec                     local, global;
45   IS                      from;
46   moab::Range::iterator   iter;
47   PetscInt                i,j,bs,gsiz,lsiz;
48   DM_Moab                *dmmoab = (DM_Moab*)dm->data;
49   PetscInt                totsize,local_min,local_max,global_min;
50   PetscSection            section;
51 
52   moab::Range adj;
53 
54   PetscFunctionBegin;
55   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
56   /* Get the local and shared vertices and cache it */
57   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.");
58 
59   /* Get the entities recursively in the current part of the mesh, if user did not set the local vertices explicitly */
60   if (dmmoab->vlocal->empty()) {
61     merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,*dmmoab->vlocal,true);MBERRNM(merr);
62 
63     /* filter based on parallel status */
64     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,dmmoab->vowned);MBERRNM(merr);
65     *dmmoab->vghost = moab::subtract(*dmmoab->vlocal, *dmmoab->vowned);
66 
67     dmmoab->nloc = dmmoab->vowned->size();
68     dmmoab->nghost = dmmoab->vghost->size();
69     ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
70 
71 #if 0
72     if(dmmoab->pcomm->rank() || dmmoab->pcomm->size()==1) {
73       PetscPrintf(PETSC_COMM_SELF, "Vertices: global: %D, local: %D", dmmoab->n, dmmoab->nloc+dmmoab->nghost);
74       dmmoab->vlocal->print(0);
75       PetscPrintf(PETSC_COMM_SELF, "Vertices: owned: %D", dmmoab->nloc);
76       dmmoab->vowned->print(0);
77       PetscPrintf(PETSC_COMM_SELF, "Vertices: ghost: %D", dmmoab->nghost);
78       dmmoab->vghost->print(0);
79     }
80 #endif
81   }
82 
83   /* get the information about the local elements in the mesh */
84   {
85     dmmoab->eghost->clear();
86 
87     /* first decipher the leading dimension */
88     for (i=3;i>0;i--) {
89       dmmoab->elocal->clear();
90       merr = dmmoab->mbiface->get_entities_by_dimension(dmmoab->fileset, i, *dmmoab->elocal, true);CHKERRQ(merr);
91 
92       /* store the current mesh dimension */
93       if (dmmoab->elocal->size()) {
94         dmmoab->dim=i;
95         break;
96       }
97     }
98 
99     *dmmoab->eghost = *dmmoab->elocal;
100     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
101     *dmmoab->eghost = moab::subtract(*dmmoab->eghost, *dmmoab->elocal);
102 
103     dmmoab->neleloc = dmmoab->elocal->size();
104     ierr = MPI_Allreduce(&dmmoab->neleloc, &dmmoab->nele, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
105   }
106 
107   bs = dmmoab->bs;
108   if (!dmmoab->ltog_tag) {
109     /* Get the global ID tag. The global ID tag is applied to each
110        vertex. It acts as an global identifier which MOAB uses to
111        assemble the individual pieces of the mesh */
112     merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr);
113   }
114 
115   totsize=dmmoab->vlocal->size();
116   ierr = PetscMalloc(totsize*sizeof(PetscInt), &dmmoab->gsindices);CHKERRQ(ierr);
117   {
118     /* first get the local indices */
119     merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vowned,&dmmoab->gsindices[0]);MBERRNM(merr);
120     /* next get the ghosted indices */
121     if (dmmoab->nghost) {
122       merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,*dmmoab->vghost,&dmmoab->gsindices[dmmoab->nloc]);MBERRNM(merr);
123     }
124 
125     /* find out the local and global minima of GLOBAL_ID */
126     local_min=local_max=dmmoab->gsindices[0];
127     for (i=0; i<totsize; ++i) {
128 //      if (dmmoab->pcomm->rank())
129 //        PetscPrintf(PETSC_COMM_SELF, "[%D] gsindices[%D] = %D\n", dmmoab->pcomm->rank(), i, dmmoab->gsindices[i]);
130       if(local_min>dmmoab->gsindices[i]) local_min=dmmoab->gsindices[i];
131       if(local_max<dmmoab->gsindices[i]) local_max=dmmoab->gsindices[i];
132     }
133 
134     ierr = MPI_Allreduce(&local_min, &global_min, 1, MPI_INT, MPI_MIN, ((PetscObject)dm)->comm);CHKERRQ(ierr);
135     PetscInfo3(dm, "GLOBAL_ID: Local minima - %D, Local maxima - %D, Global minima - %D.\n", local_min, local_max, global_min);
136 //    PetscPrintf(PETSC_COMM_SELF, "[%D] GLOBAL_ID: Local minima - %D, Local maxima - %D, Global minima - %D.\n", dmmoab->pcomm->rank(), local_min, local_max, global_min);
137   }
138 
139   {
140     ierr = PetscSectionCreate(((PetscObject)dm)->comm, &section);CHKERRQ(ierr);
141     ierr = PetscSectionSetNumFields(section, dmmoab->nfields);CHKERRQ(ierr);
142 //    ierr = PetscSectionSetChart(section, dmmoab->gsindices[0], dmmoab->gsindices[dmmoab->nloc-1]+1);CHKERRQ(ierr);
143     ierr = PetscSectionSetChart(section, local_min, local_max+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-global_min)*dmmoab->nfields+i);CHKERRQ(ierr);
150           ierr = PetscSectionSetFieldOffset(section, locgid, i, (locgid-global_min)*dmmoab->nfields);
151         }
152         else {
153           ierr = PetscSectionSetFieldDof(section, locgid, i, dmmoab->n*i+locgid-global_min);CHKERRQ(ierr);
154           ierr = PetscSectionSetFieldOffset(section, locgid, i, i*dmmoab->n);
155           PetscPrintf(PETSC_COMM_SELF, "[%D] Local_GID = %D, FDOF = %D, OFF = %D.\n", dmmoab->pcomm->rank(), locgid, dmmoab->n*i+locgid-global_min, i*dmmoab->n );
156         }
157       }
158       ierr = PetscSectionSetDof(section, locgid, dmmoab->nfields);CHKERRQ(ierr);
159     }
160     ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
161     ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr);
162   }
163 
164   {
165     for (i=0; i<totsize; ++i) {
166       dmmoab->gsindices[i]-=global_min;   /* zero based index needed for IS */
167 //      if (dmmoab->pcomm->rank())
168 //        PetscPrintf(PETSC_COMM_SELF, "[%D] modified gsindices[%D] = %D\n", dmmoab->pcomm->rank(), i, dmmoab->gsindices[i]);
169     }
170 
171     /* Create Global to Local Vector Scatter Context */
172     ierr = DMCreateGlobalVector_Moab(dm, &global);CHKERRQ(ierr);
173     ierr = DMCreateLocalVector_Moab(dm, &local);CHKERRQ(ierr);
174 
175     /* global to local must retrieve ghost points */
176     ierr = ISCreateBlock(((PetscObject)dm)->comm,bs,totsize,&dmmoab->gsindices[0],PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
177 
178     ierr = VecGetLocalSize(global,&gsiz);CHKERRQ(ierr);
179     ierr = VecGetLocalSize(local,&lsiz);CHKERRQ(ierr);
180 
181     ierr = VecScatterCreate(local,from,global,from,&dmmoab->ltog_sendrecv);CHKERRQ(ierr);
182     ierr = ISDestroy(&from);CHKERRQ(ierr);
183     ierr = VecDestroy(&local);CHKERRQ(ierr);
184     ierr = VecDestroy(&global);CHKERRQ(ierr);
185   }
186 
187   /* skin the boundary and store nodes */
188   {
189     // get the skin vertices of those faces and mark them as fixed; we don't want to fix the vertices on a
190     // part boundary, but since we exchanged a layer of ghost faces, those vertices aren't on the skin locally
191     // ok to mark non-owned skin vertices too, I won't move those anyway
192     // use MOAB's skinner class to find the skin
193     moab::Skinner skinner(dmmoab->mbiface);
194     dmmoab->bndyvtx = new moab::Range();
195     dmmoab->bndyfaces = new moab::Range();
196     dmmoab->bndyelems = new moab::Range();
197     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
198 //    merr = skinner.find_skin(dmmoab->fileset, *dmmoab->elocal, true, *dmmoab->bndyvtx, NULL, true, true, false, false);MBERRNM(merr); // 'true' param indicates we want vertices back, not faces
199 
200     if (dmmoab->dim == 3) {
201       // get the edges from faces and then do the same if needed
202     }
203 
204     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->bndyfaces,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
205     merr = dmmoab->pcomm->filter_pstatus(*dmmoab->bndyfaces,PSTATUS_SHARED,PSTATUS_NOT);MBERRNM(merr);
206 //    merr = dmmoab->pcomm->filter_pstatus(*dmmoab->bndyfaces,PSTATUS_MULTISHARED,PSTATUS_NOT);MBERRNM(merr);
207 
208     merr = dmmoab->mbiface->get_connectivity(*dmmoab->bndyfaces, *dmmoab->bndyvtx, false);MBERRNM(ierr);
209     merr = dmmoab->mbiface->get_adjacencies(*dmmoab->bndyfaces, dmmoab->dim, false, *dmmoab->bndyelems, moab::Interface::UNION);MBERRNM(ierr);
210 
211     dmmoab->bndyfaces->print(0);
212     dmmoab->bndyvtx->print(0);
213     dmmoab->bndyelems->print(0);
214 
215 //    merr = skinner.find_geometric_skin(dmmoab->fileset, *dmmoab->bndyelems, dmmoab->dim);MBERRV(dmmoab->mbiface,merr);
216 //    merr = skinner.find_skin(dmmoab->fileset, *dmmoab->elocal, dmmoab->bndyvtx, dmmoab->bndyfaces, dmmoab->bndyelems, true, true);MBERRNM(merr);
217 
218     PetscInfo3(dm, "Found %D boundary vertices, %D boundary faces and %D boundary elements.\n", dmmoab->bndyvtx->size(), dmmoab->bndyvtx->size(), dmmoab->bndyelems->size());
219   }
220   PetscFunctionReturn(0);
221 }
222 
223 
224 #undef __FUNCT__
225 #define __FUNCT__ "DMCreate_Moab"
226 PETSC_EXTERN PetscErrorCode DMCreate_Moab(DM dm)
227 {
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
232   ierr = PetscNewLog(dm,DM_Moab,&dm->data);CHKERRQ(ierr);
233 
234   ((DM_Moab*)dm->data)->bs = 1;
235   ((DM_Moab*)dm->data)->n = 0;
236   ((DM_Moab*)dm->data)->nloc = 0;
237   ((DM_Moab*)dm->data)->nele = 0;
238   ((DM_Moab*)dm->data)->neleloc = 0;
239   ((DM_Moab*)dm->data)->nghost = 0;
240   ((DM_Moab*)dm->data)->ltog_map = PETSC_NULL;
241   ((DM_Moab*)dm->data)->ltog_sendrecv = PETSC_NULL;
242 
243   ((DM_Moab*)dm->data)->vlocal = new moab::Range();
244   ((DM_Moab*)dm->data)->vowned = new moab::Range();
245   ((DM_Moab*)dm->data)->vghost = new moab::Range();
246   ((DM_Moab*)dm->data)->elocal = new moab::Range();
247   ((DM_Moab*)dm->data)->eghost = new moab::Range();
248 
249   dm->ops->createglobalvector              = DMCreateGlobalVector_Moab;
250   dm->ops->createlocalvector               = DMCreateLocalVector_Moab;
251   dm->ops->creatematrix                    = DMCreateMatrix_Moab;
252   dm->ops->setup                           = DMSetUp_Moab;
253   dm->ops->destroy                         = DMDestroy_Moab;
254   dm->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Moab;
255   dm->ops->globaltolocalend                = DMGlobalToLocalEnd_Moab;
256   dm->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Moab;
257   dm->ops->localtoglobalend                = DMLocalToGlobalEnd_Moab;
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "DMMoabCreate"
263 /*@
264   DMMoabCreate - Creates a DMMoab object, which encapsulates a moab instance
265 
266   Collective on MPI_Comm
267 
268   Input Parameter:
269 . comm - The communicator for the DMMoab object
270 
271   Output Parameter:
272 . dmb  - The DMMoab object
273 
274   Level: beginner
275 
276 .keywords: DMMoab, create
277 @*/
278 PetscErrorCode DMMoabCreate(MPI_Comm comm, DM *dmb)
279 {
280   PetscErrorCode ierr;
281 
282   PetscFunctionBegin;
283   PetscValidPointer(dmb,2);
284   ierr = DMCreate(comm, dmb);CHKERRQ(ierr);
285   ierr = DMSetType(*dmb, DMMOAB);CHKERRQ(ierr);
286   PetscFunctionReturn(0);
287 }
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "DMMoabCreateMoab"
291 /*@
292   DMMoabCreate - Creates a DMMoab object, optionally from an instance and other data
293 
294   Collective on MPI_Comm
295 
296   Input Parameter:
297 . comm - The communicator for the DMMoab object
298 . mbiface - (ptr to) the MOAB Instance; if passed in NULL, MOAB instance is created inside PETSc, and destroyed
299          along with the DMMoab
300 . pcomm - (ptr to) a ParallelComm; if NULL, creates one internally for the whole communicator
301 . ltog_tag - A tag to use to retrieve global id for an entity; if 0, will use GLOBAL_ID_TAG_NAME/tag
302 . range - If non-NULL, contains range of entities to which DOFs will be assigned
303 
304   Output Parameter:
305 . dmb  - The DMMoab object
306 
307   Level: intermediate
308 
309 .keywords: DMMoab, create
310 @*/
311 PetscErrorCode DMMoabCreateMoab(MPI_Comm comm, moab::Interface *mbiface, moab::ParallelComm *pcomm, moab::Tag *ltog_tag, moab::Range *range, DM *dmb)
312 {
313   PetscErrorCode ierr;
314   moab::ErrorCode merr;
315   moab::EntityHandle partnset;
316   PetscInt rank, nprocs;
317   DM_Moab        *dmmoab;
318 
319   PetscFunctionBegin;
320   PetscValidPointer(dmb,6);
321   ierr = DMMoabCreate(comm, dmb);CHKERRQ(ierr);
322   dmmoab = (DM_Moab*)(*dmb)->data;
323 
324   if (!mbiface) {
325     dmmoab->mbiface = new moab::Core();
326     dmmoab->icreatedinstance = PETSC_TRUE;
327   }
328   else {
329     dmmoab->mbiface = mbiface;
330     dmmoab->icreatedinstance = PETSC_FALSE;
331   }
332 
333   /* create a fileset to store the hierarchy of entities belonging to current DM */
334   merr = dmmoab->mbiface->create_meshset(moab::MESHSET_ORDERED, dmmoab->fileset);MBERR("Creating file set failed", merr);
335 
336   if (!pcomm) {
337     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
338     ierr = MPI_Comm_size(comm, &nprocs);CHKERRQ(ierr);
339 
340     /* Create root sets for each mesh.  Then pass these
341        to the load_file functions to be populated. */
342     merr = dmmoab->mbiface->create_meshset(moab::MESHSET_SET, partnset);
343     MBERR("Creating partition set failed", merr);
344 
345     /* Create the parallel communicator object with the partition handle associated with MOAB */
346     dmmoab->pcomm = moab::ParallelComm::get_pcomm(dmmoab->mbiface, partnset, &comm);
347   }
348   else {
349     ierr = DMMoabSetParallelComm(*dmb, pcomm);CHKERRQ(ierr);
350   }
351 
352   /* do the remaining initializations for DMMoab */
353   dmmoab->bs       = 1;
354 
355   /* set global ID tag handle */
356   if (!ltog_tag) {
357     merr = dmmoab->mbiface->tag_get_handle(GLOBAL_ID_TAG_NAME, dmmoab->ltog_tag);MBERRNM(merr);
358   }
359   else {
360     ierr = DMMoabSetLocalToGlobalTag(*dmb, *ltog_tag);CHKERRQ(ierr);
361   }
362 
363   /* set the local range of entities (vertices) of interest */
364   if (range) {
365     ierr = DMMoabSetLocalVertices(*dmb, range);CHKERRQ(ierr);
366   }
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "DMMoabSetParallelComm"
372 /*@
373   DMMoabSetParallelComm - Set the ParallelComm used with this DMMoab
374 
375   Collective on MPI_Comm
376 
377   Input Parameter:
378 . dm    - The DMMoab object being set
379 . pcomm - The ParallelComm being set on the DMMoab
380 
381   Level: beginner
382 
383 .keywords: DMMoab, create
384 @*/
385 PetscErrorCode DMMoabSetParallelComm(DM dm,moab::ParallelComm *pcomm)
386 {
387   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
388 
389   PetscFunctionBegin;
390   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
391   PetscValidPointer(pcomm,2);
392   dmmoab->pcomm = pcomm;
393   dmmoab->mbiface = pcomm->get_moab();
394   dmmoab->icreatedinstance = PETSC_FALSE;
395   PetscFunctionReturn(0);
396 }
397 
398 
399 #undef __FUNCT__
400 #define __FUNCT__ "DMMoabGetParallelComm"
401 /*@
402   DMMoabGetParallelComm - Get the ParallelComm used with this DMMoab
403 
404   Collective on MPI_Comm
405 
406   Input Parameter:
407 . dm    - The DMMoab object being set
408 
409   Output Parameter:
410 . pcomm - The ParallelComm for the DMMoab
411 
412   Level: beginner
413 
414 .keywords: DMMoab, create
415 @*/
416 PetscErrorCode DMMoabGetParallelComm(DM dm,moab::ParallelComm **pcomm)
417 {
418   PetscFunctionBegin;
419   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
420   *pcomm = ((DM_Moab*)(dm)->data)->pcomm;
421   PetscFunctionReturn(0);
422 }
423 
424 
425 #undef __FUNCT__
426 #define __FUNCT__ "DMMoabSetInterface"
427 /*@
428   DMMoabSetInterface - Set the MOAB instance used with this DMMoab
429 
430   Collective on MPI_Comm
431 
432   Input Parameter:
433 . dm      - The DMMoab object being set
434 . mbiface - The MOAB instance being set on this DMMoab
435 
436   Level: beginner
437 
438 .keywords: DMMoab, create
439 @*/
440 PetscErrorCode DMMoabSetInterface(DM dm,moab::Interface *mbiface)
441 {
442   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
443 
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
446   PetscValidPointer(mbiface,2);
447   dmmoab->pcomm = NULL;
448   dmmoab->mbiface = mbiface;
449   dmmoab->icreatedinstance = PETSC_FALSE;
450   PetscFunctionReturn(0);
451 }
452 
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "DMMoabGetInterface"
456 /*@
457   DMMoabGetInterface - Get the MOAB instance used with this DMMoab
458 
459   Collective on MPI_Comm
460 
461   Input Parameter:
462 . dm      - The DMMoab object being set
463 
464   Output Parameter:
465 . mbiface - The MOAB instance set on this DMMoab
466 
467   Level: beginner
468 
469 .keywords: DMMoab, create
470 @*/
471 PetscErrorCode DMMoabGetInterface(DM dm,moab::Interface **mbiface)
472 {
473   PetscErrorCode   ierr;
474   static PetscBool cite = PETSC_FALSE;
475 
476   PetscFunctionBegin;
477   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
478   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);
479   *mbiface = ((DM_Moab*)dm->data)->mbiface;
480   PetscFunctionReturn(0);
481 }
482 
483 
484 #undef __FUNCT__
485 #define __FUNCT__ "DMMoabSetLocalVertices"
486 /*@
487   DMMoabSetLocalVertices - Set the entities having DOFs on this DMMoab
488 
489   Collective on MPI_Comm
490 
491   Input Parameter:
492 . dm    - The DMMoab object being set
493 . range - The entities treated by this DMMoab
494 
495   Level: beginner
496 
497 .keywords: DMMoab, create
498 @*/
499 PetscErrorCode DMMoabSetLocalVertices(DM dm,moab::Range *range)
500 {
501   moab::ErrorCode merr;
502   PetscErrorCode  ierr;
503   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
504 
505   PetscFunctionBegin;
506   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
507   dmmoab->vlocal->clear();
508   dmmoab->vowned->clear();
509   dmmoab->vlocal->insert(range->begin(), range->end());
510   *dmmoab->vowned = *dmmoab->vlocal;
511   merr = dmmoab->pcomm->filter_pstatus(*dmmoab->vowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
512   *dmmoab->vghost = moab::subtract(*range, *dmmoab->vowned);
513   dmmoab->nloc=dmmoab->vowned->size();
514   dmmoab->nghost=dmmoab->vghost->size();
515   ierr = MPI_Allreduce(&dmmoab->nloc, &dmmoab->n, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
516   PetscFunctionReturn(0);
517 }
518 
519 
520 #undef __FUNCT__
521 #define __FUNCT__ "DMMoabGetLocalVertices"
522 /*@
523   DMMoabGetLocalVertices - Get the entities having DOFs on this DMMoab
524 
525   Collective on MPI_Comm
526 
527   Input Parameter:
528 . dm    - The DMMoab object being set
529 
530   Output Parameter:
531 . owned - The owned vertex entities in this DMMoab
532 . ghost - The ghosted entities (non-owned) stored locally in this partition
533 
534   Level: beginner
535 
536 .keywords: DMMoab, create
537 @*/
538 PetscErrorCode DMMoabGetLocalVertices(DM dm,moab::Range *owned,moab::Range *ghost)
539 {
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
542   if (owned) *owned = *((DM_Moab*)dm->data)->vowned;
543   if (ghost) *ghost = *((DM_Moab*)dm->data)->vghost;
544   PetscFunctionReturn(0);
545 }
546 
547 #undef __FUNCT__
548 #define __FUNCT__ "DMMoabGetLocalElements"
549 /*@
550   DMMoabGetLocalElements - Get the higher-dimensional entities that are locally owned
551 
552   Collective on MPI_Comm
553 
554   Input Parameter:
555 . dm    - The DMMoab object being set
556 
557   Output Parameter:
558 . range - The entities owned locally
559 
560   Level: beginner
561 
562 .keywords: DMMoab, create
563 @*/
564 PetscErrorCode DMMoabGetLocalElements(DM dm,moab::Range *range)
565 {
566   PetscFunctionBegin;
567   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
568   if (range) *range = *((DM_Moab*)dm->data)->elocal;
569   PetscFunctionReturn(0);
570 }
571 
572 
573 #undef __FUNCT__
574 #define __FUNCT__ "DMMoabSetLocalElements"
575 /*@
576   DMMoabSetLocalElements - Set the entities having DOFs on this DMMoab
577 
578   Collective on MPI_Comm
579 
580   Input Parameter:
581 . dm    - The DMMoab object being set
582 . range - The entities treated by this DMMoab
583 
584   Level: beginner
585 
586 .keywords: DMMoab, create
587 @*/
588 PetscErrorCode DMMoabSetLocalElements(DM dm,moab::Range *range)
589 {
590   moab::ErrorCode merr;
591   PetscErrorCode  ierr;
592   DM_Moab        *dmmoab = (DM_Moab*)(dm)->data;
593 
594   PetscFunctionBegin;
595   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
596   dmmoab->elocal->clear();
597   dmmoab->eghost->clear();
598   dmmoab->elocal->insert(range->begin(), range->end());
599   PetscInfo2(dm, "Range size = %D; elocal size = %D.\n", range->size(), dmmoab->elocal->size());
600   merr = dmmoab->pcomm->filter_pstatus(*dmmoab->elocal,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
601   *dmmoab->eghost = moab::subtract(*range, *dmmoab->elocal);
602   dmmoab->neleloc=dmmoab->elocal->size();
603   ierr = MPI_Allreduce(&dmmoab->nele, &dmmoab->neleloc, 1, MPI_INTEGER, MPI_SUM, ((PetscObject)dm)->comm);CHKERRQ(ierr);
604   PetscInfo2(dm, "Created %D local and %D glocal elements.\n", dmmoab->neleloc, dmmoab->nele);
605   PetscFunctionReturn(0);
606 }
607 
608 
609 #undef __FUNCT__
610 #define __FUNCT__ "DMMoabSetLocalToGlobalTag"
611 /*@
612   DMMoabSetLocalToGlobalTag - Set the tag used for local to global numbering
613 
614   Collective on MPI_Comm
615 
616   Input Parameter:
617 . dm      - The DMMoab object being set
618 . ltogtag - The MOAB tag used for local to global ids
619 
620   Level: beginner
621 
622 .keywords: DMMoab, create
623 @*/
624 PetscErrorCode DMMoabSetLocalToGlobalTag(DM dm,moab::Tag ltogtag)
625 {
626   PetscFunctionBegin;
627   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
628   ((DM_Moab*)dm->data)->ltog_tag = ltogtag;
629   PetscFunctionReturn(0);
630 }
631 
632 
633 #undef __FUNCT__
634 #define __FUNCT__ "DMMoabGetLocalToGlobalTag"
635 /*@
636   DMMoabGetLocalToGlobalTag - Get the tag used for local to global numbering
637 
638   Collective on MPI_Comm
639 
640   Input Parameter:
641 . dm      - The DMMoab object being set
642 
643   Output Parameter:
644 . ltogtag - The MOAB tag used for local to global ids
645 
646   Level: beginner
647 
648 .keywords: DMMoab, create
649 @*/
650 PetscErrorCode DMMoabGetLocalToGlobalTag(DM dm,moab::Tag *ltog_tag)
651 {
652   PetscFunctionBegin;
653   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
654   *ltog_tag = ((DM_Moab*)dm->data)->ltog_tag;
655   PetscFunctionReturn(0);
656 }
657 
658 
659 #undef __FUNCT__
660 #define __FUNCT__ "DMMoabSetBlockSize"
661 /*@
662   DMMoabSetBlockSize - Set the block size used with this DMMoab
663 
664   Collective on MPI_Comm
665 
666   Input Parameter:
667 . dm - The DMMoab object being set
668 . bs - The block size used with this DMMoab
669 
670   Level: beginner
671 
672 .keywords: DMMoab, create
673 @*/
674 PetscErrorCode DMMoabSetBlockSize(DM dm,PetscInt bs)
675 {
676   PetscFunctionBegin;
677   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
678   ((DM_Moab*)dm->data)->bs = bs;
679   PetscFunctionReturn(0);
680 }
681 
682 
683 #undef __FUNCT__
684 #define __FUNCT__ "DMMoabGetBlockSize"
685 /*@
686   DMMoabGetBlockSize - Get the block size used with this DMMoab
687 
688   Collective on MPI_Comm
689 
690   Input Parameter:
691 . dm - The DMMoab object being set
692 
693   Output Parameter:
694 . bs - The block size used with this DMMoab
695 
696   Level: beginner
697 
698 .keywords: DMMoab, create
699 @*/
700 PetscErrorCode DMMoabGetBlockSize(DM dm,PetscInt *bs)
701 {
702   PetscFunctionBegin;
703   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
704   *bs = ((DM_Moab*)dm->data)->bs;
705   PetscFunctionReturn(0);
706 }
707 
708 
709 #undef __FUNCT__
710 #define __FUNCT__ "DMMoabSetFieldVector"
711 PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec)
712 {
713   DM_Moab        *dmmoab;
714   moab::Tag     vtag,ntag;
715   PetscScalar   *varray;
716   moab::ErrorCode merr;
717   PetscErrorCode  ierr;
718   PetscSection section;
719   PetscInt doff;
720   std::string tag_name;
721   moab::Range::iterator iter;
722 
723   PetscFunctionBegin;
724   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
725   dmmoab = (DM_Moab*)(dm)->data;
726 
727   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
728 
729   /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
730   merr = dmmoab->mbiface->tag_get_handle(dmmoab->fields[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
731                                           moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT);MBERRNM(merr);
732 
733   ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr);
734 
735   merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
736   if (!tag_name.length() && merr !=moab::MB_SUCCESS) {
737     ierr = DMMoabVecGetArray(dm,fvec,&varray);CHKERRQ(ierr);
738     for(iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++) {
739       moab::EntityHandle vtx = (*iter);
740 
741       /* get field dof index */
742       ierr = PetscSectionGetFieldOffset(section, vtx, ifield, &doff);
743 
744       /* use the entity handle and the Dof index to set the right value */
745       merr = dmmoab->mbiface->tag_set_data(ntag, &vtx, 1, (const void*)&varray[doff]);MBERRNM(merr);
746     }
747     ierr = DMMoabVecRestoreArray(dm,fvec,&varray);CHKERRQ(ierr);
748   }
749   else {
750     ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&varray);CHKERRQ(ierr);
751     /* we are using a MOAB Vec - directly copy the tag data to new one */
752     merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr);
753     merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)varray);MBERRNM(merr);
754     /* make sure the parallel exchange for ghosts are done appropriately */
755     merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr);
756     ierr = PetscFree(varray);CHKERRQ(ierr);
757   }
758   PetscFunctionReturn(0);
759 }
760 
761 
762 #undef __FUNCT__
763 #define __FUNCT__ "DMMoabGetVertexCoordinates"
764 PetscErrorCode DMMoabGetVertexCoordinates(DM dm,PetscInt nconn,const moab::EntityHandle *conn,PetscScalar *vpos)
765 {
766   DM_Moab         *dmmoab;
767   PetscErrorCode  ierr;
768   moab::ErrorCode merr;
769 
770   PetscFunctionBegin;
771   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
772   PetscValidPointer(conn,3);
773   dmmoab = (DM_Moab*)(dm)->data;
774 
775   if (!vpos) {
776     ierr = PetscMalloc(sizeof(PetscScalar)*nconn*3, &vpos);CHKERRQ(ierr);
777   }
778 
779   /* Get connectivity information in MOAB canonical ordering */
780   merr = dmmoab->mbiface->get_coords(conn, nconn, vpos);MBERRNM(merr);
781   PetscFunctionReturn(0);
782 }
783 
784 
785 #undef __FUNCT__
786 #define __FUNCT__ "DMMoabGetElementConnectivity"
787 PetscErrorCode DMMoabGetElementConnectivity(DM dm,moab::EntityHandle ehandle,PetscInt* nconn,const moab::EntityHandle **conn)
788 {
789   DM_Moab        *dmmoab;
790   const moab::EntityHandle *connect;
791   moab::ErrorCode merr;
792   PetscInt nnodes;
793 
794   PetscFunctionBegin;
795   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
796   PetscValidPointer(conn,4);
797   dmmoab = (DM_Moab*)(dm)->data;
798 
799   /* Get connectivity information in MOAB canonical ordering */
800   merr = dmmoab->mbiface->get_connectivity(ehandle, connect, nnodes);MBERRNM(merr);
801   if (conn) *conn=connect;
802   if (nconn) *nconn=nnodes;
803   PetscFunctionReturn(0);
804 }
805 
806 
807 #undef __FUNCT__
808 #define __FUNCT__ "DMMoabIsEntityOnBoundary"
809 PetscErrorCode DMMoabIsEntityOnBoundary(DM dm,const moab::EntityHandle ent,PetscBool* ent_on_boundary)
810 {
811   moab::EntityType etype;
812   DM_Moab         *dmmoab;
813   PetscInt         edim;
814 
815   PetscFunctionBegin;
816   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
817   PetscValidPointer(ent_on_boundary,3);
818   dmmoab = (DM_Moab*)(dm)->data;
819 
820   /* get the entity type and handle accordingly */
821   etype=dmmoab->mbiface->type_from_handle(ent);
822   if(etype >= moab::MBPOLYHEDRON) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Entity type on the boundary skin is invalid. EntityType = %D\n",etype);
823 
824   /* get the entity dimension */
825   edim=dmmoab->mbiface->dimension_from_handle(ent);
826 
827   *ent_on_boundary=PETSC_FALSE;
828   if(etype == moab::MBVERTEX && edim == 0) {
829     moab::Range::const_iterator giter = dmmoab->bndyvtx->find(ent);
830     if (giter != dmmoab->bndyvtx->end()) *ent_on_boundary=PETSC_TRUE;
831   }
832   else {
833     if (edim == dmmoab->dim) {  /* check the higher-dimensional elements first */
834       moab::Range::const_iterator geiter = dmmoab->bndyelems->find(ent);
835       if (geiter != dmmoab->bndyelems->end()) *ent_on_boundary=PETSC_TRUE;
836     }
837     else {                      /* next check the lower-dimensional faces */
838       moab::Range::const_iterator gfiter = dmmoab->bndyfaces->find(ent);
839       if (gfiter != dmmoab->bndyfaces->end()) *ent_on_boundary=PETSC_TRUE;
840     }
841   }
842   PetscFunctionReturn(0);
843 }
844 
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "DMMoabCheckBoundaryVertices"
848 PetscErrorCode DMMoabCheckBoundaryVertices(DM dm,PetscInt nconn,const moab::EntityHandle *cnt,PetscBool* isbdvtx)
849 {
850   DM_Moab        *dmmoab;
851   PetscInt       i;
852 
853   PetscFunctionBegin;
854   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
855   PetscValidPointer(cnt,3);
856   PetscValidPointer(isbdvtx,4);
857   dmmoab = (DM_Moab*)(dm)->data;
858 
859   for (i=0; i < nconn; ++i) {
860     moab::Range::const_iterator giter = dmmoab->bndyvtx->find(cnt[i]);
861     if (giter != dmmoab->bndyvtx->end()) isbdvtx[i] = PETSC_TRUE;
862     else isbdvtx[i] = PETSC_FALSE;
863   }
864   PetscFunctionReturn(0);
865 }
866 
867 
868 #undef __FUNCT__
869 #define __FUNCT__ "DMMoabGetBoundaryEntities"
870 PetscErrorCode DMMoabGetBoundaryEntities(DM dm,moab::Range *bdvtx,moab::Range* bdfaces,moab::Range* bdelems)
871 {
872   DM_Moab        *dmmoab;
873 
874   PetscFunctionBegin;
875   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
876   dmmoab = (DM_Moab*)(dm)->data;
877 
878   if (bdvtx)  *bdvtx = *dmmoab->bndyvtx;
879   if (bdfaces)  *bdfaces = *dmmoab->bndyfaces;
880   if (bdelems)  *bdfaces = *dmmoab->bndyelems;
881   PetscFunctionReturn(0);
882 }
883 
884 
885 #undef __FUNCT__
886 #define __FUNCT__ "DMMoabSetFields"
887 PetscErrorCode DMMoabSetFields(DM dm,PetscInt nfields,const char** fields)
888 {
889   DM_Moab        *dmmoab;
890 
891   PetscFunctionBegin;
892   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
893   dmmoab = (DM_Moab*)(dm)->data;
894 
895   dmmoab->fields = fields;
896   dmmoab->nfields = nfields;
897   PetscFunctionReturn(0);
898 }
899 
900 
901 #undef __FUNCT__
902 #define __FUNCT__ "DMMoabGetFieldDof"
903 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof)
904 {
905   PetscSection section;
906   PetscInt gid;
907   PetscErrorCode ierr;
908   moab::ErrorCode merr;
909   DM_Moab        *dmmoab;
910 
911   PetscFunctionBegin;
912   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
913   dmmoab = (DM_Moab*)(dm)->data;
914 
915   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
916 
917   /* first get the global ID for the point */
918   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&point,1,&gid);MBERRNM(merr);
919 
920   /* get the dof value for the field */
921   ierr = PetscSectionGetFieldDof(section, gid, field, dof);CHKERRQ(ierr);
922   PetscFunctionReturn(0);
923 }
924 
925 
926 #undef __FUNCT__
927 #define __FUNCT__ "DMMoabGetFieldDofs"
928 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
929 {
930   PetscInt i,gid;
931   PetscSection section;
932   PetscErrorCode  ierr;
933   moab::ErrorCode merr;
934   DM_Moab        *dmmoab;
935 
936   PetscFunctionBegin;
937   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
938   dmmoab = (DM_Moab*)(dm)->data;
939 
940   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
941   if (!dof) {
942     ierr = PetscMalloc(sizeof(PetscScalar)*npoints, &dof);CHKERRQ(ierr);
943   }
944 
945   /* first get the local indices */
946   merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,points,npoints,dof);MBERRNM(merr);
947 
948   for (i=0; i<npoints; ++i) {
949     gid=dof[i];
950     ierr = PetscSectionGetFieldDof(section, gid, field, &dof[i]);CHKERRQ(ierr);
951   }
952   PetscFunctionReturn(0);
953 }
954 
955 
956 #undef __FUNCT__
957 #define __FUNCT__ "DMMoab_GetWriteOptions_Private"
958 PetscErrorCode DMMoab_GetWriteOptions_Private(PetscInt fsetid, PetscInt numproc, PetscInt dim, MoabWriteMode mode, PetscInt dbglevel, const char* extra_opts, const char** write_opts)
959 {
960   std::ostringstream str;
961 
962   PetscFunctionBegin;
963 
964   // do parallel read unless only one processor
965   if (numproc > 1) {
966     str << "PARALLEL=" << mode << ";";
967     if (fsetid>=0) str << "PARALLEL_COMM=" << fsetid << ";";
968   }
969 
970   if (dbglevel)
971     str << "CPUTIME;DEBUG_IO=" << dbglevel << ";";
972 
973   if (extra_opts)
974     str << extra_opts ;
975 
976   *write_opts = str.str().c_str();
977   PetscFunctionReturn(0);
978 }
979 
980 
981 #undef __FUNCT__
982 #define __FUNCT__ "DMMoabOutput"
983 PetscErrorCode DMMoabOutput(DM dm,const char* filename,const char* usrwriteopts)
984 {
985   DM_Moab        *dmmoab;
986   PetscInt       dbglevel=0;
987   const char *writeopts;
988 
989   PetscErrorCode ierr;
990   moab::ErrorCode merr;
991 
992   PetscFunctionBegin;
993   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
994   dmmoab = (DM_Moab*)(dm)->data;
995 
996   PetscBarrier((PetscObject)dm);
997 
998   /* TODO: Use command-line options to control by_rank, verbosity, MoabReadMode and extra options */
999   ierr  = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for reading/writing MOAB based meshes from file", "DMMoab");
1000   ierr  = PetscOptionsInt("-dmmb_rw_dbg", "The verbosity level for reading and writing MOAB meshes", "dmmbutil.cxx", dbglevel, &dbglevel, NULL);CHKERRQ(ierr);
1001   ierr  = PetscOptionsEnd();
1002 
1003   /* add mesh loading options specific to the DM */
1004   ierr = DMMoab_GetWriteOptions_Private(dmmoab->pcomm->get_id(), dmmoab->pcomm->size(), dmmoab->dim, MOAB_PARWOPTS_WRITE_PART, dbglevel, usrwriteopts, &writeopts);CHKERRQ(ierr);
1005   PetscInfo2(dm, "Writing file %s with options: %s\n",filename,writeopts);
1006 
1007   /* output file, using parallel write */
1008   merr = dmmoab->mbiface->write_file(filename, NULL, writeopts, &dmmoab->fileset, 1);MBERRVM(dmmoab->mbiface,"Writing output of DMMoab failed.",merr);
1009   PetscFunctionReturn(0);
1010 }
1011 
1012