xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision a2b725a8db0d6bf6cc2a1c6df7dd8029aadfff6e)
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 
6b117cd09SVijay Mahadevan /*@
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 
11b117cd09SVijay Mahadevan   Collective on MPI_Comm
12b117cd09SVijay Mahadevan 
13b117cd09SVijay Mahadevan   Input Parameter:
14*a2b725a8SWilliam Gropp . dmb  - The DMMoab object
15b117cd09SVijay Mahadevan 
16b117cd09SVijay Mahadevan   Output Parameter:
17b117cd09SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
18*a2b725a8SWilliam Gropp - ldegrees  - The degree of refinement at each level in the hierarchy
19b117cd09SVijay Mahadevan 
20b117cd09SVijay Mahadevan   Level: beginner
21b117cd09SVijay Mahadevan 
22b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
23b117cd09SVijay Mahadevan @*/
24b117cd09SVijay Mahadevan PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees)
25b117cd09SVijay Mahadevan {
26b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
27b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
28b117cd09SVijay Mahadevan   moab::ErrorCode merr;
292417220eSVijay Mahadevan   PetscInt *pdegrees, ilevel;
30e882eb38SVijay Mahadevan   std::vector<moab::EntityHandle> hsets;
31b117cd09SVijay Mahadevan 
32b117cd09SVijay Mahadevan   PetscFunctionBegin;
33b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
34b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
35b117cd09SVijay Mahadevan 
36b117cd09SVijay Mahadevan   if (!ldegrees) {
37b117cd09SVijay Mahadevan     ierr = PetscMalloc1(nlevels, &pdegrees);CHKERRQ(ierr);
382417220eSVijay Mahadevan     for (ilevel = 0; ilevel < nlevels; ilevel++) pdegrees[ilevel] = 2; /* default = Degree 2 refinement */
39b117cd09SVijay Mahadevan   }
40b117cd09SVijay Mahadevan   else pdegrees = ldegrees;
41b117cd09SVijay Mahadevan 
42b117cd09SVijay Mahadevan   /* initialize set level refinement data for hierarchy */
43b117cd09SVijay Mahadevan   dmmoab->nhlevels = nlevels;
44b117cd09SVijay Mahadevan 
45b117cd09SVijay Mahadevan   /* Instantiate the nested refinement class */
469daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
473f1c6e43SVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset);
489daf19fdSVijay Mahadevan #else
499daf19fdSVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), NULL, dmmoab->fileset);
509daf19fdSVijay Mahadevan #endif
51b117cd09SVijay Mahadevan 
52e882eb38SVijay Mahadevan   ierr = PetscMalloc1(nlevels + 1, &dmmoab->hsets);CHKERRQ(ierr);
53e882eb38SVijay Mahadevan 
54b117cd09SVijay Mahadevan   /* generate the mesh hierarchy */
559c368985SVijay Mahadevan   merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false);MBERRNM(merr);
56e882eb38SVijay Mahadevan 
579daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
58755f3dfbSVijay Mahadevan   if (dmmoab->pcomm->size() > 1) {
5964e1c140SVijay Mahadevan     merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);MBERRNM(merr);
60755f3dfbSVijay Mahadevan   }
619daf19fdSVijay Mahadevan #endif
6249d66b22SVijay Mahadevan 
6364e1c140SVijay Mahadevan   /* copy the mesh sets for nested refinement hierarchy */
64e92d1c7cSVijay Mahadevan   dmmoab->hsets[0] = hsets[0];
65e92d1c7cSVijay Mahadevan   for (ilevel = 1; ilevel <= nlevels; ilevel++)
662417220eSVijay Mahadevan   {
672417220eSVijay Mahadevan     dmmoab->hsets[ilevel] = hsets[ilevel];
68e92d1c7cSVijay Mahadevan 
699c368985SVijay Mahadevan #ifdef MOAB_HAVE_MPI
70e92d1c7cSVijay Mahadevan     merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false);MBERRNM(merr);
719c368985SVijay Mahadevan #endif
72e92d1c7cSVijay Mahadevan 
73e92d1c7cSVijay Mahadevan     /* Update material and other geometric tags from parent to child sets */
74e92d1c7cSVijay Mahadevan     merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]);MBERRNM(merr);
752417220eSVijay Mahadevan   }
7664e1c140SVijay Mahadevan 
7764e1c140SVijay Mahadevan   hsets.clear();
78b117cd09SVijay Mahadevan   if (!ldegrees) {
79b117cd09SVijay Mahadevan     ierr = PetscFree(pdegrees);CHKERRQ(ierr);
80b117cd09SVijay Mahadevan   }
81b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
82b117cd09SVijay Mahadevan }
83b117cd09SVijay Mahadevan 
84b117cd09SVijay Mahadevan /*@
85e882eb38SVijay Mahadevan   DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy
86e882eb38SVijay Mahadevan   by succesively refining a coarse mesh.
87b117cd09SVijay Mahadevan 
88b117cd09SVijay Mahadevan   Collective on MPI_Comm
89b117cd09SVijay Mahadevan 
90b117cd09SVijay Mahadevan   Input Parameter:
91*a2b725a8SWilliam Gropp . dm  - The DMMoab object
92b117cd09SVijay Mahadevan 
93b117cd09SVijay Mahadevan   Output Parameter:
94e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
95*a2b725a8SWilliam Gropp - dmf  - The DM objects after successive refinement of the hierarchy
96b117cd09SVijay Mahadevan 
97b117cd09SVijay Mahadevan   Level: beginner
98b117cd09SVijay Mahadevan 
99e882eb38SVijay Mahadevan .keywords: DMMoab, generate, refinement
100b117cd09SVijay Mahadevan @*/
101304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode  DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[])
102b117cd09SVijay Mahadevan {
103b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
104b117cd09SVijay Mahadevan   PetscInt        i;
105b117cd09SVijay Mahadevan 
106b117cd09SVijay Mahadevan   PetscFunctionBegin;
107b117cd09SVijay Mahadevan 
108e882eb38SVijay Mahadevan   ierr = DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]);CHKERRQ(ierr);
109e882eb38SVijay Mahadevan   for (i = 1; i < nlevels; i++) {
110e882eb38SVijay Mahadevan     ierr = DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]);CHKERRQ(ierr);
111b117cd09SVijay Mahadevan   }
112b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
113b117cd09SVijay Mahadevan }
114b117cd09SVijay Mahadevan 
115b117cd09SVijay Mahadevan /*@
116e882eb38SVijay Mahadevan   DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy
117e882eb38SVijay Mahadevan   by succesively coarsening a refined mesh.
118b117cd09SVijay Mahadevan 
119b117cd09SVijay Mahadevan   Collective on MPI_Comm
120b117cd09SVijay Mahadevan 
121b117cd09SVijay Mahadevan   Input Parameter:
122*a2b725a8SWilliam Gropp . dm  - The DMMoab object
123b117cd09SVijay Mahadevan 
124b117cd09SVijay Mahadevan   Output Parameter:
125e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
126*a2b725a8SWilliam Gropp - dmc  - The DM objects after successive coarsening of the hierarchy
127b117cd09SVijay Mahadevan 
128b117cd09SVijay Mahadevan   Level: beginner
129b117cd09SVijay Mahadevan 
130e882eb38SVijay Mahadevan .keywords: DMMoab, generate, coarsening
131b117cd09SVijay Mahadevan @*/
132304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[])
133b117cd09SVijay Mahadevan {
134b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
135b117cd09SVijay Mahadevan   PetscInt        i;
136b117cd09SVijay Mahadevan 
137b117cd09SVijay Mahadevan   PetscFunctionBegin;
138b117cd09SVijay Mahadevan 
139e882eb38SVijay Mahadevan   ierr = DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]);CHKERRQ(ierr);
140e882eb38SVijay Mahadevan   for (i = 1; i < nlevels; i++) {
141e882eb38SVijay Mahadevan     ierr = DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]);CHKERRQ(ierr);
142b117cd09SVijay Mahadevan   }
143b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
144b117cd09SVijay Mahadevan }
145b117cd09SVijay Mahadevan 
1462417220eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt*, PetscInt*, PetscInt*, PetscInt*, PetscBool);
147b117cd09SVijay Mahadevan 
148b117cd09SVijay Mahadevan /*@
149e882eb38SVijay Mahadevan   DMCreateInterpolation_Moab - Generate the interpolation operators to transform
150e882eb38SVijay Mahadevan   operators (matrices, vectors) from parent level to child level as defined by
151e882eb38SVijay Mahadevan   the DM inputs provided by the user.
152b117cd09SVijay Mahadevan 
153b117cd09SVijay Mahadevan   Collective on MPI_Comm
154b117cd09SVijay Mahadevan 
155b117cd09SVijay Mahadevan   Input Parameter:
156e882eb38SVijay Mahadevan + dm1  - The DMMoab object
157e882eb38SVijay Mahadevan - dm2  - the second, finer DMMoab object
158b117cd09SVijay Mahadevan 
159b117cd09SVijay Mahadevan   Output Parameter:
160e882eb38SVijay Mahadevan + interpl  - The interpolation operator for transferring data between the levels
161e882eb38SVijay Mahadevan - vec      - The scaling vector (optional)
162b117cd09SVijay Mahadevan 
163e882eb38SVijay Mahadevan   Level: developer
164b117cd09SVijay Mahadevan 
165b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
166b117cd09SVijay Mahadevan @*/
167304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat* interpl, Vec* vec)
168b117cd09SVijay Mahadevan {
169755f3dfbSVijay Mahadevan   DM_Moab         *dmbp, *dmbc;
170b117cd09SVijay Mahadevan   PetscErrorCode   ierr;
171b117cd09SVijay Mahadevan   moab::ErrorCode  merr;
172e882eb38SVijay Mahadevan   PetscInt         dim;
1733f1c6e43SVijay Mahadevan   PetscReal        factor;
174ce27a4eeSVijay Mahadevan   PetscInt         innz, *nnz, ionz, *onz;
175755f3dfbSVijay Mahadevan   PetscInt         nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc;
176a86ed394SVijay Mahadevan   const PetscBool  use_consistent_bases=PETSC_TRUE;
177b117cd09SVijay Mahadevan 
178b117cd09SVijay Mahadevan   PetscFunctionBegin;
179755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmp, DM_CLASSID, 1);
180755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmc, DM_CLASSID, 2);
181755f3dfbSVijay Mahadevan   dmbp = (DM_Moab*)(dmp)->data;
182755f3dfbSVijay Mahadevan   dmbc = (DM_Moab*)(dmc)->data;
183755f3dfbSVijay Mahadevan   nlsizp = dmbp->nloc;// *dmb1->numFields;
184755f3dfbSVijay Mahadevan   nlsizc = dmbc->nloc;// *dmb2->numFields;
185755f3dfbSVijay Mahadevan   ngsizp = dmbp->n; // *dmb1->numFields;
186755f3dfbSVijay Mahadevan   ngsizc = dmbc->n; // *dmb2->numFields;
187755f3dfbSVijay Mahadevan   nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields;
188b117cd09SVijay Mahadevan 
1892417220eSVijay Mahadevan   // Columns = Parent DoFs ;  Rows = Child DoFs
190755f3dfbSVijay Mahadevan   // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent))
191755f3dfbSVijay Mahadevan   // Size: nlsizc * nlghsizp
192755f3dfbSVijay Mahadevan   PetscInfo4(NULL, "Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n", ngsizc, nlghsizp, dmbp->hlevel, dmbc->hlevel);
193b117cd09SVijay Mahadevan 
194a86ed394SVijay Mahadevan   ierr = DMGetDimension(dmp, &dim);CHKERRQ(ierr);
195a86ed394SVijay Mahadevan 
196941e0cffSVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
197755f3dfbSVijay Mahadevan   ierr = PetscCalloc2(nlsizc, &nnz, nlsizc, &onz);CHKERRQ(ierr);
198941e0cffSVijay Mahadevan 
199941e0cffSVijay Mahadevan   /* Loop through the local elements and compute the relation between the current parent and the refined_level. */
2002417220eSVijay Mahadevan   for (moab::Range::iterator iter = dmbc->vowned->begin(); iter != dmbc->vowned->end(); iter++) {
201941e0cffSVijay Mahadevan 
2022417220eSVijay Mahadevan     const moab::EntityHandle vhandle = *iter;
2032417220eSVijay Mahadevan     /* define local variables */
2042417220eSVijay Mahadevan     moab::EntityHandle parent;
2052417220eSVijay Mahadevan     std::vector<moab::EntityHandle> adjs;
2062417220eSVijay Mahadevan     moab::Range     found;
2072417220eSVijay Mahadevan 
2082417220eSVijay Mahadevan     /* store the vertex DoF number */
2092417220eSVijay Mahadevan     const int ldof = dmbc->lidmap[vhandle - dmbc->seqstart];
2102417220eSVijay Mahadevan 
2112417220eSVijay Mahadevan     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
2122417220eSVijay Mahadevan        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
2132417220eSVijay Mahadevan        non-zero pattern accordingly. */
2142417220eSVijay Mahadevan     merr = dmbc->hierarchy->get_adjacencies(vhandle, dmbc->dim, adjs);MBERRNM(merr);
2152417220eSVijay Mahadevan 
2162417220eSVijay Mahadevan     /* loop over vertices and update the number of connectivity */
2172417220eSVijay Mahadevan     for (unsigned jter = 0; jter < adjs.size(); jter++) {
2182417220eSVijay Mahadevan 
2192417220eSVijay Mahadevan       const moab::EntityHandle jhandle = adjs[jter];
220941e0cffSVijay Mahadevan 
221941e0cffSVijay Mahadevan       /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
2222417220eSVijay Mahadevan       merr = dmbc->hierarchy->child_to_parent(jhandle, dmbc->hlevel, dmbp->hlevel, &parent);MBERRNM(merr);
223941e0cffSVijay Mahadevan 
2242417220eSVijay Mahadevan       /* Get connectivity information in canonical ordering for the local element */
2252417220eSVijay Mahadevan       std::vector<moab::EntityHandle> connp;
2262417220eSVijay Mahadevan       merr = dmbp->hierarchy->get_connectivity(parent, dmbp->hlevel, connp);MBERRNM(merr);
227a86ed394SVijay Mahadevan 
2282417220eSVijay Mahadevan       for (unsigned ic=0; ic < connp.size(); ++ic) {
2292417220eSVijay Mahadevan 
2302417220eSVijay Mahadevan         /* loop over each element connected to the adjacent vertex and update as needed */
2312417220eSVijay Mahadevan         /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
2322417220eSVijay Mahadevan         if (found.find(connp[ic]) != found.end()) continue; /* make sure we don't double count shared vertices */
2332417220eSVijay Mahadevan         if (dmbp->vghost->find(connp[ic]) != dmbp->vghost->end()) onz[ldof]++; /* update out-of-proc onz */
2342417220eSVijay Mahadevan         else nnz[ldof]++; /* else local vertex */
2352417220eSVijay Mahadevan         found.insert(connp[ic]);
2362417220eSVijay Mahadevan       }
2372417220eSVijay Mahadevan     }
238941e0cffSVijay Mahadevan   }
239941e0cffSVijay Mahadevan 
2402417220eSVijay Mahadevan   for (int i = 0; i < nlsizc; i++)
2412417220eSVijay Mahadevan     nnz[i] += 1; /* self count the node */
242941e0cffSVijay Mahadevan 
243941e0cffSVijay Mahadevan   ionz = onz[0];
244941e0cffSVijay Mahadevan   innz = nnz[0];
245755f3dfbSVijay Mahadevan   for (int tc = 0; tc < nlsizc; tc++) {
246941e0cffSVijay Mahadevan     // check for maximum allowed sparsity = fully dense
247755f3dfbSVijay Mahadevan     nnz[tc] = std::min(nlsizp, nnz[tc]);
2482417220eSVijay Mahadevan     onz[tc] = std::min(ngsizp - nlsizp, onz[tc]);
2492417220eSVijay Mahadevan 
2502417220eSVijay Mahadevan     PetscInfo3(NULL, "  %d: NNZ = %d, ONZ = %d\n", tc, nnz[tc], onz[tc]);
251941e0cffSVijay Mahadevan 
252941e0cffSVijay Mahadevan     innz = (innz < nnz[tc] ? nnz[tc] : innz);
253941e0cffSVijay Mahadevan     ionz = (ionz < onz[tc] ? onz[tc] : ionz);
254941e0cffSVijay Mahadevan   }
25564e1c140SVijay Mahadevan 
25664e1c140SVijay Mahadevan   /* create interpolation matrix */
257755f3dfbSVijay Mahadevan   ierr = MatCreate(PetscObjectComm((PetscObject)dmc), interpl);CHKERRQ(ierr);
258755f3dfbSVijay Mahadevan   ierr = MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp);CHKERRQ(ierr);
25964e1c140SVijay Mahadevan   ierr = MatSetType(*interpl, MATAIJ);CHKERRQ(ierr);
26064e1c140SVijay Mahadevan   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
26164e1c140SVijay Mahadevan 
262941e0cffSVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(*interpl, innz, nnz);CHKERRQ(ierr);
263941e0cffSVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz);CHKERRQ(ierr);
264941e0cffSVijay Mahadevan 
265941e0cffSVijay Mahadevan   /* clean up temporary memory */
266941e0cffSVijay Mahadevan   ierr = PetscFree2(nnz, onz);CHKERRQ(ierr);
267b117cd09SVijay Mahadevan 
268b117cd09SVijay Mahadevan   /* set up internal matrix data-structures */
269b117cd09SVijay Mahadevan   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
270b117cd09SVijay Mahadevan 
271a86ed394SVijay Mahadevan   /* Define variables for assembly */
272a86ed394SVijay Mahadevan   std::vector<moab::EntityHandle> children;
273a86ed394SVijay Mahadevan   std::vector<moab::EntityHandle> connp, connc;
274a86ed394SVijay Mahadevan   std::vector<PetscReal> pcoords, ccoords, values_phi;
275a86ed394SVijay Mahadevan 
276a86ed394SVijay Mahadevan   if (use_consistent_bases) {
2772417220eSVijay Mahadevan     const moab::EntityHandle ehandle = dmbp->elocal->front();
278a86ed394SVijay Mahadevan 
279a86ed394SVijay Mahadevan     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
280a86ed394SVijay Mahadevan 
281a86ed394SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
282a86ed394SVijay Mahadevan     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
2832417220eSVijay Mahadevan     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr);
284a86ed394SVijay Mahadevan 
285a86ed394SVijay Mahadevan     std::vector<PetscReal> natparam(3*connc.size(), 0.0);
286a86ed394SVijay Mahadevan     pcoords.resize(connp.size() * 3);
287a86ed394SVijay Mahadevan     ccoords.resize(connc.size() * 3);
288a86ed394SVijay Mahadevan     values_phi.resize(connp.size()*connc.size());
289a86ed394SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
290a86ed394SVijay Mahadevan     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr);
291a86ed394SVijay Mahadevan     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr);
292a86ed394SVijay Mahadevan 
293a86ed394SVijay Mahadevan     /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
294a86ed394SVijay Mahadevan     for (unsigned tc = 0; tc < connc.size(); tc++) {
295a86ed394SVijay Mahadevan       const PetscInt offset = tc * 3;
296a86ed394SVijay Mahadevan 
297a86ed394SVijay Mahadevan       /* Scale ccoords relative to pcoords */
298a86ed394SVijay Mahadevan       ierr = DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size()*tc]);CHKERRQ(ierr);
299a86ed394SVijay Mahadevan     }
300a86ed394SVijay Mahadevan   }
301a86ed394SVijay Mahadevan   else {
302a86ed394SVijay Mahadevan     factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0);
303a86ed394SVijay Mahadevan   }
304a86ed394SVijay Mahadevan 
305755f3dfbSVijay Mahadevan   /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */
3062417220eSVijay Mahadevan   ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
307b117cd09SVijay Mahadevan 
308e882eb38SVijay Mahadevan   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
3092417220eSVijay Mahadevan   for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) {
310b117cd09SVijay Mahadevan 
311b117cd09SVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
312b117cd09SVijay Mahadevan 
313b117cd09SVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
314a86ed394SVijay Mahadevan     children.clear();
315a86ed394SVijay Mahadevan     connc.clear();
316755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
317b117cd09SVijay Mahadevan 
318b117cd09SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
319755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
3202417220eSVijay Mahadevan     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr);
321b117cd09SVijay Mahadevan 
322a86ed394SVijay Mahadevan     pcoords.resize(connp.size() * 3);
323a86ed394SVijay Mahadevan     ccoords.resize(connc.size() * 3);
324b117cd09SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
325755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr);
326755f3dfbSVijay Mahadevan     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr);
327b117cd09SVijay Mahadevan 
328b117cd09SVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
329b117cd09SVijay Mahadevan     /* TODO: specific to scalar system - use GetDofs */
330755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr);
331755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr);
332b117cd09SVijay Mahadevan 
333a86ed394SVijay Mahadevan     /* Compute the actual interpolation weights when projecting solution/residual between levels */
334a86ed394SVijay Mahadevan     if (use_consistent_bases) {
335a86ed394SVijay Mahadevan 
336a86ed394SVijay Mahadevan       /* Use the cached values of natural parameteric coordinates and basis pre-evaluated.
337a86ed394SVijay Mahadevan          We are making an assumption here that UMR used in GMG to generate the hierarchy uses
338a86ed394SVijay Mahadevan          the same template for all elements; This will fail for mixed element meshes (TRI/QUAD).
339a86ed394SVijay Mahadevan 
3402417220eSVijay Mahadevan          TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes)
341a86ed394SVijay Mahadevan       */
342a86ed394SVijay Mahadevan 
343a86ed394SVijay Mahadevan       /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */
344a86ed394SVijay Mahadevan       for (unsigned tc = 0; tc < connc.size(); tc++) {
345a86ed394SVijay Mahadevan         /* TODO: Check if we should be using INSERT_VALUES instead */
346a86ed394SVijay Mahadevan         ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[connp.size()*tc], ADD_VALUES);CHKERRQ(ierr);
347a86ed394SVijay Mahadevan       }
348a86ed394SVijay Mahadevan     }
349a86ed394SVijay Mahadevan     else {
350e882eb38SVijay Mahadevan       /* Compute the interpolation weights by determining distance of 1-ring
351755f3dfbSVijay Mahadevan          neighbor vertices from current vertex
352755f3dfbSVijay Mahadevan 
353a86ed394SVijay Mahadevan          This should be used only when FEM basis is not used for the discretization.
354a86ed394SVijay Mahadevan          Else, the consistent interface to compute the basis function for interpolation
355a86ed394SVijay Mahadevan          between the levels should be evaluated correctly to preserve convergence of GMG.
3562417220eSVijay Mahadevan          Shephard's basis will be terrible for any unsmooth problems.
357755f3dfbSVijay Mahadevan       */
358755f3dfbSVijay Mahadevan       values_phi.resize(connp.size());
359755f3dfbSVijay Mahadevan       for (unsigned tc = 0; tc < connc.size(); tc++) {
360755f3dfbSVijay Mahadevan 
361a86ed394SVijay Mahadevan         PetscReal normsum = 0.0;
362755f3dfbSVijay Mahadevan         for (unsigned tp = 0; tp < connp.size(); tp++) {
363755f3dfbSVijay Mahadevan           values_phi[tp] = 0.0;
364e882eb38SVijay Mahadevan           for (unsigned k = 0; k < 3; k++)
365755f3dfbSVijay Mahadevan             values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim);
366755f3dfbSVijay Mahadevan           if (values_phi[tp] < 1e-12) {
367755f3dfbSVijay Mahadevan             values_phi[tp] = 1e12;
368b117cd09SVijay Mahadevan           }
369b117cd09SVijay Mahadevan           else {
370755f3dfbSVijay Mahadevan             //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim);
371755f3dfbSVijay Mahadevan             values_phi[tp] = std::pow(values_phi[tp], -1.0);
372755f3dfbSVijay Mahadevan             normsum += values_phi[tp];
373b117cd09SVijay Mahadevan           }
374b117cd09SVijay Mahadevan         }
375755f3dfbSVijay Mahadevan         for (unsigned tp = 0; tp < connp.size(); tp++) {
376755f3dfbSVijay Mahadevan           if (values_phi[tp] > 1e11)
377755f3dfbSVijay Mahadevan             values_phi[tp] = factor * 0.5 / connp.size();
378b117cd09SVijay Mahadevan           else
379755f3dfbSVijay Mahadevan             values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum);
380b117cd09SVijay Mahadevan         }
381755f3dfbSVijay Mahadevan         ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
382b117cd09SVijay Mahadevan       }
383a86ed394SVijay Mahadevan     }
384b117cd09SVijay Mahadevan   }
385304006b3SVijay Mahadevan   if (vec) *vec = NULL;
386b117cd09SVijay Mahadevan   ierr = MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
387b117cd09SVijay Mahadevan   ierr = MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
388b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
389b117cd09SVijay Mahadevan }
390b117cd09SVijay Mahadevan 
391b117cd09SVijay Mahadevan /*@
392b117cd09SVijay Mahadevan   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
393b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
394b117cd09SVijay Mahadevan   provided by the user.
395b117cd09SVijay Mahadevan 
396b117cd09SVijay Mahadevan   Collective on MPI_Comm
397b117cd09SVijay Mahadevan 
398b117cd09SVijay Mahadevan   Input Parameter:
399b117cd09SVijay Mahadevan . dmb  - The DMMoab object
400b117cd09SVijay Mahadevan 
401b117cd09SVijay Mahadevan   Output Parameter:
402*a2b725a8SWilliam Gropp + nlevels   - The number of levels of refinement needed to generate the hierarchy
403*a2b725a8SWilliam Gropp - ldegrees  - The degree of refinement at each level in the hierarchy
404b117cd09SVijay Mahadevan 
405b117cd09SVijay Mahadevan   Level: beginner
406b117cd09SVijay Mahadevan 
407b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
408b117cd09SVijay Mahadevan @*/
409304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter* ctx)
410b117cd09SVijay Mahadevan {
411e882eb38SVijay Mahadevan   //DM_Moab        *dmmoab;
412b117cd09SVijay Mahadevan 
413b117cd09SVijay Mahadevan   PetscFunctionBegin;
414b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1, DM_CLASSID, 1);
415b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2, DM_CLASSID, 2);
416e882eb38SVijay Mahadevan   //dmmoab = (DM_Moab*)(dm1)->data;
417b117cd09SVijay Mahadevan 
418b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
419b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
420b117cd09SVijay Mahadevan }
421b117cd09SVijay Mahadevan 
422470880abSPatrick Sanan static PetscErrorCode DMMoab_UMR_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref)
423b117cd09SVijay Mahadevan {
424b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
425e882eb38SVijay Mahadevan   PetscInt        i, dim;
426b117cd09SVijay Mahadevan   DM              dm2;
427b117cd09SVijay Mahadevan   moab::ErrorCode merr;
428b117cd09SVijay Mahadevan   DM_Moab        *dmb = (DM_Moab*)dm->data, *dd2;
429b117cd09SVijay Mahadevan 
430b117cd09SVijay Mahadevan   PetscFunctionBegin;
431b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
432b117cd09SVijay Mahadevan   PetscValidPointer(dmref, 4);
433b117cd09SVijay Mahadevan 
434b117cd09SVijay Mahadevan   if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) {
435e882eb38SVijay 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);
436e882eb38SVijay Mahadevan     if (dmb->hlevel - 1 < 0 && !refine) PetscInfo1(NULL, "Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n", dmb->hlevel - 1);
437e882eb38SVijay Mahadevan     *dmref = PETSC_NULL;
438b117cd09SVijay Mahadevan     PetscFunctionReturn(0);
439b117cd09SVijay Mahadevan   }
440b117cd09SVijay Mahadevan 
441b117cd09SVijay Mahadevan   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
442b117cd09SVijay Mahadevan   dd2 = (DM_Moab*)dm2->data;
443b117cd09SVijay Mahadevan 
444b117cd09SVijay Mahadevan   dd2->mbiface = dmb->mbiface;
4459daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
446b117cd09SVijay Mahadevan   dd2->pcomm = dmb->pcomm;
4479daf19fdSVijay Mahadevan #endif
448b117cd09SVijay Mahadevan   dd2->icreatedinstance = PETSC_FALSE;
44964e1c140SVijay Mahadevan   dd2->nghostrings = dmb->nghostrings;
450b117cd09SVijay Mahadevan 
451b117cd09SVijay Mahadevan   /* set the new level based on refinement/coarsening */
452b117cd09SVijay Mahadevan   if (refine) {
453b117cd09SVijay Mahadevan     dd2->hlevel = dmb->hlevel + 1;
454b117cd09SVijay Mahadevan   }
455b117cd09SVijay Mahadevan   else {
456b117cd09SVijay Mahadevan     dd2->hlevel = dmb->hlevel - 1;
457b117cd09SVijay Mahadevan   }
458b117cd09SVijay Mahadevan 
459b117cd09SVijay Mahadevan   /* Copy the multilevel hierarchy pointers in MOAB */
460b117cd09SVijay Mahadevan   dd2->hierarchy = dmb->hierarchy;
461b117cd09SVijay Mahadevan   dd2->nhlevels = dmb->nhlevels;
462b117cd09SVijay Mahadevan   ierr = PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets);CHKERRQ(ierr);
463b117cd09SVijay Mahadevan   for (i = 0; i <= dd2->nhlevels; i++) {
464b117cd09SVijay Mahadevan     dd2->hsets[i] = dmb->hsets[i];
465b117cd09SVijay Mahadevan   }
466b117cd09SVijay Mahadevan   dd2->fileset = dd2->hsets[dd2->hlevel];
467b117cd09SVijay Mahadevan 
468b117cd09SVijay Mahadevan   /* do the remaining initializations for DMMoab */
469b117cd09SVijay Mahadevan   dd2->bs = dmb->bs;
470b117cd09SVijay Mahadevan   dd2->numFields = dmb->numFields;
471b117cd09SVijay Mahadevan   dd2->rw_dbglevel = dmb->rw_dbglevel;
472b117cd09SVijay Mahadevan   dd2->partition_by_rank = dmb->partition_by_rank;
473b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
474b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
475b117cd09SVijay Mahadevan   dd2->read_mode = dmb->read_mode;
476b117cd09SVijay Mahadevan   dd2->write_mode = dmb->write_mode;
477b117cd09SVijay Mahadevan 
478b117cd09SVijay Mahadevan   /* set global ID tag handle */
479b117cd09SVijay Mahadevan   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
480b117cd09SVijay Mahadevan 
481b117cd09SVijay Mahadevan   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
482b117cd09SVijay Mahadevan 
483b117cd09SVijay Mahadevan   ierr = DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix);CHKERRQ(ierr);
484b117cd09SVijay Mahadevan   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
485b117cd09SVijay Mahadevan   ierr = DMSetDimension(dm2, dim);CHKERRQ(ierr);
486b117cd09SVijay Mahadevan 
487b117cd09SVijay Mahadevan   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
488b117cd09SVijay Mahadevan   dm2->ops->creatematrix = dm->ops->creatematrix;
489b117cd09SVijay Mahadevan 
490b117cd09SVijay Mahadevan   /* copy fill information if given */
491b117cd09SVijay Mahadevan   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
492b117cd09SVijay Mahadevan 
493b117cd09SVijay Mahadevan   /* copy vector type information */
494b117cd09SVijay Mahadevan   ierr = DMSetMatType(dm2, dm->mattype);CHKERRQ(ierr);
495b117cd09SVijay Mahadevan   ierr = DMSetVecType(dm2, dm->vectype);CHKERRQ(ierr);
4963f1c6e43SVijay Mahadevan   dd2->numFields = dmb->numFields;
4973f1c6e43SVijay Mahadevan   if (dmb->numFields) {
498b117cd09SVijay Mahadevan     ierr = DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames);CHKERRQ(ierr);
4993f1c6e43SVijay Mahadevan   }
500b117cd09SVijay Mahadevan 
501b117cd09SVijay Mahadevan   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
502b117cd09SVijay Mahadevan 
503b117cd09SVijay Mahadevan   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
504b117cd09SVijay Mahadevan   ierr = DMSetUp(dm2);CHKERRQ(ierr);
505b117cd09SVijay Mahadevan 
506b117cd09SVijay Mahadevan   *dmref = dm2;
507b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
508b117cd09SVijay Mahadevan }
509b117cd09SVijay Mahadevan 
510b117cd09SVijay Mahadevan 
511b117cd09SVijay Mahadevan /*@
512b117cd09SVijay Mahadevan   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
513b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
514b117cd09SVijay Mahadevan   provided by the user.
515b117cd09SVijay Mahadevan 
516e882eb38SVijay Mahadevan   Collective on DM
517b117cd09SVijay Mahadevan 
518b117cd09SVijay Mahadevan   Input Parameter:
519e882eb38SVijay Mahadevan + dm  - The DMMoab object
520e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
521b117cd09SVijay Mahadevan 
522b117cd09SVijay Mahadevan   Output Parameter:
523e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL
524b117cd09SVijay Mahadevan 
525e882eb38SVijay Mahadevan   Note: If no refinement was done, the return value is NULL
526e882eb38SVijay Mahadevan 
527e882eb38SVijay Mahadevan   Level: developer
528b117cd09SVijay Mahadevan 
529b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
530b117cd09SVijay Mahadevan @*/
531304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM* dmf)
532b117cd09SVijay Mahadevan {
533b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
534b117cd09SVijay Mahadevan 
535b117cd09SVijay Mahadevan   PetscFunctionBegin;
536b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
537b117cd09SVijay Mahadevan 
538470880abSPatrick Sanan   ierr = DMMoab_UMR_Private(dm, comm, PETSC_TRUE, dmf);CHKERRQ(ierr);
539b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
540b117cd09SVijay Mahadevan }
541b117cd09SVijay Mahadevan 
542b117cd09SVijay Mahadevan /*@
543b117cd09SVijay Mahadevan   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
544b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
545b117cd09SVijay Mahadevan   provided by the user.
546b117cd09SVijay Mahadevan 
547e882eb38SVijay Mahadevan   Collective on DM
548b117cd09SVijay Mahadevan 
549b117cd09SVijay Mahadevan   Input Parameter:
550*a2b725a8SWilliam Gropp + dm  - The DMMoab object
551e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
552b117cd09SVijay Mahadevan 
553b117cd09SVijay Mahadevan   Output Parameter:
554e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL
555b117cd09SVijay Mahadevan 
556e882eb38SVijay Mahadevan   Note: If no coarsening was done, the return value is NULL
557b117cd09SVijay Mahadevan 
558e882eb38SVijay Mahadevan   Level: developer
559e882eb38SVijay Mahadevan 
560e882eb38SVijay Mahadevan .keywords: DMMoab, create, coarsening
561b117cd09SVijay Mahadevan @*/
562304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM* dmc)
563b117cd09SVijay Mahadevan {
564b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
565b117cd09SVijay Mahadevan 
566b117cd09SVijay Mahadevan   PetscFunctionBegin;
567b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
568b117cd09SVijay Mahadevan 
569470880abSPatrick Sanan   ierr = DMMoab_UMR_Private(dm, comm, PETSC_FALSE, dmc);CHKERRQ(ierr);
570b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
571b117cd09SVijay Mahadevan }
572