xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1755f3dfbSVijay Mahadevan #include <petsc/private/dmmbimpl.h>
2b117cd09SVijay Mahadevan #include <petscdmmoab.h>
3b117cd09SVijay Mahadevan #include <MBTagConventions.hpp>
4b117cd09SVijay Mahadevan #include <moab/NestedRefine.hpp>
5b117cd09SVijay Mahadevan 
6cab5ea25SPierre Jolivet /*@C
7b117cd09SVijay Mahadevan   DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy
8b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
9b117cd09SVijay Mahadevan   provided by the user.
10b117cd09SVijay Mahadevan 
11d083f849SBarry Smith   Collective
12b117cd09SVijay Mahadevan 
13b117cd09SVijay Mahadevan   Input Parameter:
14a2b725a8SWilliam Gropp . dmb  - The DMMoab object
15b117cd09SVijay Mahadevan 
16d8d19677SJose E. Roman   Output Parameters:
17b117cd09SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
18a2b725a8SWilliam Gropp - ldegrees  - The degree of refinement at each level in the hierarchy
19b117cd09SVijay Mahadevan 
20b117cd09SVijay Mahadevan   Level: beginner
21b117cd09SVijay Mahadevan 
22b117cd09SVijay Mahadevan @*/
239371c9d4SSatish Balay PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees) {
24b117cd09SVijay Mahadevan   DM_Moab                        *dmmoab;
25b117cd09SVijay Mahadevan   moab::ErrorCode                 merr;
262417220eSVijay Mahadevan   PetscInt                       *pdegrees, ilevel;
27e882eb38SVijay Mahadevan   std::vector<moab::EntityHandle> hsets;
28b117cd09SVijay Mahadevan 
29b117cd09SVijay Mahadevan   PetscFunctionBegin;
30b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
31b117cd09SVijay Mahadevan   dmmoab = (DM_Moab *)(dm)->data;
32b117cd09SVijay Mahadevan 
33b117cd09SVijay Mahadevan   if (!ldegrees) {
349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nlevels, &pdegrees));
352417220eSVijay Mahadevan     for (ilevel = 0; ilevel < nlevels; ilevel++) pdegrees[ilevel] = 2; /* default = Degree 2 refinement */
369371c9d4SSatish Balay   } else pdegrees = ldegrees;
37b117cd09SVijay Mahadevan 
38b117cd09SVijay Mahadevan   /* initialize set level refinement data for hierarchy */
39b117cd09SVijay Mahadevan   dmmoab->nhlevels = nlevels;
40b117cd09SVijay Mahadevan 
41b117cd09SVijay Mahadevan   /* Instantiate the nested refinement class */
429daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
433f1c6e43SVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core *>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset);
449daf19fdSVijay Mahadevan #else
459daf19fdSVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core *>(dmmoab->mbiface), NULL, dmmoab->fileset);
469daf19fdSVijay Mahadevan #endif
47b117cd09SVijay Mahadevan 
489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nlevels + 1, &dmmoab->hsets));
49e882eb38SVijay Mahadevan 
50b117cd09SVijay Mahadevan   /* generate the mesh hierarchy */
519371c9d4SSatish Balay   merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false);
529371c9d4SSatish Balay   MBERRNM(merr);
53e882eb38SVijay Mahadevan 
549daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
55755f3dfbSVijay Mahadevan   if (dmmoab->pcomm->size() > 1) {
569371c9d4SSatish Balay     merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);
579371c9d4SSatish Balay     MBERRNM(merr);
58755f3dfbSVijay Mahadevan   }
599daf19fdSVijay Mahadevan #endif
6049d66b22SVijay Mahadevan 
6164e1c140SVijay Mahadevan   /* copy the mesh sets for nested refinement hierarchy */
62e92d1c7cSVijay Mahadevan   dmmoab->hsets[0] = hsets[0];
639371c9d4SSatish Balay   for (ilevel = 1; ilevel <= nlevels; ilevel++) {
642417220eSVijay Mahadevan     dmmoab->hsets[ilevel] = hsets[ilevel];
65e92d1c7cSVijay Mahadevan 
669c368985SVijay Mahadevan #ifdef MOAB_HAVE_MPI
679371c9d4SSatish Balay     merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false);
689371c9d4SSatish Balay     MBERRNM(merr);
699c368985SVijay Mahadevan #endif
70e92d1c7cSVijay Mahadevan 
71e92d1c7cSVijay Mahadevan     /* Update material and other geometric tags from parent to child sets */
729371c9d4SSatish Balay     merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]);
739371c9d4SSatish Balay     MBERRNM(merr);
742417220eSVijay Mahadevan   }
7564e1c140SVijay Mahadevan 
7664e1c140SVijay Mahadevan   hsets.clear();
77*48a46eb9SPierre Jolivet   if (!ldegrees) PetscCall(PetscFree(pdegrees));
78b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
79b117cd09SVijay Mahadevan }
80b117cd09SVijay Mahadevan 
81cab5ea25SPierre Jolivet /*@C
82e882eb38SVijay Mahadevan   DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy
83e882eb38SVijay Mahadevan   by succesively refining a coarse mesh.
84b117cd09SVijay Mahadevan 
85d083f849SBarry Smith   Collective
86b117cd09SVijay Mahadevan 
87b117cd09SVijay Mahadevan   Input Parameter:
88a2b725a8SWilliam Gropp . dm  - The DMMoab object
89b117cd09SVijay Mahadevan 
90d8d19677SJose E. Roman   Output Parameters:
91e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
92a2b725a8SWilliam Gropp - dmf  - The DM objects after successive refinement of the hierarchy
93b117cd09SVijay Mahadevan 
94b117cd09SVijay Mahadevan   Level: beginner
95b117cd09SVijay Mahadevan 
96b117cd09SVijay Mahadevan @*/
979371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[]) {
98b117cd09SVijay Mahadevan   PetscInt i;
99b117cd09SVijay Mahadevan 
100b117cd09SVijay Mahadevan   PetscFunctionBegin;
101b117cd09SVijay Mahadevan 
1029566063dSJacob Faibussowitsch   PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]));
103*48a46eb9SPierre Jolivet   for (i = 1; i < nlevels; i++) PetscCall(DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]));
104b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
105b117cd09SVijay Mahadevan }
106b117cd09SVijay Mahadevan 
107cab5ea25SPierre Jolivet /*@C
108e882eb38SVijay Mahadevan   DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy
109e882eb38SVijay Mahadevan   by succesively coarsening a refined mesh.
110b117cd09SVijay Mahadevan 
111d083f849SBarry Smith   Collective
112b117cd09SVijay Mahadevan 
113b117cd09SVijay Mahadevan   Input Parameter:
114a2b725a8SWilliam Gropp . dm  - The DMMoab object
115b117cd09SVijay Mahadevan 
116d8d19677SJose E. Roman   Output Parameters:
117e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
118a2b725a8SWilliam Gropp - dmc  - The DM objects after successive coarsening of the hierarchy
119b117cd09SVijay Mahadevan 
120b117cd09SVijay Mahadevan   Level: beginner
121b117cd09SVijay Mahadevan 
122b117cd09SVijay Mahadevan @*/
1239371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[]) {
124b117cd09SVijay Mahadevan   PetscInt i;
125b117cd09SVijay Mahadevan 
126b117cd09SVijay Mahadevan   PetscFunctionBegin;
127b117cd09SVijay Mahadevan 
1289566063dSJacob Faibussowitsch   PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]));
129*48a46eb9SPierre Jolivet   for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]));
130b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
131b117cd09SVijay Mahadevan }
132b117cd09SVijay Mahadevan 
1332417220eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscBool);
134b117cd09SVijay Mahadevan 
135cab5ea25SPierre Jolivet /*@C
136e882eb38SVijay Mahadevan   DMCreateInterpolation_Moab - Generate the interpolation operators to transform
137e882eb38SVijay Mahadevan   operators (matrices, vectors) from parent level to child level as defined by
138e882eb38SVijay Mahadevan   the DM inputs provided by the user.
139b117cd09SVijay Mahadevan 
140d083f849SBarry Smith   Collective
141b117cd09SVijay Mahadevan 
142d8d19677SJose E. Roman   Input Parameters:
143e882eb38SVijay Mahadevan + dm1  - The DMMoab object
144e882eb38SVijay Mahadevan - dm2  - the second, finer DMMoab object
145b117cd09SVijay Mahadevan 
146d8d19677SJose E. Roman   Output Parameters:
147e882eb38SVijay Mahadevan + interpl  - The interpolation operator for transferring data between the levels
148e882eb38SVijay Mahadevan - vec      - The scaling vector (optional)
149b117cd09SVijay Mahadevan 
150e882eb38SVijay Mahadevan   Level: developer
151b117cd09SVijay Mahadevan 
152b117cd09SVijay Mahadevan @*/
1539371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat *interpl, Vec *vec) {
154755f3dfbSVijay Mahadevan   DM_Moab        *dmbp, *dmbc;
155b117cd09SVijay Mahadevan   moab::ErrorCode merr;
156e882eb38SVijay Mahadevan   PetscInt        dim;
1573f1c6e43SVijay Mahadevan   PetscReal       factor;
158ce27a4eeSVijay Mahadevan   PetscInt        innz, *nnz, ionz, *onz;
159755f3dfbSVijay Mahadevan   PetscInt        nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc;
160a86ed394SVijay Mahadevan   const PetscBool use_consistent_bases = PETSC_TRUE;
161b117cd09SVijay Mahadevan 
162b117cd09SVijay Mahadevan   PetscFunctionBegin;
163755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmp, DM_CLASSID, 1);
164755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmc, DM_CLASSID, 2);
165755f3dfbSVijay Mahadevan   dmbp     = (DM_Moab *)(dmp)->data;
166755f3dfbSVijay Mahadevan   dmbc     = (DM_Moab *)(dmc)->data;
167755f3dfbSVijay Mahadevan   nlsizp   = dmbp->nloc;                  // *dmb1->numFields;
168755f3dfbSVijay Mahadevan   nlsizc   = dmbc->nloc;                  // *dmb2->numFields;
169755f3dfbSVijay Mahadevan   ngsizp   = dmbp->n;                     // *dmb1->numFields;
170755f3dfbSVijay Mahadevan   ngsizc   = dmbc->n;                     // *dmb2->numFields;
171755f3dfbSVijay Mahadevan   nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields;
172b117cd09SVijay Mahadevan 
1732417220eSVijay Mahadevan   // Columns = Parent DoFs ;  Rows = Child DoFs
174755f3dfbSVijay Mahadevan   // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent))
175755f3dfbSVijay Mahadevan   // Size: nlsizc * nlghsizp
17663a3b9bcSJacob Faibussowitsch   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);
177b117cd09SVijay Mahadevan 
1789566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dmp, &dim));
179a86ed394SVijay Mahadevan 
180941e0cffSVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
1819566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(nlsizc, &nnz, nlsizc, &onz));
182941e0cffSVijay Mahadevan 
183941e0cffSVijay Mahadevan   /* Loop through the local elements and compute the relation between the current parent and the refined_level. */
1842417220eSVijay Mahadevan   for (moab::Range::iterator iter = dmbc->vowned->begin(); iter != dmbc->vowned->end(); iter++) {
1852417220eSVijay Mahadevan     const moab::EntityHandle        vhandle = *iter;
1862417220eSVijay Mahadevan     /* define local variables */
1872417220eSVijay Mahadevan     moab::EntityHandle              parent;
1882417220eSVijay Mahadevan     std::vector<moab::EntityHandle> adjs;
1892417220eSVijay Mahadevan     moab::Range                     found;
1902417220eSVijay Mahadevan 
1912417220eSVijay Mahadevan     /* store the vertex DoF number */
1922417220eSVijay Mahadevan     const int ldof = dmbc->lidmap[vhandle - dmbc->seqstart];
1932417220eSVijay Mahadevan 
1942417220eSVijay Mahadevan     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
1952417220eSVijay Mahadevan        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
1962417220eSVijay Mahadevan        non-zero pattern accordingly. */
1979371c9d4SSatish Balay     merr = dmbc->hierarchy->get_adjacencies(vhandle, dmbc->dim, adjs);
1989371c9d4SSatish Balay     MBERRNM(merr);
1992417220eSVijay Mahadevan 
2002417220eSVijay Mahadevan     /* loop over vertices and update the number of connectivity */
2012417220eSVijay Mahadevan     for (unsigned jter = 0; jter < adjs.size(); jter++) {
2022417220eSVijay Mahadevan       const moab::EntityHandle jhandle = adjs[jter];
203941e0cffSVijay Mahadevan 
204941e0cffSVijay Mahadevan       /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
2059371c9d4SSatish Balay       merr = dmbc->hierarchy->child_to_parent(jhandle, dmbc->hlevel, dmbp->hlevel, &parent);
2069371c9d4SSatish Balay       MBERRNM(merr);
207941e0cffSVijay Mahadevan 
2082417220eSVijay Mahadevan       /* Get connectivity information in canonical ordering for the local element */
2092417220eSVijay Mahadevan       std::vector<moab::EntityHandle> connp;
2109371c9d4SSatish Balay       merr = dmbp->hierarchy->get_connectivity(parent, dmbp->hlevel, connp);
2119371c9d4SSatish Balay       MBERRNM(merr);
212a86ed394SVijay Mahadevan 
2132417220eSVijay Mahadevan       for (unsigned ic = 0; ic < connp.size(); ++ic) {
2142417220eSVijay Mahadevan         /* loop over each element connected to the adjacent vertex and update as needed */
2152417220eSVijay Mahadevan         /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
2162417220eSVijay Mahadevan         if (found.find(connp[ic]) != found.end()) continue;                    /* make sure we don't double count shared vertices */
2172417220eSVijay Mahadevan         if (dmbp->vghost->find(connp[ic]) != dmbp->vghost->end()) onz[ldof]++; /* update out-of-proc onz */
2182417220eSVijay Mahadevan         else nnz[ldof]++;                                                      /* else local vertex */
2192417220eSVijay Mahadevan         found.insert(connp[ic]);
2202417220eSVijay Mahadevan       }
2212417220eSVijay Mahadevan     }
222941e0cffSVijay Mahadevan   }
223941e0cffSVijay Mahadevan 
2249371c9d4SSatish Balay   for (int i = 0; i < nlsizc; i++) nnz[i] += 1; /* self count the node */
225941e0cffSVijay Mahadevan 
226941e0cffSVijay Mahadevan   ionz = onz[0];
227941e0cffSVijay Mahadevan   innz = nnz[0];
228755f3dfbSVijay Mahadevan   for (int tc = 0; tc < nlsizc; tc++) {
229941e0cffSVijay Mahadevan     // check for maximum allowed sparsity = fully dense
230755f3dfbSVijay Mahadevan     nnz[tc] = std::min(nlsizp, nnz[tc]);
2312417220eSVijay Mahadevan     onz[tc] = std::min(ngsizp - nlsizp, onz[tc]);
2322417220eSVijay Mahadevan 
2337d3de750SJacob Faibussowitsch     PetscInfo(NULL, "  %d: NNZ = %d, ONZ = %d\n", tc, nnz[tc], onz[tc]);
234941e0cffSVijay Mahadevan 
235941e0cffSVijay Mahadevan     innz = (innz < nnz[tc] ? nnz[tc] : innz);
236941e0cffSVijay Mahadevan     ionz = (ionz < onz[tc] ? onz[tc] : ionz);
237941e0cffSVijay Mahadevan   }
23864e1c140SVijay Mahadevan 
23964e1c140SVijay Mahadevan   /* create interpolation matrix */
2409566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc), interpl));
2419566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp));
2429566063dSJacob Faibussowitsch   PetscCall(MatSetType(*interpl, MATAIJ));
2439566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*interpl));
24464e1c140SVijay Mahadevan 
2459566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*interpl, innz, nnz));
2469566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz));
247941e0cffSVijay Mahadevan 
248941e0cffSVijay Mahadevan   /* clean up temporary memory */
2499566063dSJacob Faibussowitsch   PetscCall(PetscFree2(nnz, onz));
250b117cd09SVijay Mahadevan 
251b117cd09SVijay Mahadevan   /* set up internal matrix data-structures */
2529566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*interpl));
253b117cd09SVijay Mahadevan 
254a86ed394SVijay Mahadevan   /* Define variables for assembly */
255a86ed394SVijay Mahadevan   std::vector<moab::EntityHandle> children;
256a86ed394SVijay Mahadevan   std::vector<moab::EntityHandle> connp, connc;
257a86ed394SVijay Mahadevan   std::vector<PetscReal>          pcoords, ccoords, values_phi;
258a86ed394SVijay Mahadevan 
259a86ed394SVijay Mahadevan   if (use_consistent_bases) {
2602417220eSVijay Mahadevan     const moab::EntityHandle ehandle = dmbp->elocal->front();
261a86ed394SVijay Mahadevan 
2629371c9d4SSatish Balay     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);
2639371c9d4SSatish Balay     MBERRNM(merr);
264a86ed394SVijay Mahadevan 
265a86ed394SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
2669371c9d4SSatish Balay     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);
2679371c9d4SSatish Balay     MBERRNM(merr);
2689371c9d4SSatish Balay     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);
2699371c9d4SSatish Balay     MBERRNM(merr);
270a86ed394SVijay Mahadevan 
271a86ed394SVijay Mahadevan     std::vector<PetscReal> natparam(3 * connc.size(), 0.0);
272a86ed394SVijay Mahadevan     pcoords.resize(connp.size() * 3);
273a86ed394SVijay Mahadevan     ccoords.resize(connc.size() * 3);
274a86ed394SVijay Mahadevan     values_phi.resize(connp.size() * connc.size());
275a86ed394SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
2769371c9d4SSatish Balay     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);
2779371c9d4SSatish Balay     MBERRNM(merr);
2789371c9d4SSatish Balay     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);
2799371c9d4SSatish Balay     MBERRNM(merr);
280a86ed394SVijay Mahadevan 
281a86ed394SVijay Mahadevan     /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
282a86ed394SVijay Mahadevan     for (unsigned tc = 0; tc < connc.size(); tc++) {
283a86ed394SVijay Mahadevan       const PetscInt offset = tc * 3;
284a86ed394SVijay Mahadevan 
285a86ed394SVijay Mahadevan       /* Scale ccoords relative to pcoords */
2869566063dSJacob Faibussowitsch       PetscCall(DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size() * tc]));
287a86ed394SVijay Mahadevan     }
2881baa6e33SBarry Smith   } else {
289a86ed394SVijay Mahadevan     factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0);
290a86ed394SVijay Mahadevan   }
291a86ed394SVijay Mahadevan 
292755f3dfbSVijay Mahadevan   /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */
2939566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
294b117cd09SVijay Mahadevan 
295e882eb38SVijay Mahadevan   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
2962417220eSVijay Mahadevan   for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) {
297b117cd09SVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
298b117cd09SVijay Mahadevan 
299b117cd09SVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
300a86ed394SVijay Mahadevan     children.clear();
301a86ed394SVijay Mahadevan     connc.clear();
3029371c9d4SSatish Balay     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);
3039371c9d4SSatish Balay     MBERRNM(merr);
304b117cd09SVijay Mahadevan 
305b117cd09SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
3069371c9d4SSatish Balay     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);
3079371c9d4SSatish Balay     MBERRNM(merr);
3089371c9d4SSatish Balay     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);
3099371c9d4SSatish Balay     MBERRNM(merr);
310b117cd09SVijay Mahadevan 
311a86ed394SVijay Mahadevan     pcoords.resize(connp.size() * 3);
312a86ed394SVijay Mahadevan     ccoords.resize(connc.size() * 3);
313b117cd09SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
3149371c9d4SSatish Balay     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);
3159371c9d4SSatish Balay     MBERRNM(merr);
3169371c9d4SSatish Balay     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);
3179371c9d4SSatish Balay     MBERRNM(merr);
318b117cd09SVijay Mahadevan 
319b117cd09SVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
320b117cd09SVijay Mahadevan     /* TODO: specific to scalar system - use GetDofs */
3219566063dSJacob Faibussowitsch     PetscCall(DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0]));
3229566063dSJacob Faibussowitsch     PetscCall(DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0]));
323b117cd09SVijay Mahadevan 
324a86ed394SVijay Mahadevan     /* Compute the actual interpolation weights when projecting solution/residual between levels */
325a86ed394SVijay Mahadevan     if (use_consistent_bases) {
326a86ed394SVijay Mahadevan       /* Use the cached values of natural parameteric coordinates and basis pre-evaluated.
327a86ed394SVijay Mahadevan          We are making an assumption here that UMR used in GMG to generate the hierarchy uses
328a86ed394SVijay Mahadevan          the same template for all elements; This will fail for mixed element meshes (TRI/QUAD).
329a86ed394SVijay Mahadevan 
3302417220eSVijay Mahadevan          TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes)
331a86ed394SVijay Mahadevan       */
332a86ed394SVijay Mahadevan 
333a86ed394SVijay Mahadevan       /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
334a86ed394SVijay Mahadevan       for (unsigned tc = 0; tc < connc.size(); tc++) {
335a86ed394SVijay Mahadevan         /* TODO: Check if we should be using INSERT_VALUES instead */
3369566063dSJacob Faibussowitsch         PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[connp.size() * tc], ADD_VALUES));
337a86ed394SVijay Mahadevan       }
3381baa6e33SBarry Smith     } else {
339e882eb38SVijay Mahadevan       /* Compute the interpolation weights by determining distance of 1-ring
340755f3dfbSVijay Mahadevan          neighbor vertices from current vertex
341755f3dfbSVijay Mahadevan 
342a86ed394SVijay Mahadevan          This should be used only when FEM basis is not used for the discretization.
343a86ed394SVijay Mahadevan          Else, the consistent interface to compute the basis function for interpolation
344a86ed394SVijay Mahadevan          between the levels should be evaluated correctly to preserve convergence of GMG.
3452417220eSVijay Mahadevan          Shephard's basis will be terrible for any unsmooth problems.
346755f3dfbSVijay Mahadevan       */
347755f3dfbSVijay Mahadevan       values_phi.resize(connp.size());
348755f3dfbSVijay Mahadevan       for (unsigned tc = 0; tc < connc.size(); tc++) {
349a86ed394SVijay Mahadevan         PetscReal normsum = 0.0;
350755f3dfbSVijay Mahadevan         for (unsigned tp = 0; tp < connp.size(); tp++) {
351755f3dfbSVijay Mahadevan           values_phi[tp] = 0.0;
3529371c9d4SSatish Balay           for (unsigned k = 0; k < 3; k++) values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim);
353755f3dfbSVijay Mahadevan           if (values_phi[tp] < 1e-12) {
354755f3dfbSVijay Mahadevan             values_phi[tp] = 1e12;
3551baa6e33SBarry Smith           } else {
356755f3dfbSVijay Mahadevan             //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim);
357755f3dfbSVijay Mahadevan             values_phi[tp] = std::pow(values_phi[tp], -1.0);
358755f3dfbSVijay Mahadevan             normsum += values_phi[tp];
359b117cd09SVijay Mahadevan           }
360b117cd09SVijay Mahadevan         }
361755f3dfbSVijay Mahadevan         for (unsigned tp = 0; tp < connp.size(); tp++) {
3629371c9d4SSatish Balay           if (values_phi[tp] > 1e11) values_phi[tp] = factor * 0.5 / connp.size();
3639371c9d4SSatish Balay           else values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum);
364b117cd09SVijay Mahadevan         }
3659566063dSJacob Faibussowitsch         PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES));
366b117cd09SVijay Mahadevan       }
367a86ed394SVijay Mahadevan     }
368b117cd09SVijay Mahadevan   }
369304006b3SVijay Mahadevan   if (vec) *vec = NULL;
3709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY));
3719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY));
372b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
373b117cd09SVijay Mahadevan }
374b117cd09SVijay Mahadevan 
375cab5ea25SPierre Jolivet /*@C
376b117cd09SVijay Mahadevan   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
377b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
378b117cd09SVijay Mahadevan   provided by the user.
379b117cd09SVijay Mahadevan 
380d083f849SBarry Smith   Collective
381b117cd09SVijay Mahadevan 
382b117cd09SVijay Mahadevan   Input Parameter:
383b117cd09SVijay Mahadevan . dmb  - The DMMoab object
384b117cd09SVijay Mahadevan 
385d8d19677SJose E. Roman   Output Parameters:
386a2b725a8SWilliam Gropp + nlevels   - The number of levels of refinement needed to generate the hierarchy
387a2b725a8SWilliam Gropp - ldegrees  - The degree of refinement at each level in the hierarchy
388b117cd09SVijay Mahadevan 
389b117cd09SVijay Mahadevan   Level: beginner
390b117cd09SVijay Mahadevan 
391b117cd09SVijay Mahadevan @*/
3929371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter *ctx) {
393e882eb38SVijay Mahadevan   //DM_Moab        *dmmoab;
394b117cd09SVijay Mahadevan 
395b117cd09SVijay Mahadevan   PetscFunctionBegin;
396b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1, DM_CLASSID, 1);
397b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2, DM_CLASSID, 2);
398e882eb38SVijay Mahadevan   //dmmoab = (DM_Moab*)(dm1)->data;
399b117cd09SVijay Mahadevan 
400b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
401b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
402b117cd09SVijay Mahadevan }
403b117cd09SVijay Mahadevan 
4049371c9d4SSatish Balay static PetscErrorCode DMMoab_UMR_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref) {
405e882eb38SVijay Mahadevan   PetscInt        i, dim;
406b117cd09SVijay Mahadevan   DM              dm2;
407b117cd09SVijay Mahadevan   moab::ErrorCode merr;
408b117cd09SVijay Mahadevan   DM_Moab        *dmb = (DM_Moab *)dm->data, *dd2;
409b117cd09SVijay Mahadevan 
410b117cd09SVijay Mahadevan   PetscFunctionBegin;
411b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
412b117cd09SVijay Mahadevan   PetscValidPointer(dmref, 4);
413b117cd09SVijay Mahadevan 
414b117cd09SVijay Mahadevan   if ((dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine)) {
41563a3b9bcSJacob Faibussowitsch     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);
41663a3b9bcSJacob Faibussowitsch     if (dmb->hlevel - 1 < 0 && !refine) PetscInfo(NULL, "Invalid multigrid coarsen hierarchy level specified (%" PetscInt_FMT "). Creating a NULL object.\n", dmb->hlevel - 1);
417e882eb38SVijay Mahadevan     *dmref = PETSC_NULL;
418b117cd09SVijay Mahadevan     PetscFunctionReturn(0);
419b117cd09SVijay Mahadevan   }
420b117cd09SVijay Mahadevan 
4219566063dSJacob Faibussowitsch   PetscCall(DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2));
422b117cd09SVijay Mahadevan   dd2 = (DM_Moab *)dm2->data;
423b117cd09SVijay Mahadevan 
424b117cd09SVijay Mahadevan   dd2->mbiface = dmb->mbiface;
4259daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
426b117cd09SVijay Mahadevan   dd2->pcomm = dmb->pcomm;
4279daf19fdSVijay Mahadevan #endif
428b117cd09SVijay Mahadevan   dd2->icreatedinstance = PETSC_FALSE;
42964e1c140SVijay Mahadevan   dd2->nghostrings      = dmb->nghostrings;
430b117cd09SVijay Mahadevan 
431b117cd09SVijay Mahadevan   /* set the new level based on refinement/coarsening */
432b117cd09SVijay Mahadevan   if (refine) {
433b117cd09SVijay Mahadevan     dd2->hlevel = dmb->hlevel + 1;
4341baa6e33SBarry Smith   } else {
435b117cd09SVijay Mahadevan     dd2->hlevel = dmb->hlevel - 1;
436b117cd09SVijay Mahadevan   }
437b117cd09SVijay Mahadevan 
438b117cd09SVijay Mahadevan   /* Copy the multilevel hierarchy pointers in MOAB */
439b117cd09SVijay Mahadevan   dd2->hierarchy = dmb->hierarchy;
440b117cd09SVijay Mahadevan   dd2->nhlevels  = dmb->nhlevels;
4419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets));
4429371c9d4SSatish Balay   for (i = 0; i <= dd2->nhlevels; i++) { dd2->hsets[i] = dmb->hsets[i]; }
443b117cd09SVijay Mahadevan   dd2->fileset = dd2->hsets[dd2->hlevel];
444b117cd09SVijay Mahadevan 
445b117cd09SVijay Mahadevan   /* do the remaining initializations for DMMoab */
446b117cd09SVijay Mahadevan   dd2->bs                = dmb->bs;
447b117cd09SVijay Mahadevan   dd2->numFields         = dmb->numFields;
448b117cd09SVijay Mahadevan   dd2->rw_dbglevel       = dmb->rw_dbglevel;
449b117cd09SVijay Mahadevan   dd2->partition_by_rank = dmb->partition_by_rank;
4509566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options));
4519566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options));
452b117cd09SVijay Mahadevan   dd2->read_mode  = dmb->read_mode;
453b117cd09SVijay Mahadevan   dd2->write_mode = dmb->write_mode;
454b117cd09SVijay Mahadevan 
455b117cd09SVijay Mahadevan   /* set global ID tag handle */
4569566063dSJacob Faibussowitsch   PetscCall(DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag));
457b117cd09SVijay Mahadevan 
4589371c9d4SSatish Balay   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);
4599371c9d4SSatish Balay   MBERRNM(merr);
460b117cd09SVijay Mahadevan 
4619566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix));
4629566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4639566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm2, dim));
464b117cd09SVijay Mahadevan 
465b117cd09SVijay Mahadevan   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
466b117cd09SVijay Mahadevan   dm2->ops->creatematrix = dm->ops->creatematrix;
467b117cd09SVijay Mahadevan 
468b117cd09SVijay Mahadevan   /* copy fill information if given */
4699566063dSJacob Faibussowitsch   PetscCall(DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill));
470b117cd09SVijay Mahadevan 
471b117cd09SVijay Mahadevan   /* copy vector type information */
4729566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(dm2, dm->mattype));
4739566063dSJacob Faibussowitsch   PetscCall(DMSetVecType(dm2, dm->vectype));
4743f1c6e43SVijay Mahadevan   dd2->numFields = dmb->numFields;
475*48a46eb9SPierre Jolivet   if (dmb->numFields) PetscCall(DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames));
476b117cd09SVijay Mahadevan 
4779566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm2));
478b117cd09SVijay Mahadevan 
479b117cd09SVijay Mahadevan   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
4809566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm2));
481b117cd09SVijay Mahadevan 
482b117cd09SVijay Mahadevan   *dmref = dm2;
483b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
484b117cd09SVijay Mahadevan }
485b117cd09SVijay Mahadevan 
486cab5ea25SPierre Jolivet /*@C
487b117cd09SVijay Mahadevan   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
488b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
489b117cd09SVijay Mahadevan   provided by the user.
490b117cd09SVijay Mahadevan 
491d083f849SBarry Smith   Collective on dm
492b117cd09SVijay Mahadevan 
493d8d19677SJose E. Roman   Input Parameters:
494e882eb38SVijay Mahadevan + dm  - The DMMoab object
495e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
496b117cd09SVijay Mahadevan 
497b117cd09SVijay Mahadevan   Output Parameter:
498e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL
499b117cd09SVijay Mahadevan 
500e882eb38SVijay Mahadevan   Note: If no refinement was done, the return value is NULL
501e882eb38SVijay Mahadevan 
502e882eb38SVijay Mahadevan   Level: developer
503b117cd09SVijay Mahadevan 
504b117cd09SVijay Mahadevan @*/
5059371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM *dmf) {
506b117cd09SVijay Mahadevan   PetscFunctionBegin;
507b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
508b117cd09SVijay Mahadevan 
5099566063dSJacob Faibussowitsch   PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_TRUE, dmf));
510b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
511b117cd09SVijay Mahadevan }
512b117cd09SVijay Mahadevan 
513cab5ea25SPierre Jolivet /*@C
514b117cd09SVijay Mahadevan   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
515b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
516b117cd09SVijay Mahadevan   provided by the user.
517b117cd09SVijay Mahadevan 
518d083f849SBarry Smith   Collective on dm
519b117cd09SVijay Mahadevan 
520d8d19677SJose E. Roman   Input Parameters:
521a2b725a8SWilliam Gropp + dm  - The DMMoab object
522e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
523b117cd09SVijay Mahadevan 
524b117cd09SVijay Mahadevan   Output Parameter:
525e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL
526b117cd09SVijay Mahadevan 
527e882eb38SVijay Mahadevan   Note: If no coarsening was done, the return value is NULL
528b117cd09SVijay Mahadevan 
529e882eb38SVijay Mahadevan   Level: developer
530e882eb38SVijay Mahadevan 
531b117cd09SVijay Mahadevan @*/
5329371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM *dmc) {
533b117cd09SVijay Mahadevan   PetscFunctionBegin;
534b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5359566063dSJacob Faibussowitsch   PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_FALSE, dmc));
536b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
537b117cd09SVijay Mahadevan }
538