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