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: . dmb - 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: @*/
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;
31: dmmoab = (DM_Moab*)(dm)->data;
33: if (!ldegrees) {
34: PetscMalloc1(nlevels, &pdegrees);
35: for (ilevel = 0; ilevel < nlevels; ilevel++) pdegrees[ilevel] = 2; /* default = Degree 2 refinement */
36: }
37: else pdegrees = ldegrees;
39: /* initialize set level refinement data for hierarchy */
40: dmmoab->nhlevels = nlevels;
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
49: PetscMalloc1(nlevels + 1, &dmmoab->hsets);
51: /* generate the mesh hierarchy */
52: merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false);MBERRNM(merr);
54: #ifdef MOAB_HAVE_MPI
55: if (dmmoab->pcomm->size() > 1) {
56: merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);MBERRNM(merr);
57: }
58: #endif
60: /* copy the mesh sets for nested refinement hierarchy */
61: dmmoab->hsets[0] = hsets[0];
62: for (ilevel = 1; ilevel <= nlevels; ilevel++)
63: {
64: dmmoab->hsets[ilevel] = hsets[ilevel];
66: #ifdef MOAB_HAVE_MPI
67: merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false);MBERRNM(merr);
68: #endif
70: /* Update material and other geometric tags from parent to child sets */
71: merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]);MBERRNM(merr);
72: }
74: hsets.clear();
75: if (!ldegrees) {
76: PetscFree(pdegrees);
77: }
78: return 0;
79: }
81: /*@C
82: DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy
83: by succesively refining a coarse mesh.
85: Collective
87: Input Parameter:
88: . dm - The DMMoab object
90: Output Parameters:
91: + nlevels - The number of levels of refinement needed to generate the hierarchy
92: - dmf - The DM objects after successive refinement of the hierarchy
94: Level: beginner
96: @*/
97: PETSC_EXTERN PetscErrorCode DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[])
98: {
99: PetscInt i;
102: DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]);
103: for (i = 1; i < nlevels; i++) {
104: DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]);
105: }
106: return 0;
107: }
109: /*@C
110: DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy
111: by succesively coarsening a refined mesh.
113: Collective
115: Input Parameter:
116: . dm - The DMMoab object
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
122: Level: beginner
124: @*/
125: PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[])
126: {
127: PetscInt i;
130: DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]);
131: for (i = 1; i < nlevels; i++) {
132: DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]);
133: }
134: return 0;
135: }
137: PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt*, PetscInt*, PetscInt*, PetscInt*, PetscBool);
139: /*@C
140: DMCreateInterpolation_Moab - Generate the interpolation operators to transform
141: operators (matrices, vectors) from parent level to child level as defined by
142: the DM inputs provided by the user.
144: Collective
146: Input Parameters:
147: + dm1 - The DMMoab object
148: - dm2 - the second, finer DMMoab object
150: Output Parameters:
151: + interpl - The interpolation operator for transferring data between the levels
152: - vec - The scaling vector (optional)
154: 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;
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;
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 %D X %D to apply transformation between levels %D -> %D.\n", ngsizc, nlghsizp, dmbp->hlevel, dmbc->hlevel);
182: DMGetDimension(dmp, &dim);
184: /* allocate the nnz, onz arrays based on block size and local nodes */
185: PetscCalloc2(nlsizc, &nnz, nlsizc, &onz);
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++) {
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);MBERRNM(merr);
204: /* loop over vertices and update the number of connectivity */
205: 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);MBERRNM(merr);
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);MBERRNM(merr);
216: 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: }
228: for (int i = 0; i < nlsizc; i++)
229: 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: 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: MatCreate(PetscObjectComm((PetscObject)dmc), interpl);
246: MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp);
247: MatSetType(*interpl, MATAIJ);
248: MatSetFromOptions(*interpl);
250: MatSeqAIJSetPreallocation(*interpl, innz, nnz);
251: MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz);
253: /* clean up temporary memory */
254: PetscFree2(nnz, onz);
256: /* set up internal matrix data-structures */
257: 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);MBERRNM(merr);
269: /* Get connectivity and coordinates of the parent vertices */
270: merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
271: merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr);
273: std::vector<PetscReal> natparam(3*connc.size(), 0.0);
274: pcoords.resize(connp.size() * 3);
275: ccoords.resize(connc.size() * 3);
276: values_phi.resize(connp.size()*connc.size());
277: /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
278: merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr);
279: merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr);
281: /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
282: for (unsigned tc = 0; tc < connc.size(); tc++) {
283: const PetscInt offset = tc * 3;
285: /* Scale ccoords relative to pcoords */
286: DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size()*tc]);
287: }
288: }
289: else {
290: factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0);
291: }
293: /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */
294: MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
296: /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
297: for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) {
299: const moab::EntityHandle ehandle = *iter;
301: /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
302: children.clear();
303: connc.clear();
304: merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
306: /* Get connectivity and coordinates of the parent vertices */
307: merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
308: merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr);
310: pcoords.resize(connp.size() * 3);
311: ccoords.resize(connc.size() * 3);
312: /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
313: merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr);
314: merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr);
316: std::vector<int> dofsp(connp.size()), dofsc(connc.size());
317: /* TODO: specific to scalar system - use GetDofs */
318: DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0]);
319: DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0]);
321: /* Compute the actual interpolation weights when projecting solution/residual between levels */
322: if (use_consistent_bases) {
324: /* Use the cached values of natural parameteric coordinates and basis pre-evaluated.
325: We are making an assumption here that UMR used in GMG to generate the hierarchy uses
326: the same template for all elements; This will fail for mixed element meshes (TRI/QUAD).
328: TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes)
329: */
331: /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
332: for (unsigned tc = 0; tc < connc.size(); tc++) {
333: /* TODO: Check if we should be using INSERT_VALUES instead */
334: MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[connp.size()*tc], ADD_VALUES);
335: }
336: }
337: else {
338: /* Compute the interpolation weights by determining distance of 1-ring
339: neighbor vertices from current vertex
341: This should be used only when FEM basis is not used for the discretization.
342: Else, the consistent interface to compute the basis function for interpolation
343: between the levels should be evaluated correctly to preserve convergence of GMG.
344: Shephard's basis will be terrible for any unsmooth problems.
345: */
346: values_phi.resize(connp.size());
347: for (unsigned tc = 0; tc < connc.size(); tc++) {
349: PetscReal normsum = 0.0;
350: for (unsigned tp = 0; tp < connp.size(); tp++) {
351: values_phi[tp] = 0.0;
352: for (unsigned k = 0; k < 3; k++)
353: values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim);
354: if (values_phi[tp] < 1e-12) {
355: values_phi[tp] = 1e12;
356: }
357: else {
358: //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim);
359: values_phi[tp] = std::pow(values_phi[tp], -1.0);
360: normsum += values_phi[tp];
361: }
362: }
363: for (unsigned tp = 0; tp < connp.size(); tp++) {
364: if (values_phi[tp] > 1e11)
365: values_phi[tp] = factor * 0.5 / connp.size();
366: else
367: values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum);
368: }
369: MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES);
370: }
371: }
372: }
373: if (vec) *vec = NULL;
374: MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY);
375: MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY);
376: return 0;
377: }
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.
384: Collective
386: Input Parameter:
387: . dmb - The DMMoab object
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
393: Level: beginner
395: @*/
396: PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter* ctx)
397: {
398: //DM_Moab *dmmoab;
402: //dmmoab = (DM_Moab*)(dm1)->data;
404: PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
405: return 0;
406: }
408: static PetscErrorCode DMMoab_UMR_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref)
409: {
410: PetscInt i, dim;
411: DM dm2;
412: moab::ErrorCode merr;
413: DM_Moab *dmb = (DM_Moab*)dm->data, *dd2;
418: if ((dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine)) {
419: if (dmb->hlevel + 1 > dmb->nhlevels && refine) PetscInfo(NULL, "Invalid multigrid refinement hierarchy level specified (%D). MOAB UMR max levels = %D. Creating a NULL object.\n", dmb->hlevel + 1, dmb->nhlevels);
420: if (dmb->hlevel - 1 < 0 && !refine) PetscInfo(NULL, "Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n", dmb->hlevel - 1);
421: *dmref = PETSC_NULL;
422: return 0;
423: }
425: DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);
426: dd2 = (DM_Moab*)dm2->data;
428: dd2->mbiface = dmb->mbiface;
429: #ifdef MOAB_HAVE_MPI
430: dd2->pcomm = dmb->pcomm;
431: #endif
432: dd2->icreatedinstance = PETSC_FALSE;
433: dd2->nghostrings = dmb->nghostrings;
435: /* set the new level based on refinement/coarsening */
436: if (refine) {
437: dd2->hlevel = dmb->hlevel + 1;
438: }
439: else {
440: dd2->hlevel = dmb->hlevel - 1;
441: }
443: /* Copy the multilevel hierarchy pointers in MOAB */
444: dd2->hierarchy = dmb->hierarchy;
445: dd2->nhlevels = dmb->nhlevels;
446: PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets);
447: for (i = 0; i <= dd2->nhlevels; i++) {
448: dd2->hsets[i] = dmb->hsets[i];
449: }
450: dd2->fileset = dd2->hsets[dd2->hlevel];
452: /* do the remaining initializations for DMMoab */
453: dd2->bs = dmb->bs;
454: dd2->numFields = dmb->numFields;
455: dd2->rw_dbglevel = dmb->rw_dbglevel;
456: dd2->partition_by_rank = dmb->partition_by_rank;
457: PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);
458: PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);
459: dd2->read_mode = dmb->read_mode;
460: dd2->write_mode = dmb->write_mode;
462: /* set global ID tag handle */
463: DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);
465: merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
467: DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix);
468: DMGetDimension(dm, &dim);
469: DMSetDimension(dm2, dim);
471: /* allow overloaded (user replaced) operations to be inherited by refinement clones */
472: dm2->ops->creatematrix = dm->ops->creatematrix;
474: /* copy fill information if given */
475: DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);
477: /* copy vector type information */
478: DMSetMatType(dm2, dm->mattype);
479: DMSetVecType(dm2, dm->vectype);
480: dd2->numFields = dmb->numFields;
481: if (dmb->numFields) {
482: DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames);
483: }
485: DMSetFromOptions(dm2);
487: /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
488: DMSetUp(dm2);
490: *dmref = dm2;
491: return 0;
492: }
494: /*@C
495: DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
496: by succesively refining a coarse mesh, already defined in the DM object
497: provided by the user.
499: Collective on dm
501: Input Parameters:
502: + dm - The DMMoab object
503: - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
505: Output Parameter:
506: . dmf - the refined DM, or NULL
508: Note: If no refinement was done, the return value is NULL
510: Level: developer
512: @*/
513: PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM* dmf)
514: {
517: DMMoab_UMR_Private(dm, comm, PETSC_TRUE, dmf);
518: return 0;
519: }
521: /*@C
522: DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
523: by succesively refining a coarse mesh, already defined in the DM object
524: provided by the user.
526: Collective on dm
528: Input Parameters:
529: + dm - The DMMoab object
530: - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
532: Output Parameter:
533: . dmf - the coarsened DM, or NULL
535: Note: If no coarsening was done, the return value is NULL
537: Level: developer
539: @*/
540: PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM* dmc)
541: {
543: DMMoab_UMR_Private(dm, comm, PETSC_FALSE, dmc);
544: return 0;
545: }