xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision 755f3dfbc21a83a2a40e224cc9e5f876bede072f)
1*755f3dfbSVijay 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 */
483f1c6e43SVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset);
49b117cd09SVijay Mahadevan 
50e882eb38SVijay Mahadevan   ierr = PetscMalloc1(nlevels+1,&dmmoab->hsets);CHKERRQ(ierr);
51e882eb38SVijay Mahadevan 
52b117cd09SVijay Mahadevan   /* generate the mesh hierarchy */
5364e1c140SVijay Mahadevan   merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets);MBERRNM(merr);
54e882eb38SVijay Mahadevan 
55*755f3dfbSVijay Mahadevan   if (dmmoab->pcomm->size() > 1) {
5664e1c140SVijay Mahadevan     merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);MBERRNM(merr);
57*755f3dfbSVijay Mahadevan   }
5849d66b22SVijay Mahadevan 
5964e1c140SVijay Mahadevan   /* copy the mesh sets for nested refinement hierarchy */
6064e1c140SVijay Mahadevan   for (i=0; i<=nlevels; i++)
6149d66b22SVijay Mahadevan     dmmoab->hsets[i]=hsets[i];
6264e1c140SVijay Mahadevan 
6364e1c140SVijay Mahadevan   hsets.clear();
64b117cd09SVijay Mahadevan   if (!ldegrees) {
65b117cd09SVijay Mahadevan     ierr = PetscFree(pdegrees);CHKERRQ(ierr);
66b117cd09SVijay Mahadevan   }
67b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
68b117cd09SVijay Mahadevan }
69b117cd09SVijay Mahadevan 
70b117cd09SVijay Mahadevan #undef __FUNCT__
71b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefineHierarchy_Moab"
72b117cd09SVijay Mahadevan /*@
73e882eb38SVijay Mahadevan   DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy
74e882eb38SVijay Mahadevan   by succesively refining a coarse mesh.
75b117cd09SVijay Mahadevan 
76b117cd09SVijay Mahadevan   Collective on MPI_Comm
77b117cd09SVijay Mahadevan 
78b117cd09SVijay Mahadevan   Input Parameter:
79e882eb38SVijay Mahadevan + dm  - The DMMoab object
80b117cd09SVijay Mahadevan 
81b117cd09SVijay Mahadevan   Output Parameter:
82e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
83e882eb38SVijay Mahadevan . dmf  - The DM objects after successive refinement of the hierarchy
84b117cd09SVijay Mahadevan 
85b117cd09SVijay Mahadevan   Level: beginner
86b117cd09SVijay Mahadevan 
87e882eb38SVijay Mahadevan .keywords: DMMoab, generate, refinement
88b117cd09SVijay Mahadevan @*/
89b117cd09SVijay Mahadevan PetscErrorCode  DMRefineHierarchy_Moab(DM dm,PetscInt nlevels,DM dmf[])
90b117cd09SVijay Mahadevan {
91b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
92b117cd09SVijay Mahadevan   PetscInt        i;
93b117cd09SVijay Mahadevan 
94b117cd09SVijay Mahadevan   PetscFunctionBegin;
95b117cd09SVijay Mahadevan 
96e882eb38SVijay Mahadevan   ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr);
97e882eb38SVijay Mahadevan   for (i=1; i<nlevels; i++) {
98e882eb38SVijay Mahadevan     ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr);
99b117cd09SVijay Mahadevan   }
100b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
101b117cd09SVijay Mahadevan }
102b117cd09SVijay Mahadevan 
103b117cd09SVijay Mahadevan #undef __FUNCT__
104b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsenHierarchy_Moab"
105b117cd09SVijay Mahadevan /*@
106e882eb38SVijay Mahadevan   DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy
107e882eb38SVijay Mahadevan   by succesively coarsening a refined mesh.
108b117cd09SVijay Mahadevan 
109b117cd09SVijay Mahadevan   Collective on MPI_Comm
110b117cd09SVijay Mahadevan 
111b117cd09SVijay Mahadevan   Input Parameter:
112e882eb38SVijay Mahadevan + dm  - The DMMoab object
113b117cd09SVijay Mahadevan 
114b117cd09SVijay Mahadevan   Output Parameter:
115e882eb38SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
116e882eb38SVijay Mahadevan . dmc  - The DM objects after successive coarsening of the hierarchy
117b117cd09SVijay Mahadevan 
118b117cd09SVijay Mahadevan   Level: beginner
119b117cd09SVijay Mahadevan 
120e882eb38SVijay Mahadevan .keywords: DMMoab, generate, coarsening
121b117cd09SVijay Mahadevan @*/
122b117cd09SVijay Mahadevan PetscErrorCode DMCoarsenHierarchy_Moab(DM dm,PetscInt nlevels,DM dmc[])
123b117cd09SVijay Mahadevan {
124b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
125b117cd09SVijay Mahadevan   PetscInt        i;
126b117cd09SVijay Mahadevan 
127b117cd09SVijay Mahadevan   PetscFunctionBegin;
128b117cd09SVijay Mahadevan 
129e882eb38SVijay Mahadevan   ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr);
130e882eb38SVijay Mahadevan   for (i=1; i<nlevels; i++) {
131e882eb38SVijay Mahadevan     ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr);
132b117cd09SVijay Mahadevan   }
133b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
134b117cd09SVijay Mahadevan }
135b117cd09SVijay Mahadevan 
136b117cd09SVijay Mahadevan 
137b117cd09SVijay Mahadevan #undef __FUNCT__
138b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInterpolation_Moab"
139b117cd09SVijay Mahadevan /*@
140e882eb38SVijay Mahadevan   DMCreateInterpolation_Moab - Generate the interpolation operators to transform
141e882eb38SVijay Mahadevan   operators (matrices, vectors) from parent level to child level as defined by
142e882eb38SVijay Mahadevan   the DM inputs provided by the user.
143b117cd09SVijay Mahadevan 
144b117cd09SVijay Mahadevan   Collective on MPI_Comm
145b117cd09SVijay Mahadevan 
146b117cd09SVijay Mahadevan   Input Parameter:
147e882eb38SVijay Mahadevan + dm1  - The DMMoab object
148e882eb38SVijay Mahadevan - dm2  - the second, finer DMMoab object
149b117cd09SVijay Mahadevan 
150b117cd09SVijay Mahadevan   Output Parameter:
151e882eb38SVijay Mahadevan + interpl  - The interpolation operator for transferring data between the levels
152e882eb38SVijay Mahadevan - vec      - The scaling vector (optional)
153b117cd09SVijay Mahadevan 
154e882eb38SVijay Mahadevan   Level: developer
155b117cd09SVijay Mahadevan 
156b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
157b117cd09SVijay Mahadevan @*/
158*755f3dfbSVijay Mahadevan PetscErrorCode DMCreateInterpolation_Moab(DM dmp,DM dmc,Mat* interpl,Vec* vec)
159b117cd09SVijay Mahadevan {
160*755f3dfbSVijay Mahadevan   DM_Moab         *dmbp, *dmbc;
161b117cd09SVijay Mahadevan   PetscErrorCode   ierr;
162b117cd09SVijay Mahadevan   moab::ErrorCode  merr;
163e882eb38SVijay Mahadevan   PetscInt         dim;
1643f1c6e43SVijay Mahadevan   PetscReal        factor;
165ce27a4eeSVijay Mahadevan   PetscInt         innz, *nnz, ionz, *onz;
166*755f3dfbSVijay Mahadevan   PetscInt         nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc;
167*755f3dfbSVijay Mahadevan   moab::Range      eowned;
168b117cd09SVijay Mahadevan 
169b117cd09SVijay Mahadevan   PetscFunctionBegin;
170*755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmp,DM_CLASSID,1);
171*755f3dfbSVijay Mahadevan   PetscValidHeaderSpecific(dmc,DM_CLASSID,2);
172*755f3dfbSVijay Mahadevan   dmbp = (DM_Moab*)(dmp)->data;
173*755f3dfbSVijay Mahadevan   dmbc = (DM_Moab*)(dmc)->data;
174*755f3dfbSVijay Mahadevan   nlsizp = dmbp->nloc;// *dmb1->numFields;
175*755f3dfbSVijay Mahadevan   nlsizc = dmbc->nloc;// *dmb2->numFields;
176*755f3dfbSVijay Mahadevan   ngsizp = dmbp->n; // *dmb1->numFields;
177*755f3dfbSVijay Mahadevan   ngsizc = dmbc->n; // *dmb2->numFields;
178*755f3dfbSVijay Mahadevan   nlghsizp = (dmbp->nloc+dmbp->nghost); // *dmb1->numFields;
179b117cd09SVijay Mahadevan 
180*755f3dfbSVijay Mahadevan   // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent))
181*755f3dfbSVijay Mahadevan   // Size: nlsizc * nlghsizp
182*755f3dfbSVijay Mahadevan   PetscInfo4(NULL,"Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n",ngsizc,nlghsizp,dmbp->hlevel,dmbc->hlevel);
183b117cd09SVijay Mahadevan 
184941e0cffSVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
185*755f3dfbSVijay Mahadevan   ierr = PetscCalloc2(nlsizc,&nnz,nlsizc,&onz);CHKERRQ(ierr);
186941e0cffSVijay Mahadevan 
187*755f3dfbSVijay Mahadevan   eowned = *dmbp->elocal;
188*755f3dfbSVijay Mahadevan   merr = dmbp->pcomm->filter_pstatus(eowned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
189941e0cffSVijay Mahadevan   /* Loop through the local elements and compute the relation between the current parent and the refined_level. */
190*755f3dfbSVijay Mahadevan   for(moab::Range::iterator iter = eowned.begin(); iter!= eowned.end(); iter++) {
191941e0cffSVijay Mahadevan 
192941e0cffSVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
193941e0cffSVijay Mahadevan     std::vector<moab::EntityHandle> children;
194941e0cffSVijay Mahadevan     std::vector<moab::EntityHandle> connp, connc;
195*755f3dfbSVijay Mahadevan     moab::Range connc_owned;
196941e0cffSVijay Mahadevan 
197941e0cffSVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
198*755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
199941e0cffSVijay Mahadevan 
200*755f3dfbSVijay Mahadevan     /* Get connectivity of the parent elements */
201*755f3dfbSVijay Mahadevan     //merr = dmb1->mbiface->get_connectivity(ehandle,connect,vpere,false);MBERRNM(merr);
202*755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
203*755f3dfbSVijay Mahadevan     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc_owned);MBERRNM(merr);
204*755f3dfbSVijay Mahadevan     merr = dmbc->pcomm->filter_pstatus(connc_owned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
205*755f3dfbSVijay Mahadevan     for (unsigned tc=0; tc<connc_owned.size(); tc++) {
206*755f3dfbSVijay Mahadevan       connc.push_back(connc_owned[tc]);
207941e0cffSVijay Mahadevan     }
208941e0cffSVijay Mahadevan 
209941e0cffSVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
210*755f3dfbSVijay Mahadevan     /* get DoFs for both coarse and fine scale grid */
211*755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlockedLocal(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr);
212*755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlockedLocal(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr);
213941e0cffSVijay Mahadevan 
214941e0cffSVijay Mahadevan     for (unsigned tp=0;tp<connp.size(); tp++) {
215*755f3dfbSVijay Mahadevan       if (dmbp->vowned->find(connp[tp]) != dmbp->vowned->end()) {
21664e1c140SVijay Mahadevan         for (unsigned tc=0;tc<connc.size(); tc++) {
217*755f3dfbSVijay Mahadevan           if (dmbc->vowned->find(connc[tc]) != dmbc->vowned->end()) {
218*755f3dfbSVijay Mahadevan             /* FIXME: This will over-allocate since we multi-count shared nodes.
219*755f3dfbSVijay Mahadevan                Disadvantage of using element traversal for nnz computation. */
22064e1c140SVijay Mahadevan             nnz[dofsc[tc]]++;
22164e1c140SVijay Mahadevan           }
222941e0cffSVijay Mahadevan         }
223*755f3dfbSVijay Mahadevan       }
224*755f3dfbSVijay Mahadevan       else if (dmbp->vghost->find(connp[tp]) != dmbp->vghost->end()) {
22564e1c140SVijay Mahadevan         for (unsigned tc=0;tc<connc.size(); tc++) {
226*755f3dfbSVijay Mahadevan           if (dmbc->vowned->find(connc[tc]) != dmbc->vowned->end()) {
227*755f3dfbSVijay Mahadevan             /* FIXME: This will over-allocate since we multi-count shared nodes.
228*755f3dfbSVijay Mahadevan                Disadvantage of using element traversal for nnz computation. */
22964e1c140SVijay Mahadevan             onz[dofsc[tc]]++;
230941e0cffSVijay Mahadevan           }
231941e0cffSVijay Mahadevan         }
23264e1c140SVijay Mahadevan       }
233*755f3dfbSVijay Mahadevan       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in child level %D\n", connp[tp]);
23449d66b22SVijay Mahadevan     }
235941e0cffSVijay Mahadevan   }
236941e0cffSVijay Mahadevan 
237941e0cffSVijay Mahadevan   ionz=onz[0];
238941e0cffSVijay Mahadevan   innz=nnz[0];
239*755f3dfbSVijay Mahadevan   for (int tc=0; tc < nlsizc; tc++) {
240941e0cffSVijay Mahadevan     // check for maximum allowed sparsity = fully dense
241*755f3dfbSVijay Mahadevan     nnz[tc] = std::min(nlsizp,nnz[tc]);
242*755f3dfbSVijay Mahadevan     onz[tc] = std::min(nlghsizp-nlsizp,onz[tc]);
243941e0cffSVijay Mahadevan 
244941e0cffSVijay Mahadevan     innz = (innz < nnz[tc] ? nnz[tc] : innz);
245941e0cffSVijay Mahadevan     ionz = (ionz < onz[tc] ? onz[tc] : ionz);
246941e0cffSVijay Mahadevan   }
24764e1c140SVijay Mahadevan 
24864e1c140SVijay Mahadevan   /* create interpolation matrix */
249*755f3dfbSVijay Mahadevan   ierr = MatCreate(PetscObjectComm((PetscObject)dmc), interpl);CHKERRQ(ierr);
250*755f3dfbSVijay Mahadevan   ierr = MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp);CHKERRQ(ierr);
25164e1c140SVijay Mahadevan   //ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr);
25264e1c140SVijay Mahadevan   ierr = MatSetType(*interpl,MATAIJ);CHKERRQ(ierr);
25364e1c140SVijay Mahadevan   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
25464e1c140SVijay Mahadevan 
255941e0cffSVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(*interpl,innz,nnz);CHKERRQ(ierr);
256941e0cffSVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(*interpl,innz,nnz,ionz,onz);CHKERRQ(ierr);
257941e0cffSVijay Mahadevan 
258941e0cffSVijay Mahadevan   /* clean up temporary memory */
259941e0cffSVijay Mahadevan   ierr = PetscFree2(nnz,onz);CHKERRQ(ierr);
260b117cd09SVijay Mahadevan 
261b117cd09SVijay Mahadevan   /* set up internal matrix data-structures */
262b117cd09SVijay Mahadevan   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
263b117cd09SVijay Mahadevan 
264*755f3dfbSVijay Mahadevan   /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */
265*755f3dfbSVijay Mahadevan   // ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
266b117cd09SVijay Mahadevan 
267*755f3dfbSVijay Mahadevan   ierr = DMGetDimension(dmp, &dim);CHKERRQ(ierr);
268b117cd09SVijay Mahadevan 
269*755f3dfbSVijay Mahadevan   factor = std::pow(2.0 /*degree_P_for_refinement*/,(dmbc->hlevel-dmbp->hlevel)*dmbp->dim*1.0);
27049d66b22SVijay Mahadevan 
271e882eb38SVijay Mahadevan   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
272*755f3dfbSVijay Mahadevan   for(moab::Range::iterator iter = eowned.begin(); iter!= eowned.end(); iter++) {
273b117cd09SVijay Mahadevan 
274b117cd09SVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
275b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> children;
276b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> connp, connc;
277*755f3dfbSVijay Mahadevan     moab::Range connc_owned;
278b117cd09SVijay Mahadevan 
279b117cd09SVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
280*755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr);
281b117cd09SVijay Mahadevan 
282b117cd09SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
283*755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr);
284*755f3dfbSVijay Mahadevan     merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc_owned);MBERRNM(merr);
285*755f3dfbSVijay Mahadevan     merr = dmbc->pcomm->filter_pstatus(connc_owned,PSTATUS_NOT_OWNED,PSTATUS_NOT);MBERRNM(merr);
286*755f3dfbSVijay Mahadevan     for (unsigned tc=0; tc<connc_owned.size(); tc++) {
287*755f3dfbSVijay Mahadevan       connc.push_back(connc_owned[tc]);
288b117cd09SVijay Mahadevan     }
289b117cd09SVijay Mahadevan 
290b117cd09SVijay Mahadevan     std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size());
291b117cd09SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
292*755f3dfbSVijay Mahadevan     merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr);
293*755f3dfbSVijay Mahadevan     merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr);
294b117cd09SVijay Mahadevan 
295b117cd09SVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
296b117cd09SVijay Mahadevan     /* TODO: specific to scalar system - use GetDofs */
297*755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr);
298*755f3dfbSVijay Mahadevan     ierr = DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr);
299b117cd09SVijay Mahadevan 
300e882eb38SVijay Mahadevan     /* Compute the interpolation weights by determining distance of 1-ring
301*755f3dfbSVijay Mahadevan        neighbor vertices from current vertex
302*755f3dfbSVijay Mahadevan 
303*755f3dfbSVijay Mahadevan        TODO: This needs to be replaced with a consistent interface to compute
304*755f3dfbSVijay Mahadevan        the basis function interpolation between the levels evaluated correctly.
305*755f3dfbSVijay Mahadevan        RBF basis will be terrible for any unsmooth problems..
306*755f3dfbSVijay Mahadevan        -- NEEDS TO BE REMOVED --
307*755f3dfbSVijay Mahadevan     */
308*755f3dfbSVijay Mahadevan     values_phi.resize(connp.size());
309*755f3dfbSVijay Mahadevan     for (unsigned tc=0;tc<connc.size(); tc++) {
310*755f3dfbSVijay Mahadevan 
311b117cd09SVijay Mahadevan       double normsum=0.0;
312*755f3dfbSVijay Mahadevan       for (unsigned tp=0;tp<connp.size(); tp++) {
313*755f3dfbSVijay Mahadevan         values_phi[tp] = 0.0;
314e882eb38SVijay Mahadevan         for (unsigned k=0;k<3; k++)
315*755f3dfbSVijay Mahadevan           values_phi[tp] += std::pow(pcoords[tp*3+k]-ccoords[k+tc*3], dim);
316*755f3dfbSVijay Mahadevan         if (values_phi[tp] < 1e-12) {
317*755f3dfbSVijay Mahadevan           values_phi[tp] = 1e12;
318b117cd09SVijay Mahadevan         }
319b117cd09SVijay Mahadevan         else {
320*755f3dfbSVijay Mahadevan           //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim);
321*755f3dfbSVijay Mahadevan           values_phi[tp] = std::pow(values_phi[tp], -1.0);
322*755f3dfbSVijay Mahadevan           normsum += values_phi[tp];
323b117cd09SVijay Mahadevan         }
324b117cd09SVijay Mahadevan       }
325*755f3dfbSVijay Mahadevan       for (unsigned tp=0;tp<connp.size(); tp++) {
326*755f3dfbSVijay Mahadevan         if (values_phi[tp] > 1e11)
327*755f3dfbSVijay Mahadevan           values_phi[tp] = factor*0.5/connp.size();
328b117cd09SVijay Mahadevan         else
329*755f3dfbSVijay Mahadevan           values_phi[tp] = factor*values_phi[tp]*0.5/(connp.size()*normsum);
330b117cd09SVijay Mahadevan       }
331*755f3dfbSVijay Mahadevan       ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
332b117cd09SVijay Mahadevan     }
333b117cd09SVijay Mahadevan 
334b117cd09SVijay Mahadevan     //get interpolation weights
335b117cd09SVijay Mahadevan     //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr);
336e882eb38SVijay Mahadevan     // for (int j=0;j<dofs_per_element; j++)
337e882eb38SVijay Mahadevan     //  std::cout<<"values "<<values_phi[j]<<std::endl;
338b117cd09SVijay Mahadevan 
339b117cd09SVijay Mahadevan   }
340b117cd09SVijay Mahadevan 
341b117cd09SVijay Mahadevan   ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
342b117cd09SVijay Mahadevan   ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
343b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
344b117cd09SVijay Mahadevan }
345b117cd09SVijay Mahadevan 
346b117cd09SVijay Mahadevan #undef __FUNCT__
347b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInjection_Moab"
348b117cd09SVijay Mahadevan /*@
349b117cd09SVijay Mahadevan   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
350b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
351b117cd09SVijay Mahadevan   provided by the user.
352b117cd09SVijay Mahadevan 
353b117cd09SVijay Mahadevan   Collective on MPI_Comm
354b117cd09SVijay Mahadevan 
355b117cd09SVijay Mahadevan   Input Parameter:
356b117cd09SVijay Mahadevan . dmb  - The DMMoab object
357b117cd09SVijay Mahadevan 
358b117cd09SVijay Mahadevan   Output Parameter:
359b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
360b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
361b117cd09SVijay Mahadevan 
362b117cd09SVijay Mahadevan   Level: beginner
363b117cd09SVijay Mahadevan 
364b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
365b117cd09SVijay Mahadevan @*/
366b117cd09SVijay Mahadevan PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx)
367b117cd09SVijay Mahadevan {
368e882eb38SVijay Mahadevan   //DM_Moab        *dmmoab;
369b117cd09SVijay Mahadevan 
370b117cd09SVijay Mahadevan   PetscFunctionBegin;
371b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
372b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
373e882eb38SVijay Mahadevan   //dmmoab = (DM_Moab*)(dm1)->data;
374b117cd09SVijay Mahadevan 
375b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
376b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
377b117cd09SVijay Mahadevan }
378b117cd09SVijay Mahadevan 
379b117cd09SVijay Mahadevan 
380b117cd09SVijay Mahadevan #undef __FUNCT__
381b117cd09SVijay Mahadevan #define __FUNCT__ "DM_UMR_Moab_Private"
382b117cd09SVijay Mahadevan PetscErrorCode  DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref)
383b117cd09SVijay Mahadevan {
384b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
385e882eb38SVijay Mahadevan   PetscInt        i,dim;
386b117cd09SVijay Mahadevan   DM              dm2;
387b117cd09SVijay Mahadevan   moab::ErrorCode merr;
388b117cd09SVijay Mahadevan   DM_Moab        *dmb = (DM_Moab*)dm->data,*dd2;
389b117cd09SVijay Mahadevan 
390b117cd09SVijay Mahadevan   PetscFunctionBegin;
391b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
392b117cd09SVijay Mahadevan   PetscValidPointer(dmref,4);
393b117cd09SVijay Mahadevan 
394b117cd09SVijay Mahadevan   if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) {
395e882eb38SVijay 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);
396e882eb38SVijay Mahadevan     if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1);
397e882eb38SVijay Mahadevan     *dmref = PETSC_NULL;
398b117cd09SVijay Mahadevan     PetscFunctionReturn(0);
399b117cd09SVijay Mahadevan   }
400b117cd09SVijay Mahadevan 
401b117cd09SVijay Mahadevan   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
402b117cd09SVijay Mahadevan   dd2 = (DM_Moab*)dm2->data;
403b117cd09SVijay Mahadevan 
404b117cd09SVijay Mahadevan   dd2->mbiface = dmb->mbiface;
405b117cd09SVijay Mahadevan   dd2->pcomm = dmb->pcomm;
406b117cd09SVijay Mahadevan   dd2->icreatedinstance = PETSC_FALSE;
40764e1c140SVijay Mahadevan   dd2->nghostrings=dmb->nghostrings;
408b117cd09SVijay Mahadevan 
409b117cd09SVijay Mahadevan   /* set the new level based on refinement/coarsening */
410b117cd09SVijay Mahadevan   if (refine) {
411b117cd09SVijay Mahadevan     dd2->hlevel=dmb->hlevel+1;
412b117cd09SVijay Mahadevan   }
413b117cd09SVijay Mahadevan   else {
414b117cd09SVijay Mahadevan     dd2->hlevel=dmb->hlevel-1;
415b117cd09SVijay Mahadevan   }
416b117cd09SVijay Mahadevan 
417b117cd09SVijay Mahadevan   /* Copy the multilevel hierarchy pointers in MOAB */
418b117cd09SVijay Mahadevan   dd2->hierarchy = dmb->hierarchy;
419b117cd09SVijay Mahadevan   dd2->nhlevels = dmb->nhlevels;
420b117cd09SVijay Mahadevan   ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr);
421b117cd09SVijay Mahadevan   for (i=0; i<=dd2->nhlevels; i++) {
422b117cd09SVijay Mahadevan     dd2->hsets[i]=dmb->hsets[i];
423b117cd09SVijay Mahadevan   }
424b117cd09SVijay Mahadevan   dd2->fileset = dd2->hsets[dd2->hlevel];
425b117cd09SVijay Mahadevan 
426b117cd09SVijay Mahadevan   /* do the remaining initializations for DMMoab */
427b117cd09SVijay Mahadevan   dd2->bs = dmb->bs;
428b117cd09SVijay Mahadevan   dd2->numFields = dmb->numFields;
429b117cd09SVijay Mahadevan   dd2->rw_dbglevel = dmb->rw_dbglevel;
430b117cd09SVijay Mahadevan   dd2->partition_by_rank = dmb->partition_by_rank;
431b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
432b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
433b117cd09SVijay Mahadevan   dd2->read_mode = dmb->read_mode;
434b117cd09SVijay Mahadevan   dd2->write_mode = dmb->write_mode;
435b117cd09SVijay Mahadevan 
436b117cd09SVijay Mahadevan   /* set global ID tag handle */
437b117cd09SVijay Mahadevan   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
438b117cd09SVijay Mahadevan 
439b117cd09SVijay Mahadevan   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
440b117cd09SVijay Mahadevan 
441b117cd09SVijay Mahadevan   ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr);
442b117cd09SVijay Mahadevan   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
443b117cd09SVijay Mahadevan   ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr);
444b117cd09SVijay Mahadevan 
445b117cd09SVijay Mahadevan   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
446b117cd09SVijay Mahadevan   dm2->ops->creatematrix = dm->ops->creatematrix;
447b117cd09SVijay Mahadevan 
448b117cd09SVijay Mahadevan   /* copy fill information if given */
449b117cd09SVijay Mahadevan   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
450b117cd09SVijay Mahadevan 
451b117cd09SVijay Mahadevan   /* copy vector type information */
452b117cd09SVijay Mahadevan   ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr);
453b117cd09SVijay Mahadevan   ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr);
4543f1c6e43SVijay Mahadevan   dd2->numFields = dmb->numFields;
4553f1c6e43SVijay Mahadevan   if (dmb->numFields) {
456b117cd09SVijay Mahadevan     ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr);
4573f1c6e43SVijay Mahadevan   }
458b117cd09SVijay Mahadevan 
459b117cd09SVijay Mahadevan   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
460b117cd09SVijay Mahadevan 
461b117cd09SVijay Mahadevan   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
462b117cd09SVijay Mahadevan   ierr = DMSetUp(dm2);CHKERRQ(ierr);
463b117cd09SVijay Mahadevan 
464b117cd09SVijay Mahadevan   *dmref = dm2;
465b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
466b117cd09SVijay Mahadevan }
467b117cd09SVijay Mahadevan 
468b117cd09SVijay Mahadevan 
469b117cd09SVijay Mahadevan #undef __FUNCT__
470b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefine_Moab"
471b117cd09SVijay Mahadevan /*@
472b117cd09SVijay Mahadevan   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
473b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
474b117cd09SVijay Mahadevan   provided by the user.
475b117cd09SVijay Mahadevan 
476e882eb38SVijay Mahadevan   Collective on DM
477b117cd09SVijay Mahadevan 
478b117cd09SVijay Mahadevan   Input Parameter:
479e882eb38SVijay Mahadevan + dm  - The DMMoab object
480e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
481b117cd09SVijay Mahadevan 
482b117cd09SVijay Mahadevan   Output Parameter:
483e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL
484b117cd09SVijay Mahadevan 
485e882eb38SVijay Mahadevan   Note: If no refinement was done, the return value is NULL
486e882eb38SVijay Mahadevan 
487e882eb38SVijay Mahadevan   Level: developer
488b117cd09SVijay Mahadevan 
489b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
490b117cd09SVijay Mahadevan @*/
491b117cd09SVijay Mahadevan PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf)
492b117cd09SVijay Mahadevan {
493b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
494b117cd09SVijay Mahadevan 
495b117cd09SVijay Mahadevan   PetscFunctionBegin;
496b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
497b117cd09SVijay Mahadevan 
498b117cd09SVijay Mahadevan   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr);
499b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
500b117cd09SVijay Mahadevan }
501b117cd09SVijay Mahadevan 
502b117cd09SVijay Mahadevan #undef __FUNCT__
503b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsen_Moab"
504b117cd09SVijay Mahadevan /*@
505b117cd09SVijay Mahadevan   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
506b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
507b117cd09SVijay Mahadevan   provided by the user.
508b117cd09SVijay Mahadevan 
509e882eb38SVijay Mahadevan   Collective on DM
510b117cd09SVijay Mahadevan 
511b117cd09SVijay Mahadevan   Input Parameter:
512e882eb38SVijay Mahadevan . dm  - The DMMoab object
513e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
514b117cd09SVijay Mahadevan 
515b117cd09SVijay Mahadevan   Output Parameter:
516e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL
517b117cd09SVijay Mahadevan 
518e882eb38SVijay Mahadevan   Note: If no coarsening was done, the return value is NULL
519b117cd09SVijay Mahadevan 
520e882eb38SVijay Mahadevan   Level: developer
521e882eb38SVijay Mahadevan 
522e882eb38SVijay Mahadevan .keywords: DMMoab, create, coarsening
523b117cd09SVijay Mahadevan @*/
524b117cd09SVijay Mahadevan PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc)
525b117cd09SVijay Mahadevan {
526b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
527b117cd09SVijay Mahadevan 
528b117cd09SVijay Mahadevan   PetscFunctionBegin;
529b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
530b117cd09SVijay Mahadevan 
531b117cd09SVijay Mahadevan   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr);
532b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
533b117cd09SVijay Mahadevan }
534