xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision 49d66b22a367b93e538da1d54516fdaf6756fee1)
1b117cd09SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2b117cd09SVijay Mahadevan 
3b117cd09SVijay Mahadevan #include <petscdmmoab.h>
4b117cd09SVijay Mahadevan #include <MBTagConventions.hpp>
5b117cd09SVijay Mahadevan #include <moab/NestedRefine.hpp>
6b117cd09SVijay Mahadevan 
7b117cd09SVijay Mahadevan #undef __FUNCT__
8b117cd09SVijay Mahadevan #define __FUNCT__ "DMMoabGenerateHierarchy"
9b117cd09SVijay Mahadevan /*@
10b117cd09SVijay Mahadevan   DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy
11b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
12b117cd09SVijay Mahadevan   provided by the user.
13b117cd09SVijay Mahadevan 
14b117cd09SVijay Mahadevan   Collective on MPI_Comm
15b117cd09SVijay Mahadevan 
16b117cd09SVijay Mahadevan   Input Parameter:
17b117cd09SVijay Mahadevan + dmb  - The DMMoab object
18b117cd09SVijay Mahadevan 
19b117cd09SVijay Mahadevan   Output Parameter:
20b117cd09SVijay Mahadevan + nlevels   - The number of levels of refinement needed to generate the hierarchy
21b117cd09SVijay Mahadevan . ldegrees  - The degree of refinement at each level in the hierarchy
22b117cd09SVijay Mahadevan 
23b117cd09SVijay Mahadevan   Level: beginner
24b117cd09SVijay Mahadevan 
25b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
26b117cd09SVijay Mahadevan @*/
27b117cd09SVijay Mahadevan PetscErrorCode DMMoabGenerateHierarchy(DM dm,PetscInt nlevels,PetscInt *ldegrees)
28b117cd09SVijay Mahadevan {
29b117cd09SVijay Mahadevan   DM_Moab        *dmmoab;
30b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
31b117cd09SVijay Mahadevan   moab::ErrorCode merr;
32b117cd09SVijay Mahadevan   PetscInt *pdegrees,i;
33e882eb38SVijay Mahadevan   std::vector<moab::EntityHandle> hsets;
34b117cd09SVijay Mahadevan 
35b117cd09SVijay Mahadevan   PetscFunctionBegin;
36b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
37b117cd09SVijay Mahadevan   dmmoab = (DM_Moab*)(dm)->data;
38b117cd09SVijay Mahadevan 
39b117cd09SVijay Mahadevan   if (!ldegrees) {
40b117cd09SVijay Mahadevan     ierr = PetscMalloc1(nlevels,&pdegrees);CHKERRQ(ierr);
41b117cd09SVijay Mahadevan     for (i=0; i<nlevels; i++) pdegrees[i]=2; /* default = Degree 2 refinement */
42b117cd09SVijay Mahadevan   }
43b117cd09SVijay Mahadevan   else pdegrees=ldegrees;
44b117cd09SVijay Mahadevan 
45b117cd09SVijay Mahadevan   /* initialize set level refinement data for hierarchy */
46b117cd09SVijay Mahadevan   dmmoab->nhlevels=nlevels;
47b117cd09SVijay Mahadevan 
48b117cd09SVijay Mahadevan   /* Instantiate the nested refinement class */
493f1c6e43SVijay Mahadevan   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset);
50b117cd09SVijay Mahadevan 
51e882eb38SVijay Mahadevan   ierr = PetscMalloc1(nlevels+1,&dmmoab->hsets);CHKERRQ(ierr);
52e882eb38SVijay Mahadevan 
53b117cd09SVijay Mahadevan   /* generate the mesh hierarchy */
54e882eb38SVijay Mahadevan   merr = dmmoab->hierarchy->generate_mesh_hierarchy(pdegrees, nlevels, hsets);MBERRNM(merr);
55e882eb38SVijay Mahadevan 
56*49d66b22SVijay Mahadevan   PetscInfo2(NULL, "Exchanging ghost cells (dim %d) with %d rings\n",dmmoab->dim,dmmoab->nghostrings);
57*49d66b22SVijay Mahadevan   for (i=0; i<=nlevels; i++) {
58*49d66b22SVijay Mahadevan     /* resolve the shared entities by exchanging information to adjacent processors */
59*49d66b22SVijay Mahadevan     merr = dmmoab->pcomm->exchange_ghost_cells(dmmoab->dim,0,dmmoab->nghostrings,dmmoab->dim,true,false,&hsets[i]);MBERRV(dmmoab->mbiface,merr);
60*49d66b22SVijay Mahadevan 
61*49d66b22SVijay Mahadevan     dmmoab->hsets[i]=hsets[i];
62*49d66b22SVijay Mahadevan   }
63b117cd09SVijay Mahadevan 
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 @*/
158b117cd09SVijay Mahadevan PetscErrorCode DMCreateInterpolation_Moab(DM dm1,DM dm2,Mat* interpl,Vec* vec)
159b117cd09SVijay Mahadevan {
160b117cd09SVijay Mahadevan   DM_Moab         *dmb1, *dmb2;
161b117cd09SVijay Mahadevan   PetscErrorCode   ierr;
162b117cd09SVijay Mahadevan   moab::ErrorCode  merr;
163e882eb38SVijay Mahadevan   PetscInt         dim;
1643f1c6e43SVijay Mahadevan   PetscReal        factor;
1653f1c6e43SVijay Mahadevan   PetscBool        eonbnd;
166ce27a4eeSVijay Mahadevan   PetscInt         innz, *nnz, ionz, *onz;
167*49d66b22SVijay Mahadevan   PetscInt         nlsiz1, nlsiz2, nlghsiz1, nlghsiz2, ngsiz1, ngsiz2;
168b117cd09SVijay Mahadevan   std::vector<int> bndrows;
1693f1c6e43SVijay Mahadevan   std::vector<PetscBool> dbdry;
170b117cd09SVijay Mahadevan 
171b117cd09SVijay Mahadevan   PetscFunctionBegin;
172b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
173b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
174b117cd09SVijay Mahadevan   dmb1 = (DM_Moab*)(dm1)->data;
175b117cd09SVijay Mahadevan   dmb2 = (DM_Moab*)(dm2)->data;
176*49d66b22SVijay Mahadevan   nlsiz1 = dmb1->nloc*dmb1->numFields;
177941e0cffSVijay Mahadevan   nlsiz2 = dmb2->nloc*dmb2->numFields;
178*49d66b22SVijay Mahadevan   ngsiz1 = dmb1->n*dmb1->numFields;
179941e0cffSVijay Mahadevan   ngsiz2 = dmb2->n*dmb2->numFields;
180*49d66b22SVijay Mahadevan   nlghsiz1 = (dmb1->nloc+dmb1->nghost)*dmb1->numFields;
181*49d66b22SVijay Mahadevan   nlghsiz2 = (dmb2->nloc+dmb2->nghost)*dmb2->numFields;
182b117cd09SVijay Mahadevan 
1833f1c6e43SVijay Mahadevan   PetscInfo4(dm1,"Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n",dmb1->nloc,dmb2->nloc,dmb1->hlevel,dmb2->hlevel);
1843f1c6e43SVijay Mahadevan 
185b117cd09SVijay Mahadevan   ierr = MatCreate(PetscObjectComm((PetscObject)dm1), interpl);CHKERRQ(ierr);
186e882eb38SVijay Mahadevan   ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr);
187941e0cffSVijay Mahadevan   ierr = MatSetSizes(*interpl, nlsiz2, nlsiz1, ngsiz2, ngsiz1);CHKERRQ(ierr);
188b117cd09SVijay Mahadevan 
189b117cd09SVijay Mahadevan   /* TODO: This is a hack for the rectangular system - decipher NNZ pattern better */
190b117cd09SVijay Mahadevan   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
191b117cd09SVijay Mahadevan 
192941e0cffSVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
193*49d66b22SVijay Mahadevan   ierr = PetscCalloc2(nlghsiz2,&nnz,nlghsiz2,&onz);CHKERRQ(ierr);
194941e0cffSVijay Mahadevan 
195941e0cffSVijay Mahadevan   /* Loop through the local elements and compute the relation between the current parent and the refined_level. */
196941e0cffSVijay Mahadevan   for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) {
197941e0cffSVijay Mahadevan 
198941e0cffSVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
199941e0cffSVijay Mahadevan     std::vector<moab::EntityHandle> children;
200941e0cffSVijay Mahadevan     std::vector<moab::EntityHandle> connp, connc;
201941e0cffSVijay Mahadevan 
202941e0cffSVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
203941e0cffSVijay Mahadevan     merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr);
204941e0cffSVijay Mahadevan 
205941e0cffSVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
206941e0cffSVijay Mahadevan     merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr);
207941e0cffSVijay Mahadevan     for (unsigned ic=0; ic < children.size(); ic++) {
208941e0cffSVijay Mahadevan       std::vector<moab::EntityHandle> tconnc;
209941e0cffSVijay Mahadevan       /* Get coordinates of the parent vertices in canonical order */
210941e0cffSVijay Mahadevan       merr = dmb2->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr);
211941e0cffSVijay Mahadevan       for (unsigned tc=0; tc<tconnc.size(); tc++) {
212941e0cffSVijay Mahadevan         if (std::find(connc.begin(), connc.end(), tconnc[tc]) == connc.end())
213941e0cffSVijay Mahadevan           connc.push_back(tconnc[tc]);
214941e0cffSVijay Mahadevan       }
215941e0cffSVijay Mahadevan     }
216941e0cffSVijay Mahadevan 
217941e0cffSVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
218941e0cffSVijay Mahadevan     /* TODO: specific to scalar system - use GetDofs */
219941e0cffSVijay Mahadevan     //ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr);
220*49d66b22SVijay Mahadevan     ierr = DMMoabGetFieldDofsLocal(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr);
221941e0cffSVijay Mahadevan 
222941e0cffSVijay Mahadevan     for (unsigned tp=0;tp<connp.size(); tp++) {
223941e0cffSVijay Mahadevan       // ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
224941e0cffSVijay Mahadevan       if (dmb1->vowned->find(connp[tp]) != dmb1->vowned->end()) {
225941e0cffSVijay Mahadevan         for (unsigned tc=0;tc<connc.size(); tc++) nnz[dofsc[tc]]++;
226941e0cffSVijay Mahadevan       }
227941e0cffSVijay Mahadevan       else if (dmb1->vghost->find(connp[tp]) != dmb1->vghost->end()) {
228941e0cffSVijay Mahadevan         for (unsigned tc=0;tc<connc.size(); tc++) onz[dofsc[tc]]++;
229941e0cffSVijay Mahadevan       }
230941e0cffSVijay Mahadevan       else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in parent level %D\n", connc[tp]);
231941e0cffSVijay Mahadevan     }
232941e0cffSVijay Mahadevan     for(int tc = 0; tc < connc.size(); tc++) {
233941e0cffSVijay Mahadevan       if (dmb2->vowned->find(connc[tc]) != dmb2->vowned->end()) nnz[dofsc[tc]]++;
234941e0cffSVijay Mahadevan       else if (dmb2->vghost->find(connc[tc]) != dmb2->vghost->end()) onz[dofsc[tc]]++;
235941e0cffSVijay Mahadevan       else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in child level %D\n", connc[tc]);
236*49d66b22SVijay Mahadevan     }
237941e0cffSVijay Mahadevan   }
238941e0cffSVijay Mahadevan 
239*49d66b22SVijay Mahadevan   /*
240941e0cffSVijay Mahadevan   int i=0;
241941e0cffSVijay Mahadevan   std::vector<moab::EntityHandle> adjs;
242941e0cffSVijay Mahadevan   for(moab::Range::iterator iter = dmb1->vowned->begin(); iter!= dmb1->vowned->end(); iter++, i++) {
243941e0cffSVijay Mahadevan     merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr);
244941e0cffSVijay Mahadevan     nnz[i] -= adjs.size();
245941e0cffSVijay Mahadevan     adjs.clear();
246941e0cffSVijay Mahadevan   }
247*49d66b22SVijay Mahadevan   for(moab::Range::iterator iter = dmb1->vghost->begin(); iter!= dmb1->vghost->end(); iter++, i++) {
248941e0cffSVijay Mahadevan     //merr = dmb1->hierarchy->get_adjacencies(&(*iter), 1, 0, false, adjs, moab::Interface::UNION);MBERRNM(merr);
249941e0cffSVijay Mahadevan     merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr);
250941e0cffSVijay Mahadevan     onz[i] -= adjs.size();
251941e0cffSVijay Mahadevan     adjs.clear();
252941e0cffSVijay Mahadevan   }
253*49d66b22SVijay Mahadevan   */
254941e0cffSVijay Mahadevan 
255941e0cffSVijay Mahadevan   ionz=onz[0];
256941e0cffSVijay Mahadevan   innz=nnz[0];
257941e0cffSVijay Mahadevan   for (int tc=0; tc < nlsiz2; tc++) {
258941e0cffSVijay Mahadevan     // check for maximum allowed sparsity = fully dense
259941e0cffSVijay Mahadevan     nnz[tc] = std::min(nlsiz1,nnz[tc]);
260941e0cffSVijay Mahadevan     onz[tc] = std::min(nlsiz1,onz[tc]);
261941e0cffSVijay Mahadevan 
262941e0cffSVijay Mahadevan     innz = (innz < nnz[tc] ? nnz[tc] : innz);
263941e0cffSVijay Mahadevan     ionz = (ionz < onz[tc] ? onz[tc] : ionz);
264ce27a4eeSVijay Mahadevan     //PetscPrintf(PETSC_COMM_SELF, "[%D] NNZ = %D, ONZ = %D\n", tc, nnz[tc], onz[tc]);
265941e0cffSVijay Mahadevan   }
266ce27a4eeSVijay Mahadevan   //PetscPrintf(PETSC_COMM_SELF, "Final: INNZ = %D, IONZ = %D\n", innz, ionz);
267941e0cffSVijay Mahadevan 
268941e0cffSVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(*interpl,innz,nnz);CHKERRQ(ierr);
269941e0cffSVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(*interpl,innz,nnz,ionz,onz);CHKERRQ(ierr);
270941e0cffSVijay Mahadevan 
271941e0cffSVijay Mahadevan   /* clean up temporary memory */
272941e0cffSVijay Mahadevan   ierr = PetscFree2(nnz,onz);CHKERRQ(ierr);
273b117cd09SVijay Mahadevan 
274b117cd09SVijay Mahadevan   /* set up internal matrix data-structures */
275b117cd09SVijay Mahadevan   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
276b117cd09SVijay Mahadevan   ierr = MatZeroEntries(*interpl);CHKERRQ(ierr);
277b117cd09SVijay Mahadevan 
278b117cd09SVijay Mahadevan   ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr);
279b117cd09SVijay Mahadevan 
2803f1c6e43SVijay Mahadevan   factor = std::pow(2.0 /*degree_P_for_refinement*/,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0);
281b117cd09SVijay Mahadevan 
282*49d66b22SVijay Mahadevan   ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
283*49d66b22SVijay Mahadevan 
284e882eb38SVijay Mahadevan   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
285b117cd09SVijay Mahadevan   for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) {
286b117cd09SVijay Mahadevan 
287b117cd09SVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
288b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> children;
289b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> connp, connc;
290b117cd09SVijay Mahadevan 
291b117cd09SVijay Mahadevan     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
292b117cd09SVijay Mahadevan     merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr);
293b117cd09SVijay Mahadevan 
294b117cd09SVijay Mahadevan     /* Get connectivity and coordinates of the parent vertices */
295b117cd09SVijay Mahadevan     merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr);
296e882eb38SVijay Mahadevan     for (unsigned ic=0; ic < children.size(); ic++) {
297b117cd09SVijay Mahadevan       std::vector<moab::EntityHandle> tconnc;
298b117cd09SVijay Mahadevan       /* Get coordinates of the parent vertices in canonical order */
299b117cd09SVijay Mahadevan       merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr);
300e882eb38SVijay Mahadevan       for (unsigned tc=0; tc<tconnc.size(); tc++) {
301b117cd09SVijay Mahadevan         connc.push_back(tconnc[tc]);
302b117cd09SVijay Mahadevan       }
303b117cd09SVijay Mahadevan     }
304b117cd09SVijay Mahadevan 
305b117cd09SVijay Mahadevan     std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size());
306b117cd09SVijay Mahadevan     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
307b117cd09SVijay Mahadevan     merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr);
308b117cd09SVijay Mahadevan     merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr);
309b117cd09SVijay Mahadevan 
310b117cd09SVijay Mahadevan     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
311b117cd09SVijay Mahadevan     /* TODO: specific to scalar system - use GetDofs */
312*49d66b22SVijay Mahadevan     ierr = DMMoabGetFieldDofsLocal(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr);
313*49d66b22SVijay Mahadevan     ierr = DMMoabGetFieldDofsLocal(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr);
314b117cd09SVijay Mahadevan 
315e882eb38SVijay Mahadevan     /* Compute the interpolation weights by determining distance of 1-ring
316e882eb38SVijay Mahadevan        neighbor vertices from current vertex */
317e882eb38SVijay Mahadevan     for (unsigned i=0;i<connp.size(); i++) {
318b117cd09SVijay Mahadevan       double normsum=0.0;
319e882eb38SVijay Mahadevan       for (unsigned j=0;j<connc.size(); j++) {
3203f1c6e43SVijay Mahadevan         values_phi[j] = 0.0;
321e882eb38SVijay Mahadevan         for (unsigned k=0;k<3; k++)
322941e0cffSVijay Mahadevan           values_phi[j] += std::pow(pcoords[i*3+k]-ccoords[k+j*3], dim);
3233f1c6e43SVijay Mahadevan         if (values_phi[j] < 1e-12) {
3243f1c6e43SVijay Mahadevan           values_phi[j] = 1e12;
325b117cd09SVijay Mahadevan         }
326b117cd09SVijay Mahadevan         else {
327941e0cffSVijay Mahadevan           //values_phi[j] = std::pow(values_phi[j], -1.0/dim);
328941e0cffSVijay Mahadevan           values_phi[j] = std::pow(values_phi[j], -1.0);
3293f1c6e43SVijay Mahadevan           normsum += values_phi[j];
330b117cd09SVijay Mahadevan         }
331b117cd09SVijay Mahadevan       }
332e882eb38SVijay Mahadevan       for (unsigned j=0;j<connc.size(); j++) {
3333f1c6e43SVijay Mahadevan         if (values_phi[j] > 1e11)
3343f1c6e43SVijay Mahadevan           values_phi[j] = factor*0.5/connc.size();
335b117cd09SVijay Mahadevan         else
3363f1c6e43SVijay Mahadevan           values_phi[j] = factor*values_phi[j]*0.5/(connc.size()*normsum);
337b117cd09SVijay Mahadevan       }
338b117cd09SVijay Mahadevan       ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
339b117cd09SVijay Mahadevan     }
340b117cd09SVijay Mahadevan 
341b117cd09SVijay Mahadevan     /* check if element is on the boundary */
342b117cd09SVijay Mahadevan     //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr);
3433f1c6e43SVijay Mahadevan     dbdry.resize(connc.size());
3443f1c6e43SVijay Mahadevan     ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry.data());CHKERRQ(ierr);
345b117cd09SVijay Mahadevan     eonbnd=PETSC_FALSE;
346e882eb38SVijay Mahadevan     for (unsigned i=0; i< connc.size(); ++i)
347b117cd09SVijay Mahadevan       if (dbdry[i]) eonbnd=PETSC_TRUE;
348b117cd09SVijay Mahadevan 
349b117cd09SVijay Mahadevan     values_phi.clear();
350b117cd09SVijay Mahadevan     values_phi.resize(connp.size());
351b117cd09SVijay Mahadevan     /* apply dirichlet boundary conditions */
352b117cd09SVijay Mahadevan     if (eonbnd) {
353b117cd09SVijay Mahadevan 
354b117cd09SVijay Mahadevan       ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
355b117cd09SVijay Mahadevan       ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
356b117cd09SVijay Mahadevan       /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */
357b117cd09SVijay Mahadevan       //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr);
358e882eb38SVijay Mahadevan       for (unsigned i=0; i < connc.size(); i++) {
3593f1c6e43SVijay Mahadevan         if (dbdry[i]) {  /* dirichlet node */
360b117cd09SVijay Mahadevan           /* think about strongly imposing dirichlet */
361b117cd09SVijay Mahadevan           //bndrows.push_back(dofsc[i]);
362b117cd09SVijay Mahadevan 
363b117cd09SVijay Mahadevan           ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
364b117cd09SVijay Mahadevan           //values_phi[0]=1.0;
365b117cd09SVijay Mahadevan           //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
366b117cd09SVijay Mahadevan         }
367b117cd09SVijay Mahadevan       }
368b117cd09SVijay Mahadevan 
369b117cd09SVijay Mahadevan       ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
370b117cd09SVijay Mahadevan       ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
371b117cd09SVijay Mahadevan     }
372b117cd09SVijay Mahadevan 
373b117cd09SVijay Mahadevan     //get interpolation weights
374b117cd09SVijay Mahadevan     //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr);
375e882eb38SVijay Mahadevan     // for (int j=0;j<dofs_per_element; j++)
376e882eb38SVijay Mahadevan     //  std::cout<<"values "<<values_phi[j]<<std::endl;
377b117cd09SVijay Mahadevan 
378b117cd09SVijay Mahadevan     //get row and column indices, zero weights are ignored
379b117cd09SVijay Mahadevan     /*
380b117cd09SVijay Mahadevan     int nz_ind = 0;
381b117cd09SVijay Mahadevan     idx = dmb2->vowned->index(vhandle);
382b117cd09SVijay Mahadevan     for (int j=0;j<dofs_per_element; j++){
383b117cd09SVijay Mahadevan       idy[nz_ind] = dmb1->vowned->index(connectivity[j]);
384b117cd09SVijay Mahadevan       PetscPrintf(PETSC_COMM_WORLD, "Finding coarse connectivity vertex %D associated with [%D, %D] - set to %D\n", connectivity[j], parent.size(), vhandle, idy[nz_ind]);
385b117cd09SVijay Mahadevan       //values_phi[nz_ind] = values_phi[j];
386b117cd09SVijay Mahadevan       nz_ind = nz_ind+1;
387b117cd09SVijay Mahadevan     }
388b117cd09SVijay Mahadevan     */
389b117cd09SVijay Mahadevan 
390b117cd09SVijay Mahadevan     //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr);
391b117cd09SVijay Mahadevan     //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
392b117cd09SVijay Mahadevan   }
393b117cd09SVijay Mahadevan 
394b117cd09SVijay Mahadevan   //PetscPrintf(PETSC_COMM_WORLD, "[Boundary vertices = %D] :: A few: %D %D %D %D \n", bndrows.size(), bndrows[0], bndrows[1], bndrows[2], bndrows[3]);
395b117cd09SVijay Mahadevan   //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr);
396b117cd09SVijay Mahadevan 
397b117cd09SVijay Mahadevan   ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
398b117cd09SVijay Mahadevan   ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
399b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
400b117cd09SVijay Mahadevan }
401b117cd09SVijay Mahadevan 
402b117cd09SVijay Mahadevan #undef __FUNCT__
403b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInjection_Moab"
404b117cd09SVijay Mahadevan /*@
405b117cd09SVijay Mahadevan   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
406b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
407b117cd09SVijay Mahadevan   provided by the user.
408b117cd09SVijay Mahadevan 
409b117cd09SVijay Mahadevan   Collective on MPI_Comm
410b117cd09SVijay Mahadevan 
411b117cd09SVijay Mahadevan   Input Parameter:
412b117cd09SVijay Mahadevan . dmb  - The DMMoab object
413b117cd09SVijay Mahadevan 
414b117cd09SVijay Mahadevan   Output Parameter:
415b117cd09SVijay Mahadevan . nlevels   - The number of levels of refinement needed to generate the hierarchy
416b117cd09SVijay Mahadevan + ldegrees  - The degree of refinement at each level in the hierarchy
417b117cd09SVijay Mahadevan 
418b117cd09SVijay Mahadevan   Level: beginner
419b117cd09SVijay Mahadevan 
420b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
421b117cd09SVijay Mahadevan @*/
422b117cd09SVijay Mahadevan PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx)
423b117cd09SVijay Mahadevan {
424e882eb38SVijay Mahadevan   //DM_Moab        *dmmoab;
425b117cd09SVijay Mahadevan 
426b117cd09SVijay Mahadevan   PetscFunctionBegin;
427b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
428b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
429e882eb38SVijay Mahadevan   //dmmoab = (DM_Moab*)(dm1)->data;
430b117cd09SVijay Mahadevan 
431b117cd09SVijay Mahadevan   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
432b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
433b117cd09SVijay Mahadevan }
434b117cd09SVijay Mahadevan 
435b117cd09SVijay Mahadevan 
436b117cd09SVijay Mahadevan #undef __FUNCT__
437b117cd09SVijay Mahadevan #define __FUNCT__ "DM_UMR_Moab_Private"
438b117cd09SVijay Mahadevan PetscErrorCode  DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref)
439b117cd09SVijay Mahadevan {
440b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
441e882eb38SVijay Mahadevan   PetscInt        i,dim;
442b117cd09SVijay Mahadevan   DM              dm2;
443b117cd09SVijay Mahadevan   moab::ErrorCode merr;
444b117cd09SVijay Mahadevan   DM_Moab        *dmb = (DM_Moab*)dm->data,*dd2;
445b117cd09SVijay Mahadevan 
446b117cd09SVijay Mahadevan   PetscFunctionBegin;
447b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
448b117cd09SVijay Mahadevan   PetscValidPointer(dmref,4);
449b117cd09SVijay Mahadevan 
450b117cd09SVijay Mahadevan   if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) {
451e882eb38SVijay 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);
452e882eb38SVijay Mahadevan     if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1);
453e882eb38SVijay Mahadevan     *dmref = PETSC_NULL;
454b117cd09SVijay Mahadevan     PetscFunctionReturn(0);
455b117cd09SVijay Mahadevan   }
456b117cd09SVijay Mahadevan 
457b117cd09SVijay Mahadevan   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
458b117cd09SVijay Mahadevan   dd2 = (DM_Moab*)dm2->data;
459b117cd09SVijay Mahadevan 
460b117cd09SVijay Mahadevan   dd2->mbiface = dmb->mbiface;
461b117cd09SVijay Mahadevan   dd2->pcomm = dmb->pcomm;
462b117cd09SVijay Mahadevan   dd2->icreatedinstance = PETSC_FALSE;
463b117cd09SVijay Mahadevan 
464b117cd09SVijay Mahadevan   /* set the new level based on refinement/coarsening */
465b117cd09SVijay Mahadevan   if (refine) {
466b117cd09SVijay Mahadevan     dd2->hlevel=dmb->hlevel+1;
467b117cd09SVijay Mahadevan   }
468b117cd09SVijay Mahadevan   else {
469b117cd09SVijay Mahadevan     dd2->hlevel=dmb->hlevel-1;
470b117cd09SVijay Mahadevan   }
471b117cd09SVijay Mahadevan 
472b117cd09SVijay Mahadevan   /* Copy the multilevel hierarchy pointers in MOAB */
473b117cd09SVijay Mahadevan   dd2->hierarchy = dmb->hierarchy;
474b117cd09SVijay Mahadevan   dd2->nhlevels = dmb->nhlevels;
475b117cd09SVijay Mahadevan   ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr);
476b117cd09SVijay Mahadevan   for (i=0; i<=dd2->nhlevels; i++) {
477b117cd09SVijay Mahadevan     dd2->hsets[i]=dmb->hsets[i];
478b117cd09SVijay Mahadevan   }
479b117cd09SVijay Mahadevan   dd2->fileset = dd2->hsets[dd2->hlevel];
480b117cd09SVijay Mahadevan 
481b117cd09SVijay Mahadevan   /* do the remaining initializations for DMMoab */
482b117cd09SVijay Mahadevan   dd2->bs = dmb->bs;
483b117cd09SVijay Mahadevan   dd2->numFields = dmb->numFields;
484b117cd09SVijay Mahadevan   dd2->rw_dbglevel = dmb->rw_dbglevel;
485b117cd09SVijay Mahadevan   dd2->partition_by_rank = dmb->partition_by_rank;
486b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
487b117cd09SVijay Mahadevan   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
488b117cd09SVijay Mahadevan   dd2->read_mode = dmb->read_mode;
489b117cd09SVijay Mahadevan   dd2->write_mode = dmb->write_mode;
490b117cd09SVijay Mahadevan 
491b117cd09SVijay Mahadevan   /* set global ID tag handle */
492b117cd09SVijay Mahadevan   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
493b117cd09SVijay Mahadevan 
494b117cd09SVijay Mahadevan   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
495b117cd09SVijay Mahadevan 
496b117cd09SVijay Mahadevan   ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr);
497b117cd09SVijay Mahadevan   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
498b117cd09SVijay Mahadevan   ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr);
499b117cd09SVijay Mahadevan 
500b117cd09SVijay Mahadevan   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
501b117cd09SVijay Mahadevan   dm2->ops->creatematrix = dm->ops->creatematrix;
502b117cd09SVijay Mahadevan 
503b117cd09SVijay Mahadevan   /* copy fill information if given */
504b117cd09SVijay Mahadevan   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
505b117cd09SVijay Mahadevan 
506b117cd09SVijay Mahadevan   /* copy vector type information */
507b117cd09SVijay Mahadevan   ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr);
508b117cd09SVijay Mahadevan   ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr);
5093f1c6e43SVijay Mahadevan   dd2->numFields = dmb->numFields;
5103f1c6e43SVijay Mahadevan   if (dmb->numFields) {
511b117cd09SVijay Mahadevan     ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr);
5123f1c6e43SVijay Mahadevan   }
513b117cd09SVijay Mahadevan 
514b117cd09SVijay Mahadevan   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
515b117cd09SVijay Mahadevan 
516b117cd09SVijay Mahadevan   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
517b117cd09SVijay Mahadevan   ierr = DMSetUp(dm2);CHKERRQ(ierr);
518b117cd09SVijay Mahadevan 
519b117cd09SVijay Mahadevan   *dmref = dm2;
520b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
521b117cd09SVijay Mahadevan }
522b117cd09SVijay Mahadevan 
523b117cd09SVijay Mahadevan 
524b117cd09SVijay Mahadevan #undef __FUNCT__
525b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefine_Moab"
526b117cd09SVijay Mahadevan /*@
527b117cd09SVijay Mahadevan   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
528b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
529b117cd09SVijay Mahadevan   provided by the user.
530b117cd09SVijay Mahadevan 
531e882eb38SVijay Mahadevan   Collective on DM
532b117cd09SVijay Mahadevan 
533b117cd09SVijay Mahadevan   Input Parameter:
534e882eb38SVijay Mahadevan + dm  - The DMMoab object
535e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
536b117cd09SVijay Mahadevan 
537b117cd09SVijay Mahadevan   Output Parameter:
538e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL
539b117cd09SVijay Mahadevan 
540e882eb38SVijay Mahadevan   Note: If no refinement was done, the return value is NULL
541e882eb38SVijay Mahadevan 
542e882eb38SVijay Mahadevan   Level: developer
543b117cd09SVijay Mahadevan 
544b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement
545b117cd09SVijay Mahadevan @*/
546b117cd09SVijay Mahadevan PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf)
547b117cd09SVijay Mahadevan {
548b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
549b117cd09SVijay Mahadevan 
550b117cd09SVijay Mahadevan   PetscFunctionBegin;
551b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
552b117cd09SVijay Mahadevan 
553b117cd09SVijay Mahadevan   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr);
554b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
555b117cd09SVijay Mahadevan }
556b117cd09SVijay Mahadevan 
557b117cd09SVijay Mahadevan #undef __FUNCT__
558b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsen_Moab"
559b117cd09SVijay Mahadevan /*@
560b117cd09SVijay Mahadevan   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
561b117cd09SVijay Mahadevan   by succesively refining a coarse mesh, already defined in the DM object
562b117cd09SVijay Mahadevan   provided by the user.
563b117cd09SVijay Mahadevan 
564e882eb38SVijay Mahadevan   Collective on DM
565b117cd09SVijay Mahadevan 
566b117cd09SVijay Mahadevan   Input Parameter:
567e882eb38SVijay Mahadevan . dm  - The DMMoab object
568e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
569b117cd09SVijay Mahadevan 
570b117cd09SVijay Mahadevan   Output Parameter:
571e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL
572b117cd09SVijay Mahadevan 
573e882eb38SVijay Mahadevan   Note: If no coarsening was done, the return value is NULL
574b117cd09SVijay Mahadevan 
575e882eb38SVijay Mahadevan   Level: developer
576e882eb38SVijay Mahadevan 
577e882eb38SVijay Mahadevan .keywords: DMMoab, create, coarsening
578b117cd09SVijay Mahadevan @*/
579b117cd09SVijay Mahadevan PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc)
580b117cd09SVijay Mahadevan {
581b117cd09SVijay Mahadevan   PetscErrorCode  ierr;
582b117cd09SVijay Mahadevan 
583b117cd09SVijay Mahadevan   PetscFunctionBegin;
584b117cd09SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
585b117cd09SVijay Mahadevan 
586b117cd09SVijay Mahadevan   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr);
587b117cd09SVijay Mahadevan   PetscFunctionReturn(0);
588b117cd09SVijay Mahadevan }
589