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