xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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 
16*d8d19677SJose 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 @*/
23b117cd09SVijay Mahadevan PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees)
24b117cd09SVijay Mahadevan {
25b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
26b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
27b117cd09SVijay Mahadevan   moab::ErrorCode merr;
282417220eSVijay Mahadevan   PetscInt *pdegrees, ilevel;
29e882eb38SVijay Mahadevan   std::vector<moab::EntityHandle> hsets;
30b117cd09SVijay Mahadevan 
31b117cd09SVijay Mahadevan   PetscFunctionBegin;
32b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
33b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
34b117cd09SVijay Mahadevan 
35b117cd09SVijay Mahadevan   if (!ldegrees) {
36b117cd09SVijay Mahadevan     ierr = PetscMalloc1(nlevels, &pdegrees);CHKERRQ(ierr);
372417220eSVijay Mahadevan     for (ilevel = 0; ilevel < nlevels; ilevel++) pdegrees[ilevel] = 2; /* default = Degree 2 refinement */
38b117cd09SVijay Mahadevan   }
39b117cd09SVijay Mahadevan   else pdegrees = ldegrees;
40b117cd09SVijay Mahadevan 
41b117cd09SVijay Mahadevan   /* initialize set level refinement data for hierarchy */
42b117cd09SVijay Mahadevan   dmmoab->nhlevels = nlevels;
43b117cd09SVijay Mahadevan 
44b117cd09SVijay Mahadevan   /* Instantiate the nested refinement class */
459daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
463f1c6e43SVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset);
479daf19fdSVijay Mahadevan #else
489daf19fdSVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), NULL, dmmoab->fileset);
499daf19fdSVijay Mahadevan #endif
50b117cd09SVijay Mahadevan 
51e882eb38SVijay Mahadevan   ierr = PetscMalloc1(nlevels + 1, &dmmoab->hsets);CHKERRQ(ierr);
52e882eb38SVijay Mahadevan 
53b117cd09SVijay Mahadevan   /* generate the mesh hierarchy */
549c368985SVijay Mahadevan   merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false);MBERRNM(merr);
55e882eb38SVijay Mahadevan 
569daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
57755f3dfbSVijay Mahadevan   if (dmmoab->pcomm->size() > 1) {
5864e1c140SVijay Mahadevan     merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);MBERRNM(merr);
59755f3dfbSVijay Mahadevan   }
609daf19fdSVijay Mahadevan #endif
6149d66b22SVijay Mahadevan 
6264e1c140SVijay Mahadevan   /* copy the mesh sets for nested refinement hierarchy */
63e92d1c7cSVijay Mahadevan   dmmoab->hsets[0] = hsets[0];
64e92d1c7cSVijay Mahadevan   for (ilevel = 1; ilevel <= nlevels; ilevel++)
652417220eSVijay Mahadevan   {
662417220eSVijay Mahadevan     dmmoab->hsets[ilevel] = hsets[ilevel];
67e92d1c7cSVijay Mahadevan 
689c368985SVijay Mahadevan #ifdef MOAB_HAVE_MPI
69e92d1c7cSVijay Mahadevan     merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false);MBERRNM(merr);
709c368985SVijay Mahadevan #endif
71e92d1c7cSVijay Mahadevan 
72e92d1c7cSVijay Mahadevan     /* Update material and other geometric tags from parent to child sets */
73e92d1c7cSVijay Mahadevan     merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]);MBERRNM(merr);
742417220eSVijay Mahadevan   }
7564e1c140SVijay Mahadevan 
7664e1c140SVijay Mahadevan   hsets.clear();
77b117cd09SVijay Mahadevan   if (!ldegrees) {
78b117cd09SVijay Mahadevan     ierr = PetscFree(pdegrees);CHKERRQ(ierr);
79b117cd09SVijay Mahadevan   }
80b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
81b117cd09SVijay Mahadevan }
82b117cd09SVijay Mahadevan 
83cab5ea25SPierre Jolivet /*@C
84e882eb38SVijay Mahadevan   DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy
85e882eb38SVijay Mahadevan   by succesively refining a coarse mesh.
86b117cd09SVijay Mahadevan 
87d083f849SBarry Smith   Collective
88b117cd09SVijay Mahadevan 
89b117cd09SVijay Mahadevan   Input Parameter:
90a2b725a8SWilliam Gropp . dm  - The DMMoab object
91b117cd09SVijay Mahadevan 
92*d8d19677SJose E. Roman   Output Parameters:
93e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
94a2b725a8SWilliam Gropp - dmf  - The DM objects after successive refinement of the hierarchy
95b117cd09SVijay Mahadevan 
96b117cd09SVijay Mahadevan   Level: beginner
97b117cd09SVijay Mahadevan 
98b117cd09SVijay Mahadevan @*/
99304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode  DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[])
100b117cd09SVijay Mahadevan {
101b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
102b117cd09SVijay Mahadevan   PetscInt        i;
103b117cd09SVijay Mahadevan 
104b117cd09SVijay Mahadevan   PetscFunctionBegin;
105b117cd09SVijay Mahadevan 
106e882eb38SVijay Mahadevan   ierr = DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]);CHKERRQ(ierr);
107e882eb38SVijay Mahadevan   for (i = 1; i < nlevels; i++) {
108e882eb38SVijay Mahadevan     ierr = DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]);CHKERRQ(ierr);
109b117cd09SVijay Mahadevan   }
110b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
111b117cd09SVijay Mahadevan }
112b117cd09SVijay Mahadevan 
113cab5ea25SPierre Jolivet /*@C
114e882eb38SVijay Mahadevan   DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy
115e882eb38SVijay Mahadevan   by succesively coarsening a refined mesh.
116b117cd09SVijay Mahadevan 
117d083f849SBarry Smith   Collective
118b117cd09SVijay Mahadevan 
119b117cd09SVijay Mahadevan   Input Parameter:
120a2b725a8SWilliam Gropp . dm  - The DMMoab object
121b117cd09SVijay Mahadevan 
122*d8d19677SJose E. Roman   Output Parameters:
123e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
124a2b725a8SWilliam Gropp - dmc  - The DM objects after successive coarsening of the hierarchy
125b117cd09SVijay Mahadevan 
126b117cd09SVijay Mahadevan   Level: beginner
127b117cd09SVijay Mahadevan 
128b117cd09SVijay Mahadevan @*/
129304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[])
130b117cd09SVijay Mahadevan {
131b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
132b117cd09SVijay Mahadevan   PetscInt        i;
133b117cd09SVijay Mahadevan 
134b117cd09SVijay Mahadevan   PetscFunctionBegin;
135b117cd09SVijay Mahadevan 
136e882eb38SVijay Mahadevan   ierr = DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]);CHKERRQ(ierr);
137e882eb38SVijay Mahadevan   for (i = 1; i < nlevels; i++) {
138e882eb38SVijay Mahadevan     ierr = DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]);CHKERRQ(ierr);
139b117cd09SVijay Mahadevan   }
140b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
141b117cd09SVijay Mahadevan }
142b117cd09SVijay Mahadevan 
1432417220eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt*, PetscInt*, PetscInt*, PetscInt*, PetscBool);
144b117cd09SVijay Mahadevan 
145cab5ea25SPierre Jolivet /*@C
146e882eb38SVijay Mahadevan   DMCreateInterpolation_Moab - Generate the interpolation operators to transform
147e882eb38SVijay Mahadevan   operators (matrices, vectors) from parent level to child level as defined by
148e882eb38SVijay Mahadevan   the DM inputs provided by the user.
149b117cd09SVijay Mahadevan 
150d083f849SBarry Smith   Collective
151b117cd09SVijay Mahadevan 
152*d8d19677SJose E. Roman   Input Parameters:
153e882eb38SVijay Mahadevan + dm1  - The DMMoab object
154e882eb38SVijay Mahadevan - dm2  - the second, finer DMMoab object
155b117cd09SVijay Mahadevan 
156*d8d19677SJose E. Roman   Output Parameters:
157e882eb38SVijay Mahadevan + interpl  - The interpolation operator for transferring data between the levels
158e882eb38SVijay Mahadevan - vec      - The scaling vector (optional)
159b117cd09SVijay Mahadevan 
160e882eb38SVijay Mahadevan   Level: developer
161b117cd09SVijay Mahadevan 
162b117cd09SVijay Mahadevan @*/
163304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat* interpl, Vec* vec)
164b117cd09SVijay Mahadevan {
165755f3dfbSVijay Mahadevan   DM_Moab         *dmbp, *dmbc;
166b117cd09SVijay Mahadevan   PetscErrorCode   ierr;
167b117cd09SVijay Mahadevan   moab::ErrorCode  merr;
168e882eb38SVijay Mahadevan   PetscInt         dim;
1693f1c6e43SVijay Mahadevan   PetscReal        factor;
170ce27a4eeSVijay Mahadevan   PetscInt         innz, *nnz, ionz, *onz;
171755f3dfbSVijay Mahadevan   PetscInt         nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc;
172a86ed394SVijay Mahadevan   const PetscBool  use_consistent_bases=PETSC_TRUE;
173b117cd09SVijay Mahadevan 
174b117cd09SVijay Mahadevan   PetscFunctionBegin;
175755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmp, DM_CLASSID, 1);
176755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmc, DM_CLASSID, 2);
177755f3dfbSVijay Mahadevan   dmbp = (DM_Moab*)(dmp)->data;
178755f3dfbSVijay Mahadevan   dmbc = (DM_Moab*)(dmc)->data;
179755f3dfbSVijay Mahadevan   nlsizp = dmbp->nloc;// *dmb1->numFields;
180755f3dfbSVijay Mahadevan   nlsizc = dmbc->nloc;// *dmb2->numFields;
181755f3dfbSVijay Mahadevan   ngsizp = dmbp->n; // *dmb1->numFields;
182755f3dfbSVijay Mahadevan   ngsizc = dmbc->n; // *dmb2->numFields;
183755f3dfbSVijay Mahadevan   nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields;
184b117cd09SVijay Mahadevan 
1852417220eSVijay Mahadevan   // Columns = Parent DoFs ;  Rows = Child DoFs
186755f3dfbSVijay Mahadevan   // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent))
187755f3dfbSVijay Mahadevan   // Size: nlsizc * nlghsizp
188755f3dfbSVijay Mahadevan   PetscInfo4(NULL, "Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n", ngsizc, nlghsizp, dmbp->hlevel, dmbc->hlevel);
189b117cd09SVijay Mahadevan 
190a86ed394SVijay Mahadevan   ierr = DMGetDimension(dmp, &dim);CHKERRQ(ierr);
191a86ed394SVijay Mahadevan 
192941e0cffSVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
193755f3dfbSVijay Mahadevan   ierr = PetscCalloc2(nlsizc, &nnz, nlsizc, &onz);CHKERRQ(ierr);
194941e0cffSVijay Mahadevan 
195941e0cffSVijay Mahadevan   /* Loop through the local elements and compute the relation between the current parent and the refined_level. */
1962417220eSVijay Mahadevan   for (moab::Range::iterator iter = dmbc->vowned->begin(); iter != dmbc->vowned->end(); iter++) {
197941e0cffSVijay Mahadevan 
1982417220eSVijay Mahadevan     const moab::EntityHandle vhandle = *iter;
1992417220eSVijay Mahadevan     /* define local variables */
2002417220eSVijay Mahadevan     moab::EntityHandle parent;
2012417220eSVijay Mahadevan     std::vector<moab::EntityHandle> adjs;
2022417220eSVijay Mahadevan     moab::Range     found;
2032417220eSVijay Mahadevan 
2042417220eSVijay Mahadevan     /* store the vertex DoF number */
2052417220eSVijay Mahadevan     const int ldof = dmbc->lidmap[vhandle - dmbc->seqstart];
2062417220eSVijay Mahadevan 
2072417220eSVijay Mahadevan     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
2082417220eSVijay Mahadevan        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
2092417220eSVijay Mahadevan        non-zero pattern accordingly. */
2102417220eSVijay Mahadevan     merr = dmbc->hierarchy->get_adjacencies(vhandle, dmbc->dim, adjs);MBERRNM(merr);
2112417220eSVijay Mahadevan 
2122417220eSVijay Mahadevan     /* loop over vertices and update the number of connectivity */
2132417220eSVijay Mahadevan     for (unsigned jter = 0; jter < adjs.size(); jter++) {
2142417220eSVijay Mahadevan 
2152417220eSVijay Mahadevan       const moab::EntityHandle jhandle = adjs[jter];
216941e0cffSVijay Mahadevan 
217941e0cffSVijay Mahadevan       /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
2182417220eSVijay Mahadevan       merr = dmbc->hierarchy->child_to_parent(jhandle, dmbc->hlevel, dmbp->hlevel, &parent);MBERRNM(merr);
219941e0cffSVijay Mahadevan 
2202417220eSVijay Mahadevan       /* Get connectivity information in canonical ordering for the local element */
2212417220eSVijay Mahadevan       std::vector<moab::EntityHandle> connp;
2222417220eSVijay Mahadevan       merr = dmbp->hierarchy->get_connectivity(parent, dmbp->hlevel, connp);MBERRNM(merr);
223a86ed394SVijay Mahadevan 
2242417220eSVijay Mahadevan       for (unsigned ic=0; ic < connp.size(); ++ic) {
2252417220eSVijay Mahadevan 
2262417220eSVijay Mahadevan         /* loop over each element connected to the adjacent vertex and update as needed */
2272417220eSVijay Mahadevan         /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
2282417220eSVijay Mahadevan         if (found.find(connp[ic]) != found.end()) continue; /* make sure we don't double count shared vertices */
2292417220eSVijay Mahadevan         if (dmbp->vghost->find(connp[ic]) != dmbp->vghost->end()) onz[ldof]++; /* update out-of-proc onz */
2302417220eSVijay Mahadevan         else nnz[ldof]++; /* else local vertex */
2312417220eSVijay Mahadevan         found.insert(connp[ic]);
2322417220eSVijay Mahadevan       }
2332417220eSVijay Mahadevan     }
234941e0cffSVijay Mahadevan   }
235941e0cffSVijay Mahadevan 
2362417220eSVijay Mahadevan   for (int i = 0; i < nlsizc; i++)
2372417220eSVijay Mahadevan     nnz[i] += 1; /* self count the node */
238941e0cffSVijay Mahadevan 
239941e0cffSVijay Mahadevan   ionz = onz[0];
240941e0cffSVijay Mahadevan   innz = nnz[0];
241755f3dfbSVijay Mahadevan   for (int tc = 0; tc < nlsizc; tc++) {
242941e0cffSVijay Mahadevan     // check for maximum allowed sparsity = fully dense
243755f3dfbSVijay Mahadevan     nnz[tc] = std::min(nlsizp, nnz[tc]);
2442417220eSVijay Mahadevan     onz[tc] = std::min(ngsizp - nlsizp, onz[tc]);
2452417220eSVijay Mahadevan 
2462417220eSVijay Mahadevan     PetscInfo3(NULL, "  %d: NNZ = %d, ONZ = %d\n", tc, nnz[tc], onz[tc]);
247941e0cffSVijay Mahadevan 
248941e0cffSVijay Mahadevan     innz = (innz < nnz[tc] ? nnz[tc] : innz);
249941e0cffSVijay Mahadevan     ionz = (ionz < onz[tc] ? onz[tc] : ionz);
250941e0cffSVijay Mahadevan   }
25164e1c140SVijay Mahadevan 
25264e1c140SVijay Mahadevan   /* create interpolation matrix */
253755f3dfbSVijay Mahadevan   ierr = MatCreate(PetscObjectComm((PetscObject)dmc), interpl);CHKERRQ(ierr);
254755f3dfbSVijay Mahadevan   ierr = MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp);CHKERRQ(ierr);
25564e1c140SVijay Mahadevan   ierr = MatSetType(*interpl, MATAIJ);CHKERRQ(ierr);
25664e1c140SVijay Mahadevan   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
25764e1c140SVijay Mahadevan 
258941e0cffSVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(*interpl, innz, nnz);CHKERRQ(ierr);
259941e0cffSVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz);CHKERRQ(ierr);
260941e0cffSVijay Mahadevan 
261941e0cffSVijay Mahadevan   /* clean up temporary memory */
262941e0cffSVijay Mahadevan   ierr = PetscFree2(nnz, onz);CHKERRQ(ierr);
263b117cd09SVijay Mahadevan 
264b117cd09SVijay Mahadevan   /* set up internal matrix data-structures */
265b117cd09SVijay Mahadevan   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
266b117cd09SVijay Mahadevan 
267a86ed394SVijay Mahadevan   /* Define variables for assembly */
268a86ed394SVijay Mahadevan   std::vector<moab::EntityHandle> children;
269a86ed394SVijay Mahadevan   std::vector<moab::EntityHandle> connp, connc;
270a86ed394SVijay Mahadevan   std::vector<PetscReal> pcoords, ccoords, values_phi;
271a86ed394SVijay Mahadevan 
272a86ed394SVijay Mahadevan   if (use_consistent_bases) {
2732417220eSVijay Mahadevan     const moab::EntityHandle ehandle = dmbp->elocal->front();
274a86ed394SVijay Mahadevan 
275a86ed394SVijay Mahadevan     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
276a86ed394SVijay Mahadevan 
277a86ed394SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
278a86ed394SVijay Mahadevan     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
2792417220eSVijay Mahadevan     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr);
280a86ed394SVijay Mahadevan 
281a86ed394SVijay Mahadevan     std::vector<PetscReal> natparam(3*connc.size(), 0.0);
282a86ed394SVijay Mahadevan     pcoords.resize(connp.size() * 3);
283a86ed394SVijay Mahadevan     ccoords.resize(connc.size() * 3);
284a86ed394SVijay Mahadevan     values_phi.resize(connp.size()*connc.size());
285a86ed394SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
286a86ed394SVijay Mahadevan     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr);
287a86ed394SVijay Mahadevan     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr);
288a86ed394SVijay Mahadevan 
289a86ed394SVijay Mahadevan     /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
290a86ed394SVijay Mahadevan     for (unsigned tc = 0; tc < connc.size(); tc++) {
291a86ed394SVijay Mahadevan       const PetscInt offset = tc * 3;
292a86ed394SVijay Mahadevan 
293a86ed394SVijay Mahadevan       /* Scale ccoords relative to pcoords */
294a86ed394SVijay Mahadevan       ierr = DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size()*tc]);CHKERRQ(ierr);
295a86ed394SVijay Mahadevan     }
296a86ed394SVijay Mahadevan   }
297a86ed394SVijay Mahadevan   else {
298a86ed394SVijay Mahadevan     factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0);
299a86ed394SVijay Mahadevan   }
300a86ed394SVijay Mahadevan 
301755f3dfbSVijay Mahadevan   /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */
3022417220eSVijay Mahadevan   ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
303b117cd09SVijay Mahadevan 
304e882eb38SVijay Mahadevan   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
3052417220eSVijay Mahadevan   for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) {
306b117cd09SVijay Mahadevan 
307b117cd09SVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
308b117cd09SVijay Mahadevan 
309b117cd09SVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
310a86ed394SVijay Mahadevan     children.clear();
311a86ed394SVijay Mahadevan     connc.clear();
312755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
313b117cd09SVijay Mahadevan 
314b117cd09SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
315755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
3162417220eSVijay Mahadevan     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr);
317b117cd09SVijay Mahadevan 
318a86ed394SVijay Mahadevan     pcoords.resize(connp.size() * 3);
319a86ed394SVijay Mahadevan     ccoords.resize(connc.size() * 3);
320b117cd09SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
321755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr);
322755f3dfbSVijay Mahadevan     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr);
323b117cd09SVijay Mahadevan 
324b117cd09SVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
325b117cd09SVijay Mahadevan     /* TODO: specific to scalar system - use GetDofs */
326755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr);
327755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr);
328b117cd09SVijay Mahadevan 
329a86ed394SVijay Mahadevan     /* Compute the actual interpolation weights when projecting solution/residual between levels */
330a86ed394SVijay Mahadevan     if (use_consistent_bases) {
331a86ed394SVijay Mahadevan 
332a86ed394SVijay Mahadevan       /* Use the cached values of natural parameteric coordinates and basis pre-evaluated.
333a86ed394SVijay Mahadevan          We are making an assumption here that UMR used in GMG to generate the hierarchy uses
334a86ed394SVijay Mahadevan          the same template for all elements; This will fail for mixed element meshes (TRI/QUAD).
335a86ed394SVijay Mahadevan 
3362417220eSVijay Mahadevan          TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes)
337a86ed394SVijay Mahadevan       */
338a86ed394SVijay Mahadevan 
339a86ed394SVijay Mahadevan       /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
340a86ed394SVijay Mahadevan       for (unsigned tc = 0; tc < connc.size(); tc++) {
341a86ed394SVijay Mahadevan         /* TODO: Check if we should be using INSERT_VALUES instead */
342a86ed394SVijay Mahadevan         ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[connp.size()*tc], ADD_VALUES);CHKERRQ(ierr);
343a86ed394SVijay Mahadevan       }
344a86ed394SVijay Mahadevan     }
345a86ed394SVijay Mahadevan     else {
346e882eb38SVijay Mahadevan       /* Compute the interpolation weights by determining distance of 1-ring
347755f3dfbSVijay Mahadevan          neighbor vertices from current vertex
348755f3dfbSVijay Mahadevan 
349a86ed394SVijay Mahadevan          This should be used only when FEM basis is not used for the discretization.
350a86ed394SVijay Mahadevan          Else, the consistent interface to compute the basis function for interpolation
351a86ed394SVijay Mahadevan          between the levels should be evaluated correctly to preserve convergence of GMG.
3522417220eSVijay Mahadevan          Shephard's basis will be terrible for any unsmooth problems.
353755f3dfbSVijay Mahadevan       */
354755f3dfbSVijay Mahadevan       values_phi.resize(connp.size());
355755f3dfbSVijay Mahadevan       for (unsigned tc = 0; tc < connc.size(); tc++) {
356755f3dfbSVijay Mahadevan 
357a86ed394SVijay Mahadevan         PetscReal normsum = 0.0;
358755f3dfbSVijay Mahadevan         for (unsigned tp = 0; tp < connp.size(); tp++) {
359755f3dfbSVijay Mahadevan           values_phi[tp] = 0.0;
360e882eb38SVijay Mahadevan           for (unsigned k = 0; k < 3; k++)
361755f3dfbSVijay Mahadevan             values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim);
362755f3dfbSVijay Mahadevan           if (values_phi[tp] < 1e-12) {
363755f3dfbSVijay Mahadevan             values_phi[tp] = 1e12;
364b117cd09SVijay Mahadevan           }
365b117cd09SVijay Mahadevan           else {
366755f3dfbSVijay Mahadevan             //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim);
367755f3dfbSVijay Mahadevan             values_phi[tp] = std::pow(values_phi[tp], -1.0);
368755f3dfbSVijay Mahadevan             normsum += values_phi[tp];
369b117cd09SVijay Mahadevan           }
370b117cd09SVijay Mahadevan         }
371755f3dfbSVijay Mahadevan         for (unsigned tp = 0; tp < connp.size(); tp++) {
372755f3dfbSVijay Mahadevan           if (values_phi[tp] > 1e11)
373755f3dfbSVijay Mahadevan             values_phi[tp] = factor * 0.5 / connp.size();
374b117cd09SVijay Mahadevan           else
375755f3dfbSVijay Mahadevan             values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum);
376b117cd09SVijay Mahadevan         }
377755f3dfbSVijay Mahadevan         ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
378b117cd09SVijay Mahadevan       }
379a86ed394SVijay Mahadevan     }
380b117cd09SVijay Mahadevan   }
381304006b3SVijay Mahadevan   if (vec) *vec = NULL;
382b117cd09SVijay Mahadevan   ierr = MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
383b117cd09SVijay Mahadevan   ierr = MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
384b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
385b117cd09SVijay Mahadevan }
386b117cd09SVijay Mahadevan 
387cab5ea25SPierre Jolivet /*@C
388b117cd09SVijay Mahadevan   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
389b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
390b117cd09SVijay Mahadevan   provided by the user.
391b117cd09SVijay Mahadevan 
392d083f849SBarry Smith   Collective
393b117cd09SVijay Mahadevan 
394b117cd09SVijay Mahadevan   Input Parameter:
395b117cd09SVijay Mahadevan . dmb  - The DMMoab object
396b117cd09SVijay Mahadevan 
397*d8d19677SJose E. Roman   Output Parameters:
398a2b725a8SWilliam Gropp + nlevels   - The number of levels of refinement needed to generate the hierarchy
399a2b725a8SWilliam Gropp - ldegrees  - The degree of refinement at each level in the hierarchy
400b117cd09SVijay Mahadevan 
401b117cd09SVijay Mahadevan   Level: beginner
402b117cd09SVijay Mahadevan 
403b117cd09SVijay Mahadevan @*/
404304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter* ctx)
405b117cd09SVijay Mahadevan {
406e882eb38SVijay Mahadevan   //DM_Moab        *dmmoab;
407b117cd09SVijay Mahadevan 
408b117cd09SVijay Mahadevan   PetscFunctionBegin;
409b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1, DM_CLASSID, 1);
410b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2, DM_CLASSID, 2);
411e882eb38SVijay Mahadevan   //dmmoab = (DM_Moab*)(dm1)->data;
412b117cd09SVijay Mahadevan 
413b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
414b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
415b117cd09SVijay Mahadevan }
416b117cd09SVijay Mahadevan 
417470880abSPatrick Sanan static PetscErrorCode DMMoab_UMR_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref)
418b117cd09SVijay Mahadevan {
419b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
420e882eb38SVijay Mahadevan   PetscInt        i, dim;
421b117cd09SVijay Mahadevan   DM              dm2;
422b117cd09SVijay Mahadevan   moab::ErrorCode merr;
423b117cd09SVijay Mahadevan   DM_Moab        *dmb = (DM_Moab*)dm->data, *dd2;
424b117cd09SVijay Mahadevan 
425b117cd09SVijay Mahadevan   PetscFunctionBegin;
426b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
427b117cd09SVijay Mahadevan   PetscValidPointer(dmref, 4);
428b117cd09SVijay Mahadevan 
429b117cd09SVijay Mahadevan   if ((dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine)) {
430e882eb38SVijay Mahadevan     if (dmb->hlevel + 1 > dmb->nhlevels && refine) PetscInfo2(NULL, "Invalid multigrid refinement hierarchy level specified (%D). MOAB UMR max levels = %D. Creating a NULL object.\n", dmb->hlevel + 1, dmb->nhlevels);
431e882eb38SVijay Mahadevan     if (dmb->hlevel - 1 < 0 && !refine) PetscInfo1(NULL, "Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n", dmb->hlevel - 1);
432e882eb38SVijay Mahadevan     *dmref = PETSC_NULL;
433b117cd09SVijay Mahadevan     PetscFunctionReturn(0);
434b117cd09SVijay Mahadevan   }
435b117cd09SVijay Mahadevan 
436b117cd09SVijay Mahadevan   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
437b117cd09SVijay Mahadevan   dd2 = (DM_Moab*)dm2->data;
438b117cd09SVijay Mahadevan 
439b117cd09SVijay Mahadevan   dd2->mbiface = dmb->mbiface;
4409daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
441b117cd09SVijay Mahadevan   dd2->pcomm = dmb->pcomm;
4429daf19fdSVijay Mahadevan #endif
443b117cd09SVijay Mahadevan   dd2->icreatedinstance = PETSC_FALSE;
44464e1c140SVijay Mahadevan   dd2->nghostrings = dmb->nghostrings;
445b117cd09SVijay Mahadevan 
446b117cd09SVijay Mahadevan   /* set the new level based on refinement/coarsening */
447b117cd09SVijay Mahadevan   if (refine) {
448b117cd09SVijay Mahadevan     dd2->hlevel = dmb->hlevel + 1;
449b117cd09SVijay Mahadevan   }
450b117cd09SVijay Mahadevan   else {
451b117cd09SVijay Mahadevan     dd2->hlevel = dmb->hlevel - 1;
452b117cd09SVijay Mahadevan   }
453b117cd09SVijay Mahadevan 
454b117cd09SVijay Mahadevan   /* Copy the multilevel hierarchy pointers in MOAB */
455b117cd09SVijay Mahadevan   dd2->hierarchy = dmb->hierarchy;
456b117cd09SVijay Mahadevan   dd2->nhlevels = dmb->nhlevels;
457b117cd09SVijay Mahadevan   ierr = PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets);CHKERRQ(ierr);
458b117cd09SVijay Mahadevan   for (i = 0; i <= dd2->nhlevels; i++) {
459b117cd09SVijay Mahadevan     dd2->hsets[i] = dmb->hsets[i];
460b117cd09SVijay Mahadevan   }
461b117cd09SVijay Mahadevan   dd2->fileset = dd2->hsets[dd2->hlevel];
462b117cd09SVijay Mahadevan 
463b117cd09SVijay Mahadevan   /* do the remaining initializations for DMMoab */
464b117cd09SVijay Mahadevan   dd2->bs = dmb->bs;
465b117cd09SVijay Mahadevan   dd2->numFields = dmb->numFields;
466b117cd09SVijay Mahadevan   dd2->rw_dbglevel = dmb->rw_dbglevel;
467b117cd09SVijay Mahadevan   dd2->partition_by_rank = dmb->partition_by_rank;
468b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
469b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
470b117cd09SVijay Mahadevan   dd2->read_mode = dmb->read_mode;
471b117cd09SVijay Mahadevan   dd2->write_mode = dmb->write_mode;
472b117cd09SVijay Mahadevan 
473b117cd09SVijay Mahadevan   /* set global ID tag handle */
474b117cd09SVijay Mahadevan   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
475b117cd09SVijay Mahadevan 
476b117cd09SVijay Mahadevan   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
477b117cd09SVijay Mahadevan 
478b117cd09SVijay Mahadevan   ierr = DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix);CHKERRQ(ierr);
479b117cd09SVijay Mahadevan   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
480b117cd09SVijay Mahadevan   ierr = DMSetDimension(dm2, dim);CHKERRQ(ierr);
481b117cd09SVijay Mahadevan 
482b117cd09SVijay Mahadevan   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
483b117cd09SVijay Mahadevan   dm2->ops->creatematrix = dm->ops->creatematrix;
484b117cd09SVijay Mahadevan 
485b117cd09SVijay Mahadevan   /* copy fill information if given */
486b117cd09SVijay Mahadevan   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
487b117cd09SVijay Mahadevan 
488b117cd09SVijay Mahadevan   /* copy vector type information */
489b117cd09SVijay Mahadevan   ierr = DMSetMatType(dm2, dm->mattype);CHKERRQ(ierr);
490b117cd09SVijay Mahadevan   ierr = DMSetVecType(dm2, dm->vectype);CHKERRQ(ierr);
4913f1c6e43SVijay Mahadevan   dd2->numFields = dmb->numFields;
4923f1c6e43SVijay Mahadevan   if (dmb->numFields) {
493b117cd09SVijay Mahadevan     ierr = DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames);CHKERRQ(ierr);
4943f1c6e43SVijay Mahadevan   }
495b117cd09SVijay Mahadevan 
496b117cd09SVijay Mahadevan   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
497b117cd09SVijay Mahadevan 
498b117cd09SVijay Mahadevan   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
499b117cd09SVijay Mahadevan   ierr = DMSetUp(dm2);CHKERRQ(ierr);
500b117cd09SVijay Mahadevan 
501b117cd09SVijay Mahadevan   *dmref = dm2;
502b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
503b117cd09SVijay Mahadevan }
504b117cd09SVijay Mahadevan 
505cab5ea25SPierre Jolivet /*@C
506b117cd09SVijay Mahadevan   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
507b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
508b117cd09SVijay Mahadevan   provided by the user.
509b117cd09SVijay Mahadevan 
510d083f849SBarry Smith   Collective on dm
511b117cd09SVijay Mahadevan 
512*d8d19677SJose E. Roman   Input Parameters:
513e882eb38SVijay Mahadevan + dm  - The DMMoab object
514e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
515b117cd09SVijay Mahadevan 
516b117cd09SVijay Mahadevan   Output Parameter:
517e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL
518b117cd09SVijay Mahadevan 
519e882eb38SVijay Mahadevan   Note: If no refinement was done, the return value is NULL
520e882eb38SVijay Mahadevan 
521e882eb38SVijay Mahadevan   Level: developer
522b117cd09SVijay Mahadevan 
523b117cd09SVijay Mahadevan @*/
524304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM* dmf)
525b117cd09SVijay Mahadevan {
526b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
527b117cd09SVijay Mahadevan 
528b117cd09SVijay Mahadevan   PetscFunctionBegin;
529b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
530b117cd09SVijay Mahadevan 
531470880abSPatrick Sanan   ierr = DMMoab_UMR_Private(dm, comm, PETSC_TRUE, dmf);CHKERRQ(ierr);
532b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
533b117cd09SVijay Mahadevan }
534b117cd09SVijay Mahadevan 
535cab5ea25SPierre Jolivet /*@C
536b117cd09SVijay Mahadevan   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
537b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
538b117cd09SVijay Mahadevan   provided by the user.
539b117cd09SVijay Mahadevan 
540d083f849SBarry Smith   Collective on dm
541b117cd09SVijay Mahadevan 
542*d8d19677SJose E. Roman   Input Parameters:
543a2b725a8SWilliam Gropp + dm  - The DMMoab object
544e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
545b117cd09SVijay Mahadevan 
546b117cd09SVijay Mahadevan   Output Parameter:
547e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL
548b117cd09SVijay Mahadevan 
549e882eb38SVijay Mahadevan   Note: If no coarsening was done, the return value is NULL
550b117cd09SVijay Mahadevan 
551e882eb38SVijay Mahadevan   Level: developer
552e882eb38SVijay Mahadevan 
553b117cd09SVijay Mahadevan @*/
554304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM* dmc)
555b117cd09SVijay Mahadevan {
556b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
557b117cd09SVijay Mahadevan 
558b117cd09SVijay Mahadevan   PetscFunctionBegin;
559b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
560b117cd09SVijay Mahadevan 
561470880abSPatrick Sanan   ierr = DMMoab_UMR_Private(dm, comm, PETSC_FALSE, dmc);CHKERRQ(ierr);
562b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
563b117cd09SVijay Mahadevan }
564