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