xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision a86ed394d0a6627fd6ef9fd16cfa3d0db6a8bc8a)
1 #include <petsc/private/dmmbimpl.h>
2 #include <petscdmmoab.h>
3 #include <MBTagConventions.hpp>
4 #include <moab/NestedRefine.hpp>
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "DMMoabGenerateHierarchy"
8 /*@
9   DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy
10   by succesively refining a coarse mesh, already defined in the DM object
11   provided by the user.
12 
13   Collective on MPI_Comm
14 
15   Input Parameter:
16 + dmb  - The DMMoab object
17 
18   Output Parameter:
19 + nlevels   - The number of levels of refinement needed to generate the hierarchy
20 . ldegrees  - The degree of refinement at each level in the hierarchy
21 
22   Level: beginner
23 
24 .keywords: DMMoab, create, refinement
25 @*/
26 PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees)
27 {
28   DM_Moab        *dmmoab;
29   PetscErrorCode  ierr;
30   moab::ErrorCode merr;
31   PetscInt *pdegrees, i;
32   std::vector<moab::EntityHandle> hsets;
33 
34   PetscFunctionBegin;
35   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36   dmmoab = (DM_Moab*)(dm)->data;
37 
38   if (!ldegrees) {
39     ierr = PetscMalloc1(nlevels, &pdegrees);CHKERRQ(ierr);
40     for (i = 0; i < nlevels; i++) pdegrees[i] = 2; /* default = Degree 2 refinement */
41   }
42   else pdegrees = ldegrees;
43 
44   /* initialize set level refinement data for hierarchy */
45   dmmoab->nhlevels = nlevels;
46 
47   /* Instantiate the nested refinement class */
48 #ifdef MOAB_HAVE_MPI
49   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset);
50 #else
51   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), NULL, dmmoab->fileset);
52 #endif
53 
54   ierr = PetscMalloc1(nlevels + 1, &dmmoab->hsets);CHKERRQ(ierr);
55 
56   /* generate the mesh hierarchy */
57   merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets); MBERRNM(merr);
58 
59 #ifdef MOAB_HAVE_MPI
60   if (dmmoab->pcomm->size() > 1) {
61     merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings); MBERRNM(merr);
62   }
63 #endif
64 
65   /* copy the mesh sets for nested refinement hierarchy */
66   for (i = 0; i <= nlevels; i++)
67     dmmoab->hsets[i] = hsets[i];
68 
69   hsets.clear();
70   if (!ldegrees) {
71     ierr = PetscFree(pdegrees);CHKERRQ(ierr);
72   }
73   PetscFunctionReturn(0);
74 }
75 
76 #undef __FUNCT__
77 #define __FUNCT__ "DMRefineHierarchy_Moab"
78 /*@
79   DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy
80   by succesively refining a coarse mesh.
81 
82   Collective on MPI_Comm
83 
84   Input Parameter:
85 + dm  - The DMMoab object
86 
87   Output Parameter:
88 + nlevels   - The number of levels of refinement needed to generate the hierarchy
89 . dmf  - The DM objects after successive refinement of the hierarchy
90 
91   Level: beginner
92 
93 .keywords: DMMoab, generate, refinement
94 @*/
95 PETSC_EXTERN PetscErrorCode  DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[])
96 {
97   PetscErrorCode  ierr;
98   PetscInt        i;
99 
100   PetscFunctionBegin;
101 
102   ierr = DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]);CHKERRQ(ierr);
103   for (i = 1; i < nlevels; i++) {
104     ierr = DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]);CHKERRQ(ierr);
105   }
106   PetscFunctionReturn(0);
107 }
108 
109 #undef __FUNCT__
110 #define __FUNCT__ "DMCoarsenHierarchy_Moab"
111 /*@
112   DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy
113   by succesively coarsening a refined mesh.
114 
115   Collective on MPI_Comm
116 
117   Input Parameter:
118 + dm  - The DMMoab object
119 
120   Output Parameter:
121 + nlevels   - The number of levels of refinement needed to generate the hierarchy
122 . dmc  - The DM objects after successive coarsening of the hierarchy
123 
124   Level: beginner
125 
126 .keywords: DMMoab, generate, coarsening
127 @*/
128 PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[])
129 {
130   PetscErrorCode  ierr;
131   PetscInt        i;
132 
133   PetscFunctionBegin;
134 
135   ierr = DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]);CHKERRQ(ierr);
136   for (i = 1; i < nlevels; i++) {
137     ierr = DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]);CHKERRQ(ierr);
138   }
139   PetscFunctionReturn(0);
140 }
141 
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "DMCreateInterpolation_Moab"
145 /*@
146   DMCreateInterpolation_Moab - Generate the interpolation operators to transform
147   operators (matrices, vectors) from parent level to child level as defined by
148   the DM inputs provided by the user.
149 
150   Collective on MPI_Comm
151 
152   Input Parameter:
153 + dm1  - The DMMoab object
154 - dm2  - the second, finer DMMoab object
155 
156   Output Parameter:
157 + interpl  - The interpolation operator for transferring data between the levels
158 - vec      - The scaling vector (optional)
159 
160   Level: developer
161 
162 .keywords: DMMoab, create, refinement
163 @*/
164 PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat* interpl, Vec* vec)
165 {
166   DM_Moab         *dmbp, *dmbc;
167   PetscErrorCode   ierr;
168   moab::ErrorCode  merr;
169   PetscInt         dim;
170   PetscReal        factor;
171   PetscInt         innz, *nnz, ionz, *onz;
172   PetscInt         nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc;
173   moab::Range      eowned;
174   const PetscBool  use_consistent_bases=PETSC_TRUE;
175 
176   PetscFunctionBegin;
177   PetscValidHeaderSpecific(dmp, DM_CLASSID, 1);
178   PetscValidHeaderSpecific(dmc, DM_CLASSID, 2);
179   dmbp = (DM_Moab*)(dmp)->data;
180   dmbc = (DM_Moab*)(dmc)->data;
181   nlsizp = dmbp->nloc;// *dmb1->numFields;
182   nlsizc = dmbc->nloc;// *dmb2->numFields;
183   ngsizp = dmbp->n; // *dmb1->numFields;
184   ngsizc = dmbc->n; // *dmb2->numFields;
185   nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields;
186 
187   // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent))
188   // Size: nlsizc * nlghsizp
189   PetscInfo4(NULL, "Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n", ngsizc, nlghsizp, dmbp->hlevel, dmbc->hlevel);
190 
191   ierr = DMGetDimension(dmp, &dim);CHKERRQ(ierr);
192 
193   /* allocate the nnz, onz arrays based on block size and local nodes */
194   ierr = PetscCalloc2(nlsizc, &nnz, nlsizc, &onz);CHKERRQ(ierr);
195 
196   eowned = *dmbp->elocal;
197 #ifdef MOAB_HAVE_MPI
198   merr = dmbp->pcomm->filter_pstatus(eowned, PSTATUS_NOT_OWNED, PSTATUS_NOT); MBERRNM(merr);
199 #endif
200   /* Loop through the local elements and compute the relation between the current parent and the refined_level. */
201   for (moab::Range::iterator iter = eowned.begin(); iter != eowned.end(); iter++) {
202 
203     const moab::EntityHandle ehandle = *iter;
204     std::vector<moab::EntityHandle> children;
205     std::vector<moab::EntityHandle> connp, connc;
206     moab::Range connc_owned;
207 
208     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
209     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); MBERRNM(merr);
210 
211     /* Get connectivity of the parent elements */
212     //merr = dmb1->mbiface->get_connectivity(ehandle,connect,vpere,false);MBERRNM(merr);
213     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); MBERRNM(merr);
214     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc_owned); MBERRNM(merr);
215 
216     for (unsigned tc = 0; tc < connc_owned.size(); tc++) {
217       connc.push_back(connc_owned[tc]);
218     }
219 
220     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
221     /* get DoFs for both coarse and fine scale grid */
222     ierr = DMMoabGetDofsBlockedLocal(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr);
223     ierr = DMMoabGetDofsBlockedLocal(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr);
224 
225     for (unsigned tp = 0; tp < connp.size(); tp++) {
226       if (dmbp->vowned->find(connp[tp]) != dmbp->vowned->end()) {
227         for (unsigned tc = 0; tc < connc.size(); tc++) {
228           if (dmbc->vowned->find(connc[tc]) != dmbc->vowned->end()) {
229             /* FIXME: This will over-allocate since we multi-count shared nodes.
230                Disadvantage of using element traversal for nnz computation. */
231             nnz[dofsc[tc]]++;
232           }
233         }
234       }
235       else if (dmbp->vghost->find(connp[tp]) != dmbp->vghost->end()) {
236         for (unsigned tc = 0; tc < connc.size(); tc++) {
237           if (dmbc->vowned->find(connc[tc]) != dmbc->vowned->end()) {
238             /* FIXME: This will over-allocate since we multi-count shared nodes.
239                Disadvantage of using element traversal for nnz computation. */
240             onz[dofsc[tc]]++;
241           }
242         }
243       }
244       else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid entity in child level %D\n", connp[tp]);
245     }
246   }
247 
248   ionz = onz[0];
249   innz = nnz[0];
250   for (int tc = 0; tc < nlsizc; tc++) {
251     // check for maximum allowed sparsity = fully dense
252     nnz[tc] = std::min(nlsizp, nnz[tc]);
253     onz[tc] = std::min(nlghsizp - nlsizp, onz[tc]);
254 
255     innz = (innz < nnz[tc] ? nnz[tc] : innz);
256     ionz = (ionz < onz[tc] ? onz[tc] : ionz);
257   }
258 
259   /* create interpolation matrix */
260   ierr = MatCreate(PetscObjectComm((PetscObject)dmc), interpl);CHKERRQ(ierr);
261   ierr = MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp);CHKERRQ(ierr);
262   //ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr);
263   ierr = MatSetType(*interpl, MATAIJ);CHKERRQ(ierr);
264   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
265 
266   ierr = MatSeqAIJSetPreallocation(*interpl, innz, nnz);CHKERRQ(ierr);
267   ierr = MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz);CHKERRQ(ierr);
268 
269   /* clean up temporary memory */
270   ierr = PetscFree2(nnz, onz);CHKERRQ(ierr);
271 
272   /* set up internal matrix data-structures */
273   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
274 
275   /* Define variables for assembly */
276   std::vector<moab::EntityHandle> children;
277   std::vector<moab::EntityHandle> connp, connc;
278   std::vector<PetscReal> pcoords, ccoords, values_phi;
279   moab::Range connc_owned;
280 
281   if (use_consistent_bases) {
282     const moab::EntityHandle ehandle = eowned.front();
283 
284     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); MBERRNM(merr);
285 
286     /* Get connectivity and coordinates of the parent vertices */
287     connc_owned.clear();
288     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); MBERRNM(merr);
289     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc_owned); MBERRNM(merr);
290 
291     for (unsigned tc = 0; tc < connc_owned.size(); tc++) {
292       connc.push_back(connc_owned[tc]);
293     }
294 
295     std::vector<PetscReal> natparam(3*connc.size(), 0.0);
296     pcoords.resize(connp.size() * 3);
297     ccoords.resize(connc.size() * 3);
298     values_phi.resize(connp.size()*connc.size());
299     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
300     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]); MBERRNM(merr);
301     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]); MBERRNM(merr);
302 
303     /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
304     for (unsigned tc = 0; tc < connc.size(); tc++) {
305       const PetscInt offset = tc * 3;
306 
307       /* Scale ccoords relative to pcoords */
308       ierr = DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size()*tc]);CHKERRQ(ierr);
309     }
310   }
311   else {
312     factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0);
313   }
314 
315   /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */
316   // ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
317 
318   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
319   for (moab::Range::iterator iter = eowned.begin(); iter != eowned.end(); iter++) {
320 
321     const moab::EntityHandle ehandle = *iter;
322 
323     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
324     children.clear();
325     connc.clear();
326     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); MBERRNM(merr);
327 
328     /* Get connectivity and coordinates of the parent vertices */
329     connc_owned.clear();
330     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); MBERRNM(merr);
331     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc_owned); MBERRNM(merr);
332     for (unsigned tc = 0; tc < connc_owned.size(); tc++) {
333       connc.push_back(connc_owned[tc]);
334     }
335 
336     pcoords.resize(connp.size() * 3);
337     ccoords.resize(connc.size() * 3);
338     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
339     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]); MBERRNM(merr);
340     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]); MBERRNM(merr);
341 
342     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
343     /* TODO: specific to scalar system - use GetDofs */
344     ierr = DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr);
345     ierr = DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr);
346 
347     /* Compute the actual interpolation weights when projecting solution/residual between levels */
348     if (use_consistent_bases) {
349 
350       /* Use the cached values of natural parameteric coordinates and basis pre-evaluated.
351          We are making an assumption here that UMR used in GMG to generate the hierarchy uses
352          the same template for all elements; This will fail for mixed element meshes (TRI/QUAD).
353 
354          TODO: Fix the above assumption by caching data for families
355       */
356 
357       /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
358       for (unsigned tc = 0; tc < connc.size(); tc++) {
359         /* TODO: Check if we should be using INSERT_VALUES instead */
360         ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[connp.size()*tc], ADD_VALUES);CHKERRQ(ierr);
361       }
362     }
363     else {
364       /* Compute the interpolation weights by determining distance of 1-ring
365          neighbor vertices from current vertex
366 
367          This should be used only when FEM basis is not used for the discretization.
368          Else, the consistent interface to compute the basis function for interpolation
369          between the levels should be evaluated correctly to preserve convergence of GMG.
370          RBF basis will be terrible for any unsmooth problems..
371       */
372       values_phi.resize(connp.size());
373       for (unsigned tc = 0; tc < connc.size(); tc++) {
374 
375         PetscReal normsum = 0.0;
376         for (unsigned tp = 0; tp < connp.size(); tp++) {
377           values_phi[tp] = 0.0;
378           for (unsigned k = 0; k < 3; k++)
379             values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim);
380           if (values_phi[tp] < 1e-12) {
381             values_phi[tp] = 1e12;
382           }
383           else {
384             //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim);
385             values_phi[tp] = std::pow(values_phi[tp], -1.0);
386             normsum += values_phi[tp];
387           }
388         }
389         for (unsigned tp = 0; tp < connp.size(); tp++) {
390           if (values_phi[tp] > 1e11)
391             values_phi[tp] = factor * 0.5 / connp.size();
392           else
393             values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum);
394         }
395         ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
396       }
397     }
398   }
399   if (vec) *vec = NULL;
400   ierr = MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
401   ierr = MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
402   PetscFunctionReturn(0);
403 }
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "DMCreateInjection_Moab"
407 /*@
408   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
409   by succesively refining a coarse mesh, already defined in the DM object
410   provided by the user.
411 
412   Collective on MPI_Comm
413 
414   Input Parameter:
415 . dmb  - The DMMoab object
416 
417   Output Parameter:
418 . nlevels   - The number of levels of refinement needed to generate the hierarchy
419 + ldegrees  - The degree of refinement at each level in the hierarchy
420 
421   Level: beginner
422 
423 .keywords: DMMoab, create, refinement
424 @*/
425 PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter* ctx)
426 {
427   //DM_Moab        *dmmoab;
428 
429   PetscFunctionBegin;
430   PetscValidHeaderSpecific(dm1, DM_CLASSID, 1);
431   PetscValidHeaderSpecific(dm2, DM_CLASSID, 2);
432   //dmmoab = (DM_Moab*)(dm1)->data;
433 
434   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
435   PetscFunctionReturn(0);
436 }
437 
438 
439 #undef __FUNCT__
440 #define __FUNCT__ "DM_UMR_Moab_Private"
441 PetscErrorCode  DM_UMR_Moab_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref)
442 {
443   PetscErrorCode  ierr;
444   PetscInt        i, dim;
445   DM              dm2;
446   moab::ErrorCode merr;
447   DM_Moab        *dmb = (DM_Moab*)dm->data, *dd2;
448 
449   PetscFunctionBegin;
450   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
451   PetscValidPointer(dmref, 4);
452 
453   if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) {
454     if (dmb->hlevel + 1 > dmb->nhlevels && refine) PetscInfo2(NULL, "Invalid multigrid refinement hierarchy level specified (%D). MOAB UMR max levels = %D. Creating a NULL object.\n", dmb->hlevel + 1, dmb->nhlevels);
455     if (dmb->hlevel - 1 < 0 && !refine) PetscInfo1(NULL, "Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n", dmb->hlevel - 1);
456     *dmref = PETSC_NULL;
457     PetscFunctionReturn(0);
458   }
459 
460   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
461   dd2 = (DM_Moab*)dm2->data;
462 
463   dd2->mbiface = dmb->mbiface;
464 #ifdef MOAB_HAVE_MPI
465   dd2->pcomm = dmb->pcomm;
466 #endif
467   dd2->icreatedinstance = PETSC_FALSE;
468   dd2->nghostrings = dmb->nghostrings;
469 
470   /* set the new level based on refinement/coarsening */
471   if (refine) {
472     dd2->hlevel = dmb->hlevel + 1;
473   }
474   else {
475     dd2->hlevel = dmb->hlevel - 1;
476   }
477 
478   /* Copy the multilevel hierarchy pointers in MOAB */
479   dd2->hierarchy = dmb->hierarchy;
480   dd2->nhlevels = dmb->nhlevels;
481   ierr = PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets);CHKERRQ(ierr);
482   for (i = 0; i <= dd2->nhlevels; i++) {
483     dd2->hsets[i] = dmb->hsets[i];
484   }
485   dd2->fileset = dd2->hsets[dd2->hlevel];
486 
487   /* do the remaining initializations for DMMoab */
488   dd2->bs = dmb->bs;
489   dd2->numFields = dmb->numFields;
490   dd2->rw_dbglevel = dmb->rw_dbglevel;
491   dd2->partition_by_rank = dmb->partition_by_rank;
492   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
493   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
494   dd2->read_mode = dmb->read_mode;
495   dd2->write_mode = dmb->write_mode;
496 
497   /* set global ID tag handle */
498   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
499 
500   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag); MBERRNM(merr);
501 
502   ierr = DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix);CHKERRQ(ierr);
503   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
504   ierr = DMSetDimension(dm2, dim);CHKERRQ(ierr);
505 
506   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
507   dm2->ops->creatematrix = dm->ops->creatematrix;
508 
509   /* copy fill information if given */
510   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
511 
512   /* copy vector type information */
513   ierr = DMSetMatType(dm2, dm->mattype);CHKERRQ(ierr);
514   ierr = DMSetVecType(dm2, dm->vectype);CHKERRQ(ierr);
515   dd2->numFields = dmb->numFields;
516   if (dmb->numFields) {
517     ierr = DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames);CHKERRQ(ierr);
518   }
519 
520   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
521 
522   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
523   ierr = DMSetUp(dm2);CHKERRQ(ierr);
524 
525   *dmref = dm2;
526   PetscFunctionReturn(0);
527 }
528 
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "DMRefine_Moab"
532 /*@
533   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
534   by succesively refining a coarse mesh, already defined in the DM object
535   provided by the user.
536 
537   Collective on DM
538 
539   Input Parameter:
540 + dm  - The DMMoab object
541 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
542 
543   Output Parameter:
544 . dmf - the refined DM, or NULL
545 
546   Note: If no refinement was done, the return value is NULL
547 
548   Level: developer
549 
550 .keywords: DMMoab, create, refinement
551 @*/
552 PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM* dmf)
553 {
554   PetscErrorCode  ierr;
555 
556   PetscFunctionBegin;
557   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
558 
559   ierr = DM_UMR_Moab_Private(dm, comm, PETSC_TRUE, dmf);CHKERRQ(ierr);
560   PetscFunctionReturn(0);
561 }
562 
563 #undef __FUNCT__
564 #define __FUNCT__ "DMCoarsen_Moab"
565 /*@
566   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
567   by succesively refining a coarse mesh, already defined in the DM object
568   provided by the user.
569 
570   Collective on DM
571 
572   Input Parameter:
573 . dm  - The DMMoab object
574 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
575 
576   Output Parameter:
577 . dmf - the coarsened DM, or NULL
578 
579   Note: If no coarsening was done, the return value is NULL
580 
581   Level: developer
582 
583 .keywords: DMMoab, create, coarsening
584 @*/
585 PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM* dmc)
586 {
587   PetscErrorCode  ierr;
588 
589   PetscFunctionBegin;
590   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
591 
592   ierr = DM_UMR_Moab_Private(dm, comm, PETSC_FALSE, dmc);CHKERRQ(ierr);
593   PetscFunctionReturn(0);
594 }
595