xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision 9daf19fd0d22011a060902ed53673b48ac0c4e7d)
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 #undef __FUNCT__
7b117cd09SVijay Mahadevan #define __FUNCT__ "DMMoabGenerateHierarchy"
8b117cd09SVijay Mahadevan /*@
9b117cd09SVijay Mahadevan   DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy
10b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
11b117cd09SVijay Mahadevan   provided by the user.
12b117cd09SVijay Mahadevan 
13b117cd09SVijay Mahadevan   Collective on MPI_Comm
14b117cd09SVijay Mahadevan 
15b117cd09SVijay Mahadevan   Input Parameter:
16b117cd09SVijay Mahadevan + dmb  - The DMMoab object
17b117cd09SVijay Mahadevan 
18b117cd09SVijay Mahadevan   Output Parameter:
19b117cd09SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
20b117cd09SVijay Mahadevan . ldegrees  - The degree of refinement at each level in the hierarchy
21b117cd09SVijay Mahadevan 
22b117cd09SVijay Mahadevan   Level: beginner
23b117cd09SVijay Mahadevan 
24b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
25b117cd09SVijay Mahadevan @*/
26b117cd09SVijay Mahadevan PetscErrorCode DMMoabGenerateHierarchy(DM dm,PetscInt nlevels,PetscInt *ldegrees)
27b117cd09SVijay Mahadevan {
28b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
29b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
30b117cd09SVijay Mahadevan   moab::ErrorCode merr;
31b117cd09SVijay Mahadevan   PetscInt *pdegrees,i;
32e882eb38SVijay Mahadevan   std::vector<moab::EntityHandle> hsets;
33b117cd09SVijay Mahadevan 
34b117cd09SVijay Mahadevan   PetscFunctionBegin;
35b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
37b117cd09SVijay Mahadevan 
38b117cd09SVijay Mahadevan   if (!ldegrees) {
39b117cd09SVijay Mahadevan     ierr = PetscMalloc1(nlevels,&pdegrees);CHKERRQ(ierr);
40b117cd09SVijay Mahadevan     for (i=0; i<nlevels; i++) pdegrees[i]=2; /* default = Degree 2 refinement */
41b117cd09SVijay Mahadevan   }
42b117cd09SVijay Mahadevan   else pdegrees=ldegrees;
43b117cd09SVijay Mahadevan 
44b117cd09SVijay Mahadevan   /* initialize set level refinement data for hierarchy */
45b117cd09SVijay Mahadevan   dmmoab->nhlevels=nlevels;
46b117cd09SVijay Mahadevan 
47b117cd09SVijay Mahadevan   /* Instantiate the nested refinement class */
48*9daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
493f1c6e43SVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset);
50*9daf19fdSVijay Mahadevan #else
51*9daf19fdSVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), NULL, dmmoab->fileset);
52*9daf19fdSVijay Mahadevan #endif
53b117cd09SVijay Mahadevan 
54e882eb38SVijay Mahadevan   ierr = PetscMalloc1(nlevels+1,&dmmoab->hsets);CHKERRQ(ierr);
55e882eb38SVijay Mahadevan 
56b117cd09SVijay Mahadevan   /* generate the mesh hierarchy */
5764e1c140SVijay Mahadevan   merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets);MBERRNM(merr);
58e882eb38SVijay Mahadevan 
59*9daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
60755f3dfbSVijay Mahadevan   if (dmmoab->pcomm->size() > 1) {
6164e1c140SVijay Mahadevan     merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);MBERRNM(merr);
62755f3dfbSVijay Mahadevan   }
63*9daf19fdSVijay Mahadevan #endif
6449d66b22SVijay Mahadevan 
6564e1c140SVijay Mahadevan   /* copy the mesh sets for nested refinement hierarchy */
6664e1c140SVijay Mahadevan   for (i=0; i<=nlevels; i++)
6749d66b22SVijay Mahadevan     dmmoab->hsets[i]=hsets[i];
6864e1c140SVijay Mahadevan 
6964e1c140SVijay Mahadevan   hsets.clear();
70b117cd09SVijay Mahadevan   if (!ldegrees) {
71b117cd09SVijay Mahadevan     ierr = PetscFree(pdegrees);CHKERRQ(ierr);
72b117cd09SVijay Mahadevan   }
73b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
74b117cd09SVijay Mahadevan }
75b117cd09SVijay Mahadevan 
76b117cd09SVijay Mahadevan #undef __FUNCT__
77b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefineHierarchy_Moab"
78b117cd09SVijay Mahadevan /*@
79e882eb38SVijay Mahadevan   DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy
80e882eb38SVijay Mahadevan   by succesively refining a coarse mesh.
81b117cd09SVijay Mahadevan 
82b117cd09SVijay Mahadevan   Collective on MPI_Comm
83b117cd09SVijay Mahadevan 
84b117cd09SVijay Mahadevan   Input Parameter:
85e882eb38SVijay Mahadevan + dm  - The DMMoab object
86b117cd09SVijay Mahadevan 
87b117cd09SVijay Mahadevan   Output Parameter:
88e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
89e882eb38SVijay Mahadevan . dmf  - The DM objects after successive refinement of the hierarchy
90b117cd09SVijay Mahadevan 
91b117cd09SVijay Mahadevan   Level: beginner
92b117cd09SVijay Mahadevan 
93e882eb38SVijay Mahadevan .keywords: DMMoab, generate, refinement
94b117cd09SVijay Mahadevan @*/
95304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode  DMRefineHierarchy_Moab(DM dm,PetscInt nlevels,DM dmf[])
96b117cd09SVijay Mahadevan {
97b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
98b117cd09SVijay Mahadevan   PetscInt        i;
99b117cd09SVijay Mahadevan 
100b117cd09SVijay Mahadevan   PetscFunctionBegin;
101b117cd09SVijay Mahadevan 
102e882eb38SVijay Mahadevan   ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr);
103e882eb38SVijay Mahadevan   for (i=1; i<nlevels; i++) {
104e882eb38SVijay Mahadevan     ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr);
105b117cd09SVijay Mahadevan   }
106b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
107b117cd09SVijay Mahadevan }
108b117cd09SVijay Mahadevan 
109b117cd09SVijay Mahadevan #undef __FUNCT__
110b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsenHierarchy_Moab"
111b117cd09SVijay Mahadevan /*@
112e882eb38SVijay Mahadevan   DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy
113e882eb38SVijay Mahadevan   by succesively coarsening a refined mesh.
114b117cd09SVijay Mahadevan 
115b117cd09SVijay Mahadevan   Collective on MPI_Comm
116b117cd09SVijay Mahadevan 
117b117cd09SVijay Mahadevan   Input Parameter:
118e882eb38SVijay Mahadevan + dm  - The DMMoab object
119b117cd09SVijay Mahadevan 
120b117cd09SVijay Mahadevan   Output Parameter:
121e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
122e882eb38SVijay Mahadevan . dmc  - The DM objects after successive coarsening of the hierarchy
123b117cd09SVijay Mahadevan 
124b117cd09SVijay Mahadevan   Level: beginner
125b117cd09SVijay Mahadevan 
126e882eb38SVijay Mahadevan .keywords: DMMoab, generate, coarsening
127b117cd09SVijay Mahadevan @*/
128304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm,PetscInt nlevels,DM dmc[])
129b117cd09SVijay Mahadevan {
130b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
131b117cd09SVijay Mahadevan   PetscInt        i;
132b117cd09SVijay Mahadevan 
133b117cd09SVijay Mahadevan   PetscFunctionBegin;
134b117cd09SVijay Mahadevan 
135e882eb38SVijay Mahadevan   ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr);
136e882eb38SVijay Mahadevan   for (i=1; i<nlevels; i++) {
137e882eb38SVijay Mahadevan     ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr);
138b117cd09SVijay Mahadevan   }
139b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
140b117cd09SVijay Mahadevan }
141b117cd09SVijay Mahadevan 
142b117cd09SVijay Mahadevan 
143b117cd09SVijay Mahadevan #undef __FUNCT__
144b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInterpolation_Moab"
145b117cd09SVijay Mahadevan /*@
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 
150b117cd09SVijay Mahadevan   Collective on MPI_Comm
151b117cd09SVijay Mahadevan 
152b117cd09SVijay Mahadevan   Input Parameter:
153e882eb38SVijay Mahadevan + dm1  - The DMMoab object
154e882eb38SVijay Mahadevan - dm2  - the second, finer DMMoab object
155b117cd09SVijay Mahadevan 
156b117cd09SVijay Mahadevan   Output Parameter:
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 .keywords: DMMoab, create, refinement
163b117cd09SVijay Mahadevan @*/
164304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp,DM dmc,Mat* interpl,Vec* vec)
165b117cd09SVijay Mahadevan {
166755f3dfbSVijay Mahadevan   DM_Moab         *dmbp, *dmbc;
167b117cd09SVijay Mahadevan   PetscErrorCode   ierr;
168b117cd09SVijay Mahadevan   moab::ErrorCode  merr;
169e882eb38SVijay Mahadevan   PetscInt         dim;
1703f1c6e43SVijay Mahadevan   PetscReal        factor;
171ce27a4eeSVijay Mahadevan   PetscInt         innz, *nnz, ionz, *onz;
172755f3dfbSVijay Mahadevan   PetscInt         nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc;
173755f3dfbSVijay Mahadevan   moab::Range      eowned;
174b117cd09SVijay Mahadevan 
175b117cd09SVijay Mahadevan   PetscFunctionBegin;
176755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmp,DM_CLASSID,1);
177755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmc,DM_CLASSID,2);
178755f3dfbSVijay Mahadevan   dmbp = (DM_Moab*)(dmp)->data;
179755f3dfbSVijay Mahadevan   dmbc = (DM_Moab*)(dmc)->data;
180755f3dfbSVijay Mahadevan   nlsizp = dmbp->nloc;// *dmb1->numFields;
181755f3dfbSVijay Mahadevan   nlsizc = dmbc->nloc;// *dmb2->numFields;
182755f3dfbSVijay Mahadevan   ngsizp = dmbp->n; // *dmb1->numFields;
183755f3dfbSVijay Mahadevan   ngsizc = dmbc->n; // *dmb2->numFields;
184755f3dfbSVijay Mahadevan   nlghsizp = (dmbp->nloc+dmbp->nghost); // *dmb1->numFields;
185b117cd09SVijay Mahadevan 
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 
190941e0cffSVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
191755f3dfbSVijay Mahadevan   ierr = PetscCalloc2(nlsizc,&nnz,nlsizc,&onz);CHKERRQ(ierr);
192941e0cffSVijay Mahadevan 
193755f3dfbSVijay Mahadevan   eowned = *dmbp->elocal;
194*9daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
195755f3dfbSVijay Mahadevan   merr = dmbp->pcomm->filter_pstatus(eowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
196*9daf19fdSVijay Mahadevan #endif
197941e0cffSVijay Mahadevan   /* Loop through the local elements and compute the relation between the current parent and the refined_level. */
198755f3dfbSVijay Mahadevan   for(moab::Range::iterator iter = eowned.begin(); iter!= eowned.end(); iter++) {
199941e0cffSVijay Mahadevan 
200941e0cffSVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
201941e0cffSVijay Mahadevan     std::vector<moab::EntityHandle> children;
202941e0cffSVijay Mahadevan     std::vector<moab::EntityHandle> connp, connc;
203755f3dfbSVijay Mahadevan     moab::Range connc_owned;
204941e0cffSVijay Mahadevan 
205941e0cffSVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
206755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
207941e0cffSVijay Mahadevan 
208755f3dfbSVijay Mahadevan     /* Get connectivity of the parent elements */
209755f3dfbSVijay Mahadevan     //merr = dmb1->mbiface->get_connectivity(ehandle,connect,vpere,false);MBERRNM(merr);
210755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
211755f3dfbSVijay Mahadevan     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc_owned);MBERRNM(merr);
212*9daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
213755f3dfbSVijay Mahadevan     merr = dmbc->pcomm->filter_pstatus(connc_owned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
214*9daf19fdSVijay Mahadevan #endif
215755f3dfbSVijay Mahadevan     for (unsigned tc=0; tc<connc_owned.size(); tc++) {
216755f3dfbSVijay Mahadevan       connc.push_back(connc_owned[tc]);
217941e0cffSVijay Mahadevan     }
218941e0cffSVijay Mahadevan 
219941e0cffSVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
220755f3dfbSVijay Mahadevan     /* get DoFs for both coarse and fine scale grid */
221755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlockedLocal(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr);
222755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlockedLocal(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr);
223941e0cffSVijay Mahadevan 
224941e0cffSVijay Mahadevan     for (unsigned tp=0;tp<connp.size(); tp++) {
225755f3dfbSVijay Mahadevan       if (dmbp->vowned->find(connp[tp]) != dmbp->vowned->end()) {
22664e1c140SVijay Mahadevan         for (unsigned tc=0;tc<connc.size(); tc++) {
227755f3dfbSVijay Mahadevan           if (dmbc->vowned->find(connc[tc]) != dmbc->vowned->end()) {
228755f3dfbSVijay Mahadevan             /* FIXME: This will over-allocate since we multi-count shared nodes.
229755f3dfbSVijay Mahadevan                Disadvantage of using element traversal for nnz computation. */
23064e1c140SVijay Mahadevan             nnz[dofsc[tc]]++;
23164e1c140SVijay Mahadevan           }
232941e0cffSVijay Mahadevan         }
233755f3dfbSVijay Mahadevan       }
234755f3dfbSVijay Mahadevan       else if (dmbp->vghost->find(connp[tp]) != dmbp->vghost->end()) {
23564e1c140SVijay Mahadevan         for (unsigned tc=0;tc<connc.size(); tc++) {
236755f3dfbSVijay Mahadevan           if (dmbc->vowned->find(connc[tc]) != dmbc->vowned->end()) {
237755f3dfbSVijay Mahadevan             /* FIXME: This will over-allocate since we multi-count shared nodes.
238755f3dfbSVijay Mahadevan                Disadvantage of using element traversal for nnz computation. */
23964e1c140SVijay Mahadevan             onz[dofsc[tc]]++;
240941e0cffSVijay Mahadevan           }
241941e0cffSVijay Mahadevan         }
24264e1c140SVijay Mahadevan       }
243755f3dfbSVijay Mahadevan       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in child level %D\n", connp[tp]);
24449d66b22SVijay Mahadevan     }
245941e0cffSVijay Mahadevan   }
246941e0cffSVijay Mahadevan 
247941e0cffSVijay Mahadevan   ionz=onz[0];
248941e0cffSVijay Mahadevan   innz=nnz[0];
249755f3dfbSVijay Mahadevan   for (int tc=0; tc < nlsizc; tc++) {
250941e0cffSVijay Mahadevan     // check for maximum allowed sparsity = fully dense
251755f3dfbSVijay Mahadevan     nnz[tc] = std::min(nlsizp,nnz[tc]);
252755f3dfbSVijay Mahadevan     onz[tc] = std::min(nlghsizp-nlsizp,onz[tc]);
253941e0cffSVijay Mahadevan 
254941e0cffSVijay Mahadevan     innz = (innz < nnz[tc] ? nnz[tc] : innz);
255941e0cffSVijay Mahadevan     ionz = (ionz < onz[tc] ? onz[tc] : ionz);
256941e0cffSVijay Mahadevan   }
25764e1c140SVijay Mahadevan 
25864e1c140SVijay Mahadevan   /* create interpolation matrix */
259755f3dfbSVijay Mahadevan   ierr = MatCreate(PetscObjectComm((PetscObject)dmc), interpl);CHKERRQ(ierr);
260755f3dfbSVijay Mahadevan   ierr = MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp);CHKERRQ(ierr);
26164e1c140SVijay Mahadevan   //ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr);
26264e1c140SVijay Mahadevan   ierr = MatSetType(*interpl,MATAIJ);CHKERRQ(ierr);
26364e1c140SVijay Mahadevan   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
26464e1c140SVijay Mahadevan 
265941e0cffSVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(*interpl,innz,nnz);CHKERRQ(ierr);
266941e0cffSVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(*interpl,innz,nnz,ionz,onz);CHKERRQ(ierr);
267941e0cffSVijay Mahadevan 
268941e0cffSVijay Mahadevan   /* clean up temporary memory */
269941e0cffSVijay Mahadevan   ierr = PetscFree2(nnz,onz);CHKERRQ(ierr);
270b117cd09SVijay Mahadevan 
271b117cd09SVijay Mahadevan   /* set up internal matrix data-structures */
272b117cd09SVijay Mahadevan   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
273b117cd09SVijay Mahadevan 
274755f3dfbSVijay Mahadevan   /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */
275755f3dfbSVijay Mahadevan   // ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
276b117cd09SVijay Mahadevan 
277755f3dfbSVijay Mahadevan   ierr = DMGetDimension(dmp, &dim);CHKERRQ(ierr);
278b117cd09SVijay Mahadevan 
279755f3dfbSVijay Mahadevan   factor = std::pow(2.0 /*degree_P_for_refinement*/,(dmbc->hlevel-dmbp->hlevel)*dmbp->dim*1.0);
28049d66b22SVijay Mahadevan 
281e882eb38SVijay Mahadevan   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
282755f3dfbSVijay Mahadevan   for(moab::Range::iterator iter = eowned.begin(); iter!= eowned.end(); iter++) {
283b117cd09SVijay Mahadevan 
284b117cd09SVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
285b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> children;
286b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> connp, connc;
287755f3dfbSVijay Mahadevan     moab::Range connc_owned;
288b117cd09SVijay Mahadevan 
289b117cd09SVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
290755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
291b117cd09SVijay Mahadevan 
292b117cd09SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
293755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
294755f3dfbSVijay Mahadevan     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc_owned);MBERRNM(merr);
295*9daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
296755f3dfbSVijay Mahadevan     merr = dmbc->pcomm->filter_pstatus(connc_owned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
297*9daf19fdSVijay Mahadevan #endif
298755f3dfbSVijay Mahadevan     for (unsigned tc=0; tc<connc_owned.size(); tc++) {
299755f3dfbSVijay Mahadevan       connc.push_back(connc_owned[tc]);
300b117cd09SVijay Mahadevan     }
301b117cd09SVijay Mahadevan 
302b117cd09SVijay Mahadevan     std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size());
303b117cd09SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
304755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr);
305755f3dfbSVijay Mahadevan     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr);
306b117cd09SVijay Mahadevan 
307b117cd09SVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
308b117cd09SVijay Mahadevan     /* TODO: specific to scalar system - use GetDofs */
309755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr);
310755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr);
311b117cd09SVijay Mahadevan 
312e882eb38SVijay Mahadevan     /* Compute the interpolation weights by determining distance of 1-ring
313755f3dfbSVijay Mahadevan        neighbor vertices from current vertex
314755f3dfbSVijay Mahadevan 
315755f3dfbSVijay Mahadevan        TODO: This needs to be replaced with a consistent interface to compute
316755f3dfbSVijay Mahadevan        the basis function interpolation between the levels evaluated correctly.
317755f3dfbSVijay Mahadevan        RBF basis will be terrible for any unsmooth problems..
318755f3dfbSVijay Mahadevan        -- NEEDS TO BE REMOVED --
319755f3dfbSVijay Mahadevan     */
320755f3dfbSVijay Mahadevan     values_phi.resize(connp.size());
321755f3dfbSVijay Mahadevan     for (unsigned tc=0;tc<connc.size(); tc++) {
322755f3dfbSVijay Mahadevan 
323b117cd09SVijay Mahadevan       double normsum=0.0;
324755f3dfbSVijay Mahadevan       for (unsigned tp=0;tp<connp.size(); tp++) {
325755f3dfbSVijay Mahadevan         values_phi[tp] = 0.0;
326e882eb38SVijay Mahadevan         for (unsigned k=0;k<3; k++)
327755f3dfbSVijay Mahadevan           values_phi[tp] += std::pow(pcoords[tp*3+k]-ccoords[k+tc*3], dim);
328755f3dfbSVijay Mahadevan         if (values_phi[tp] < 1e-12) {
329755f3dfbSVijay Mahadevan           values_phi[tp] = 1e12;
330b117cd09SVijay Mahadevan         }
331b117cd09SVijay Mahadevan         else {
332755f3dfbSVijay Mahadevan           //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim);
333755f3dfbSVijay Mahadevan           values_phi[tp] = std::pow(values_phi[tp], -1.0);
334755f3dfbSVijay Mahadevan           normsum += values_phi[tp];
335b117cd09SVijay Mahadevan         }
336b117cd09SVijay Mahadevan       }
337755f3dfbSVijay Mahadevan       for (unsigned tp=0;tp<connp.size(); tp++) {
338755f3dfbSVijay Mahadevan         if (values_phi[tp] > 1e11)
339755f3dfbSVijay Mahadevan           values_phi[tp] = factor*0.5/connp.size();
340b117cd09SVijay Mahadevan         else
341755f3dfbSVijay Mahadevan           values_phi[tp] = factor*values_phi[tp]*0.5/(connp.size()*normsum);
342b117cd09SVijay Mahadevan       }
343755f3dfbSVijay Mahadevan       ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
344b117cd09SVijay Mahadevan     }
345b117cd09SVijay Mahadevan 
346b117cd09SVijay Mahadevan     //get interpolation weights
347b117cd09SVijay Mahadevan     //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr);
348e882eb38SVijay Mahadevan     // for (int j=0;j<dofs_per_element; j++)
349e882eb38SVijay Mahadevan     //  std::cout<<"values "<<values_phi[j]<<std::endl;
350b117cd09SVijay Mahadevan 
351b117cd09SVijay Mahadevan   }
352304006b3SVijay Mahadevan   if (vec) *vec=NULL;
353b117cd09SVijay Mahadevan   ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
354b117cd09SVijay Mahadevan   ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
355b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
356b117cd09SVijay Mahadevan }
357b117cd09SVijay Mahadevan 
358b117cd09SVijay Mahadevan #undef __FUNCT__
359b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInjection_Moab"
360b117cd09SVijay Mahadevan /*@
361b117cd09SVijay Mahadevan   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
362b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
363b117cd09SVijay Mahadevan   provided by the user.
364b117cd09SVijay Mahadevan 
365b117cd09SVijay Mahadevan   Collective on MPI_Comm
366b117cd09SVijay Mahadevan 
367b117cd09SVijay Mahadevan   Input Parameter:
368b117cd09SVijay Mahadevan . dmb  - The DMMoab object
369b117cd09SVijay Mahadevan 
370b117cd09SVijay Mahadevan   Output Parameter:
371b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
372b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
373b117cd09SVijay Mahadevan 
374b117cd09SVijay Mahadevan   Level: beginner
375b117cd09SVijay Mahadevan 
376b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
377b117cd09SVijay Mahadevan @*/
378304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx)
379b117cd09SVijay Mahadevan {
380e882eb38SVijay Mahadevan   //DM_Moab        *dmmoab;
381b117cd09SVijay Mahadevan 
382b117cd09SVijay Mahadevan   PetscFunctionBegin;
383b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
384b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
385e882eb38SVijay Mahadevan   //dmmoab = (DM_Moab*)(dm1)->data;
386b117cd09SVijay Mahadevan 
387b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
388b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
389b117cd09SVijay Mahadevan }
390b117cd09SVijay Mahadevan 
391b117cd09SVijay Mahadevan 
392b117cd09SVijay Mahadevan #undef __FUNCT__
393b117cd09SVijay Mahadevan #define __FUNCT__ "DM_UMR_Moab_Private"
394b117cd09SVijay Mahadevan PetscErrorCode  DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref)
395b117cd09SVijay Mahadevan {
396b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
397e882eb38SVijay Mahadevan   PetscInt        i,dim;
398b117cd09SVijay Mahadevan   DM              dm2;
399b117cd09SVijay Mahadevan   moab::ErrorCode merr;
400b117cd09SVijay Mahadevan   DM_Moab        *dmb = (DM_Moab*)dm->data,*dd2;
401b117cd09SVijay Mahadevan 
402b117cd09SVijay Mahadevan   PetscFunctionBegin;
403b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
404b117cd09SVijay Mahadevan   PetscValidPointer(dmref,4);
405b117cd09SVijay Mahadevan 
406b117cd09SVijay Mahadevan   if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) {
407e882eb38SVijay 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);
408e882eb38SVijay Mahadevan     if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1);
409e882eb38SVijay Mahadevan     *dmref = PETSC_NULL;
410b117cd09SVijay Mahadevan     PetscFunctionReturn(0);
411b117cd09SVijay Mahadevan   }
412b117cd09SVijay Mahadevan 
413b117cd09SVijay Mahadevan   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
414b117cd09SVijay Mahadevan   dd2 = (DM_Moab*)dm2->data;
415b117cd09SVijay Mahadevan 
416b117cd09SVijay Mahadevan   dd2->mbiface = dmb->mbiface;
417*9daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI
418b117cd09SVijay Mahadevan   dd2->pcomm = dmb->pcomm;
419*9daf19fdSVijay Mahadevan #endif
420b117cd09SVijay Mahadevan   dd2->icreatedinstance = PETSC_FALSE;
42164e1c140SVijay Mahadevan   dd2->nghostrings=dmb->nghostrings;
422b117cd09SVijay Mahadevan 
423b117cd09SVijay Mahadevan   /* set the new level based on refinement/coarsening */
424b117cd09SVijay Mahadevan   if (refine) {
425b117cd09SVijay Mahadevan     dd2->hlevel=dmb->hlevel+1;
426b117cd09SVijay Mahadevan   }
427b117cd09SVijay Mahadevan   else {
428b117cd09SVijay Mahadevan     dd2->hlevel=dmb->hlevel-1;
429b117cd09SVijay Mahadevan   }
430b117cd09SVijay Mahadevan 
431b117cd09SVijay Mahadevan   /* Copy the multilevel hierarchy pointers in MOAB */
432b117cd09SVijay Mahadevan   dd2->hierarchy = dmb->hierarchy;
433b117cd09SVijay Mahadevan   dd2->nhlevels = dmb->nhlevels;
434b117cd09SVijay Mahadevan   ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr);
435b117cd09SVijay Mahadevan   for (i=0; i<=dd2->nhlevels; i++) {
436b117cd09SVijay Mahadevan     dd2->hsets[i]=dmb->hsets[i];
437b117cd09SVijay Mahadevan   }
438b117cd09SVijay Mahadevan   dd2->fileset = dd2->hsets[dd2->hlevel];
439b117cd09SVijay Mahadevan 
440b117cd09SVijay Mahadevan   /* do the remaining initializations for DMMoab */
441b117cd09SVijay Mahadevan   dd2->bs = dmb->bs;
442b117cd09SVijay Mahadevan   dd2->numFields = dmb->numFields;
443b117cd09SVijay Mahadevan   dd2->rw_dbglevel = dmb->rw_dbglevel;
444b117cd09SVijay Mahadevan   dd2->partition_by_rank = dmb->partition_by_rank;
445b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
446b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
447b117cd09SVijay Mahadevan   dd2->read_mode = dmb->read_mode;
448b117cd09SVijay Mahadevan   dd2->write_mode = dmb->write_mode;
449b117cd09SVijay Mahadevan 
450b117cd09SVijay Mahadevan   /* set global ID tag handle */
451b117cd09SVijay Mahadevan   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
452b117cd09SVijay Mahadevan 
453b117cd09SVijay Mahadevan   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
454b117cd09SVijay Mahadevan 
455b117cd09SVijay Mahadevan   ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr);
456b117cd09SVijay Mahadevan   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
457b117cd09SVijay Mahadevan   ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr);
458b117cd09SVijay Mahadevan 
459b117cd09SVijay Mahadevan   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
460b117cd09SVijay Mahadevan   dm2->ops->creatematrix = dm->ops->creatematrix;
461b117cd09SVijay Mahadevan 
462b117cd09SVijay Mahadevan   /* copy fill information if given */
463b117cd09SVijay Mahadevan   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
464b117cd09SVijay Mahadevan 
465b117cd09SVijay Mahadevan   /* copy vector type information */
466b117cd09SVijay Mahadevan   ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr);
467b117cd09SVijay Mahadevan   ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr);
4683f1c6e43SVijay Mahadevan   dd2->numFields = dmb->numFields;
4693f1c6e43SVijay Mahadevan   if (dmb->numFields) {
470b117cd09SVijay Mahadevan     ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr);
4713f1c6e43SVijay Mahadevan   }
472b117cd09SVijay Mahadevan 
473b117cd09SVijay Mahadevan   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
474b117cd09SVijay Mahadevan 
475b117cd09SVijay Mahadevan   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
476b117cd09SVijay Mahadevan   ierr = DMSetUp(dm2);CHKERRQ(ierr);
477b117cd09SVijay Mahadevan 
478b117cd09SVijay Mahadevan   *dmref = dm2;
479b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
480b117cd09SVijay Mahadevan }
481b117cd09SVijay Mahadevan 
482b117cd09SVijay Mahadevan 
483b117cd09SVijay Mahadevan #undef __FUNCT__
484b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefine_Moab"
485b117cd09SVijay Mahadevan /*@
486b117cd09SVijay Mahadevan   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
487b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
488b117cd09SVijay Mahadevan   provided by the user.
489b117cd09SVijay Mahadevan 
490e882eb38SVijay Mahadevan   Collective on DM
491b117cd09SVijay Mahadevan 
492b117cd09SVijay Mahadevan   Input Parameter:
493e882eb38SVijay Mahadevan + dm  - The DMMoab object
494e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
495b117cd09SVijay Mahadevan 
496b117cd09SVijay Mahadevan   Output Parameter:
497e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL
498b117cd09SVijay Mahadevan 
499e882eb38SVijay Mahadevan   Note: If no refinement was done, the return value is NULL
500e882eb38SVijay Mahadevan 
501e882eb38SVijay Mahadevan   Level: developer
502b117cd09SVijay Mahadevan 
503b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
504b117cd09SVijay Mahadevan @*/
505304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf)
506b117cd09SVijay Mahadevan {
507b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
508b117cd09SVijay Mahadevan 
509b117cd09SVijay Mahadevan   PetscFunctionBegin;
510b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
511b117cd09SVijay Mahadevan 
512b117cd09SVijay Mahadevan   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr);
513b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
514b117cd09SVijay Mahadevan }
515b117cd09SVijay Mahadevan 
516b117cd09SVijay Mahadevan #undef __FUNCT__
517b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsen_Moab"
518b117cd09SVijay Mahadevan /*@
519b117cd09SVijay Mahadevan   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
520b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
521b117cd09SVijay Mahadevan   provided by the user.
522b117cd09SVijay Mahadevan 
523e882eb38SVijay Mahadevan   Collective on DM
524b117cd09SVijay Mahadevan 
525b117cd09SVijay Mahadevan   Input Parameter:
526e882eb38SVijay Mahadevan . dm  - The DMMoab object
527e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
528b117cd09SVijay Mahadevan 
529b117cd09SVijay Mahadevan   Output Parameter:
530e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL
531b117cd09SVijay Mahadevan 
532e882eb38SVijay Mahadevan   Note: If no coarsening was done, the return value is NULL
533b117cd09SVijay Mahadevan 
534e882eb38SVijay Mahadevan   Level: developer
535e882eb38SVijay Mahadevan 
536e882eb38SVijay Mahadevan .keywords: DMMoab, create, coarsening
537b117cd09SVijay Mahadevan @*/
538304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc)
539b117cd09SVijay Mahadevan {
540b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
541b117cd09SVijay Mahadevan 
542b117cd09SVijay Mahadevan   PetscFunctionBegin;
543b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
544b117cd09SVijay Mahadevan 
545b117cd09SVijay Mahadevan   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr);
546b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
547b117cd09SVijay Mahadevan }
548