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