xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision 64e1c140f280b9f14e7e1bfe948bce0ef936fcec)
1 #include <petsc-private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2 
3 #include <petscdmmoab.h>
4 #include <MBTagConventions.hpp>
5 #include <moab/NestedRefine.hpp>
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "DMMoabGenerateHierarchy"
9 /*@
10   DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy
11   by succesively refining a coarse mesh, already defined in the DM object
12   provided by the user.
13 
14   Collective on MPI_Comm
15 
16   Input Parameter:
17 + dmb  - The DMMoab object
18 
19   Output Parameter:
20 + nlevels   - The number of levels of refinement needed to generate the hierarchy
21 . ldegrees  - The degree of refinement at each level in the hierarchy
22 
23   Level: beginner
24 
25 .keywords: DMMoab, create, refinement
26 @*/
27 PetscErrorCode DMMoabGenerateHierarchy(DM dm,PetscInt nlevels,PetscInt *ldegrees)
28 {
29   DM_Moab        *dmmoab;
30   PetscErrorCode  ierr;
31   moab::ErrorCode merr;
32   PetscInt *pdegrees,i;
33   std::vector<moab::EntityHandle> hsets;
34 
35   PetscFunctionBegin;
36   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
37   dmmoab = (DM_Moab*)(dm)->data;
38 
39   if (!ldegrees) {
40     ierr = PetscMalloc1(nlevels,&pdegrees);CHKERRQ(ierr);
41     for (i=0; i<nlevels; i++) pdegrees[i]=2; /* default = Degree 2 refinement */
42   }
43   else pdegrees=ldegrees;
44 
45   /* initialize set level refinement data for hierarchy */
46   dmmoab->nhlevels=nlevels;
47 
48   /* Instantiate the nested refinement class */
49   dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset);
50 
51   ierr = PetscMalloc1(nlevels+1,&dmmoab->hsets);CHKERRQ(ierr);
52   hsets.resize(nlevels+1);
53 
54   /* generate the mesh hierarchy */
55   merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets);MBERRNM(merr);
56 
57   merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);MBERRNM(merr);
58 
59   /* copy the mesh sets for nested refinement hierarchy */
60   for (i=0; i<=nlevels; i++)
61       dmmoab->hsets[i]=hsets[i];
62 
63   if (dmmoab->nghostrings && false) {
64     PetscInfo2(NULL, "Exchanging ghost cells (dim %d) with %d rings\n",dmmoab->dim,dmmoab->nghostrings);
65     // for (i=1; i<=nlevels; i++) {
66     //   /* resolve the shared entities by exchanging information to adjacent processors */
67     //   merr = dmmoab->pcomm->exchange_ghost_cells(dmmoab->dim,0,dmmoab->nghostrings,dmmoab->dim,true,false,&hsets[i]);MBERRV(dmmoab->mbiface,merr);
68     // }
69     merr = dmmoab->pcomm->exchange_ghost_cells(dmmoab->dim,0,dmmoab->nghostrings,dmmoab->dim,true,false);MBERRV(dmmoab->mbiface,merr);
70 
71     // moab::Range vtxall, elmsall;
72     // merr = dmmoab->mbiface->get_entities_by_dimension(0, 0, vtxall, true);MBERRNM(merr);
73     // merr = dmmoab->mbiface->get_entities_by_dimension(0, dmmoab->dim, elmsall, true);MBERRNM(merr);
74 
75   }
76 
77   hsets.clear();
78   if (!ldegrees) {
79     ierr = PetscFree(pdegrees);CHKERRQ(ierr);
80   }
81   PetscFunctionReturn(0);
82 }
83 
84 #undef __FUNCT__
85 #define __FUNCT__ "DMRefineHierarchy_Moab"
86 /*@
87   DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy
88   by succesively refining a coarse mesh.
89 
90   Collective on MPI_Comm
91 
92   Input Parameter:
93 + dm  - The DMMoab object
94 
95   Output Parameter:
96 + nlevels   - The number of levels of refinement needed to generate the hierarchy
97 . dmf  - The DM objects after successive refinement of the hierarchy
98 
99   Level: beginner
100 
101 .keywords: DMMoab, generate, refinement
102 @*/
103 PetscErrorCode  DMRefineHierarchy_Moab(DM dm,PetscInt nlevels,DM dmf[])
104 {
105   PetscErrorCode  ierr;
106   PetscInt        i;
107 
108   PetscFunctionBegin;
109 
110   ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr);
111   for (i=1; i<nlevels; i++) {
112     ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr);
113   }
114   PetscFunctionReturn(0);
115 }
116 
117 #undef __FUNCT__
118 #define __FUNCT__ "DMCoarsenHierarchy_Moab"
119 /*@
120   DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy
121   by succesively coarsening a refined mesh.
122 
123   Collective on MPI_Comm
124 
125   Input Parameter:
126 + dm  - The DMMoab object
127 
128   Output Parameter:
129 + nlevels   - The number of levels of refinement needed to generate the hierarchy
130 . dmc  - The DM objects after successive coarsening of the hierarchy
131 
132   Level: beginner
133 
134 .keywords: DMMoab, generate, coarsening
135 @*/
136 PetscErrorCode DMCoarsenHierarchy_Moab(DM dm,PetscInt nlevels,DM dmc[])
137 {
138   PetscErrorCode  ierr;
139   PetscInt        i;
140 
141   PetscFunctionBegin;
142 
143   ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr);
144   for (i=1; i<nlevels; i++) {
145     ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr);
146   }
147   PetscFunctionReturn(0);
148 }
149 
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "DMCreateInterpolation_Moab"
153 /*@
154   DMCreateInterpolation_Moab - Generate the interpolation operators to transform
155   operators (matrices, vectors) from parent level to child level as defined by
156   the DM inputs provided by the user.
157 
158   Collective on MPI_Comm
159 
160   Input Parameter:
161 + dm1  - The DMMoab object
162 - dm2  - the second, finer DMMoab object
163 
164   Output Parameter:
165 + interpl  - The interpolation operator for transferring data between the levels
166 - vec      - The scaling vector (optional)
167 
168   Level: developer
169 
170 .keywords: DMMoab, create, refinement
171 @*/
172 PetscErrorCode DMCreateInterpolation_Moab(DM dm1,DM dm2,Mat* interpl,Vec* vec)
173 {
174   DM_Moab         *dmb1, *dmb2;
175   PetscErrorCode   ierr;
176   moab::ErrorCode  merr;
177   PetscInt         dim;
178   PetscReal        factor;
179   PetscBool        eonbnd;
180   PetscInt         innz, *nnz, ionz, *onz;
181   PetscInt         nlsiz1, nlsiz2, nlghsiz1, nlghsiz2, ngsiz1, ngsiz2;
182   std::vector<int> bndrows;
183   std::vector<PetscBool> dbdry;
184 
185   PetscFunctionBegin;
186   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
187   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
188   dmb1 = (DM_Moab*)(dm1)->data;
189   dmb2 = (DM_Moab*)(dm2)->data;
190   nlsiz1 = dmb1->nloc*dmb1->numFields;
191   nlsiz2 = dmb2->nloc*dmb2->numFields;
192   ngsiz1 = dmb1->n*dmb1->numFields;
193   ngsiz2 = dmb2->n*dmb2->numFields;
194   nlghsiz1 = (dmb1->nloc+dmb1->nghost)*dmb1->numFields;
195   nlghsiz2 = (dmb2->nloc+dmb2->nghost)*dmb2->numFields;
196 
197   int rank = dmb1->pcomm->rank();
198 
199   PetscInfo4(dm1,"Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n",ngsiz2,ngsiz1,dmb1->hlevel,dmb2->hlevel);
200 
201   PetscPrintf(PETSC_COMM_SELF, "[%d] Local matrix: %D X %D\n", rank, nlsiz2, nlsiz1);
202 
203   /* allocate the nnz, onz arrays based on block size and local nodes */
204   ierr = PetscCalloc2(nlghsiz2,&nnz,nlghsiz2,&onz);CHKERRQ(ierr);
205 
206   PetscPrintf(PETSC_COMM_SELF, "[%d] Coarse Local elements: %D, vertices: %D\n", rank, dmb1->elocal->size(), dmb1->vlocal->size());
207   PetscPrintf(PETSC_COMM_SELF, "[%d] Fine   Local elements: %D, vertices: %D\n", rank, dmb2->elocal->size(), dmb2->vlocal->size());
208   PetscPrintf(PETSC_COMM_SELF, "[%d] Sequence Start: %D, End: %D\n", rank, dmb2->seqstart, dmb2->seqend);
209 
210   /* Loop through the local elements and compute the relation between the current parent and the refined_level. */
211   for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) {
212 
213     const moab::EntityHandle ehandle = *iter;
214     std::vector<moab::EntityHandle> children;
215     std::vector<moab::EntityHandle> connp, connc;
216 
217     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
218     merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr);
219 
220     /* Get connectivity and coordinates of the parent vertices */
221     merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr);
222     for (unsigned ic=0; ic < children.size(); ic++) {
223       std::vector<moab::EntityHandle> tconnc;
224       /* Get handles of the parent vertices in canonical order and intersect */
225       merr = dmb2->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr);
226       for (unsigned tc=0; tc<tconnc.size(); tc++) {
227         if (std::find(connc.begin(), connc.end(), tconnc[tc]) == connc.end())
228           connc.push_back(tconnc[tc]);
229       }
230     }
231     //PetscPrintf(PETSC_COMM_SELF, "[%d] EntityHandle %d, children = %d, total intersection = %d\n", rank, ehandle, children.size(), connc.size());
232 
233     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
234     /* TODO: specific to scalar system - use GetDofs */
235     //ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr);
236     ierr = DMMoabGetFieldDofsLocal(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr);
237 
238     //PetscPrintf(PETSC_COMM_SELF, "[%d] EntityHandle %d, dofs [4] = %d, %d, %d, %d\n", rank, ehandle, dofsc[0], dofsc[1], dofsc[2], dofsc[3]);
239     // if (rank == 1) {
240     //   for (unsigned tp=0;tp<connc.size(); tp++)
241     //     PetscPrintf(PETSC_COMM_SELF, "[%d] EntityHandle %d \t %d  -- dofs [%d] = %d, %d\n", rank, ehandle, connc[tp]-dmb2->seqstart, tp, dmb2->gidmap[(PetscInt)connc[tp]-dmb2->seqstart], dofsc[tp]);
242     // }
243     for (unsigned tp=0;tp<connp.size(); tp++) {
244       // ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
245       if (dmb1->vowned->find(connp[tp]) != dmb1->vowned->end()) {
246         for (unsigned tc=0;tc<connc.size(); tc++) {
247           nnz[dofsc[tc]]++;
248           //PetscPrintf(PETSC_COMM_SELF, "[%d] Found nnz coupling for %D = %d\n", rank, connp[tp], dofsc[tc]);
249         }
250       }
251       else if (dmb1->vghost->find(connp[tp]) != dmb1->vghost->end()) {
252         for (unsigned tc=0;tc<connc.size(); tc++) {
253           onz[dofsc[tc]]++;
254           //PetscPrintf(PETSC_COMM_SELF, "[%d] Found onz coupling for %D = %d\n", rank, connp[tp], dofsc[tc]);
255         }
256       }
257       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in parent level %D\n", connc[tp]);
258     }
259     for(unsigned tc = 0; tc < connc.size(); tc++) {
260       if (dmb2->vowned->find(connc[tc]) != dmb2->vowned->end()) nnz[dofsc[tc]]++;
261       else if (dmb2->vghost->find(connc[tc]) != dmb2->vghost->end()) onz[dofsc[tc]]++;
262       else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in child level %D\n", connc[tc]);
263     }
264   }
265 
266 /*
267   int i=0;
268   std::vector<moab::EntityHandle> adjs;
269   for(moab::Range::iterator iter = dmb1->vowned->begin(); iter!= dmb1->vowned->end(); iter++, i++) {
270     merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr);
271     nnz[i] -= adjs.size();
272     adjs.clear();
273   }
274   i=0;
275   for(moab::Range::iterator iter = dmb1->vghost->begin(); iter!= dmb1->vghost->end(); iter++, i++) {
276     //merr = dmb1->hierarchy->get_adjacencies(&(*iter), 1, 0, false, adjs, moab::Interface::UNION);MBERRNM(merr);
277     merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr);
278     onz[i] -= adjs.size();
279     adjs.clear();
280   }
281 */
282 
283   PetscInt* ldofs = dmb2->lidmap;
284   PetscInt* gdofs = dmb2->gidmap;
285   ionz=onz[0];
286   innz=nnz[0];
287   for (int tc=0; tc < nlsiz2; tc++) {
288     // check for maximum allowed sparsity = fully dense
289     nnz[tc] = std::min(nlsiz1,nnz[tc]);
290     onz[tc] = std::min(nlsiz1,onz[tc]);
291 
292     innz = (innz < nnz[tc] ? nnz[tc] : innz);
293     ionz = (ionz < onz[tc] ? onz[tc] : ionz);
294     PetscPrintf(PETSC_COMM_SELF, "[%D]: %D NNZ = %D, ONZ = %D\n", rank, gdofs[ldofs[tc]], nnz[tc], onz[tc]);
295   }
296   PetscPrintf(PETSC_COMM_SELF, "[%D]: Final: INNZ = %D, IONZ = %D\n", rank, innz, ionz);
297 
298   MPI_Barrier(PETSC_COMM_WORLD);
299 
300   /* create interpolation matrix */
301   PetscPrintf(PETSC_COMM_SELF, "[%D]: Creating matrix\n", rank);
302   ierr = MatCreate(PetscObjectComm((PetscObject)dm2), interpl);CHKERRQ(ierr);
303   PetscPrintf(PETSC_COMM_SELF, "[%D]: set sizes\n", rank);
304   ierr = MatSetSizes(*interpl, nlsiz2, ngsiz1, ngsiz2, ngsiz1);CHKERRQ(ierr);
305   //ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr);
306   PetscPrintf(PETSC_COMM_SELF, "[%D]: set type\n", rank);
307   ierr = MatSetType(*interpl,MATAIJ);CHKERRQ(ierr);
308   PetscPrintf(PETSC_COMM_SELF, "[%D]: set from opts\n", rank);
309   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
310 
311   PetscPrintf(PETSC_COMM_SELF, "[%D]: Setting prealloc\n", rank);
312   ierr = MatSeqAIJSetPreallocation(*interpl,innz,nnz);CHKERRQ(ierr);
313   ierr = MatMPIAIJSetPreallocation(*interpl,innz,nnz,ionz,onz);CHKERRQ(ierr);
314   //ierr = MatMPIAIJSetPreallocation(*interpl,innz,0,ionz,0);CHKERRQ(ierr);
315 
316   PetscPrintf(PETSC_COMM_SELF, "[%D]: Cleaning arrays\n", rank);
317   /* clean up temporary memory */
318   ierr = PetscFree2(nnz,onz);CHKERRQ(ierr);
319 
320   /* set up internal matrix data-structures */
321   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
322   //ierr = MatZeroEntries(*interpl);CHKERRQ(ierr);
323 
324   ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr);
325 
326   factor = std::pow(2.0 /*degree_P_for_refinement*/,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0);
327 
328   ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
329 
330   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
331   for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) {
332 
333     const moab::EntityHandle ehandle = *iter;
334     std::vector<moab::EntityHandle> children;
335     std::vector<moab::EntityHandle> connp, connc;
336 
337     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
338     merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr);
339 
340     /* Get connectivity and coordinates of the parent vertices */
341     merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr);
342     for (unsigned ic=0; ic < children.size(); ic++) {
343       std::vector<moab::EntityHandle> tconnc;
344       /* Get coordinates of the parent vertices in canonical order */
345       merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr);
346       for (unsigned tc=0; tc<tconnc.size(); tc++) {
347         connc.push_back(tconnc[tc]);
348       }
349     }
350 
351     std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size());
352     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
353     merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr);
354     merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr);
355 
356     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
357     /* TODO: specific to scalar system - use GetDofs */
358     ierr = DMMoabGetFieldDofsLocal(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr);
359     ierr = DMMoabGetFieldDofsLocal(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr);
360 
361     /* Compute the interpolation weights by determining distance of 1-ring
362        neighbor vertices from current vertex */
363     for (unsigned i=0;i<connp.size(); i++) {
364       double normsum=0.0;
365       for (unsigned j=0;j<connc.size(); j++) {
366         values_phi[j] = 0.0;
367         for (unsigned k=0;k<3; k++)
368           values_phi[j] += std::pow(pcoords[i*3+k]-ccoords[k+j*3], dim);
369         if (values_phi[j] < 1e-12) {
370           values_phi[j] = 1e12;
371         }
372         else {
373           //values_phi[j] = std::pow(values_phi[j], -1.0/dim);
374           values_phi[j] = std::pow(values_phi[j], -1.0);
375           normsum += values_phi[j];
376         }
377       }
378       for (unsigned j=0;j<connc.size(); j++) {
379         if (values_phi[j] > 1e11)
380           values_phi[j] = factor*0.5/connc.size();
381         else
382           values_phi[j] = factor*values_phi[j]*0.5/(connc.size()*normsum);
383       }
384       ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
385     }
386 
387     /* check if element is on the boundary */
388     //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr);
389     dbdry.resize(connc.size());
390     ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry.data());CHKERRQ(ierr);
391     eonbnd=PETSC_FALSE;
392     for (unsigned i=0; i< connc.size(); ++i)
393       if (dbdry[i]) eonbnd=PETSC_TRUE;
394 
395     values_phi.clear();
396     values_phi.resize(connp.size());
397     /* apply dirichlet boundary conditions */
398     if (eonbnd) {
399 
400       ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
401       ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
402       /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */
403       //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr);
404       for (unsigned i=0; i < connc.size(); i++) {
405         if (dbdry[i]) {  /* dirichlet node */
406           /* think about strongly imposing dirichlet */
407           //bndrows.push_back(dofsc[i]);
408 
409           ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
410           //values_phi[0]=1.0;
411           //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
412         }
413       }
414 
415       ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
416       ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
417     }
418 
419     //get interpolation weights
420     //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr);
421     // for (int j=0;j<dofs_per_element; j++)
422     //  std::cout<<"values "<<values_phi[j]<<std::endl;
423 
424     //get row and column indices, zero weights are ignored
425     /*
426     int nz_ind = 0;
427     idx = dmb2->vowned->index(vhandle);
428     for (int j=0;j<dofs_per_element; j++){
429       idy[nz_ind] = dmb1->vowned->index(connectivity[j]);
430       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]);
431       //values_phi[nz_ind] = values_phi[j];
432       nz_ind = nz_ind+1;
433     }
434     */
435 
436     //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr);
437     //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
438   }
439 
440   //PetscPrintf(PETSC_COMM_WORLD, "[Boundary vertices = %D] :: A few: %D %D %D %D \n", bndrows.size(), bndrows[0], bndrows[1], bndrows[2], bndrows[3]);
441   //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr);
442 
443   ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
444   ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
445   PetscFunctionReturn(0);
446 }
447 
448 #undef __FUNCT__
449 #define __FUNCT__ "DMCreateInjection_Moab"
450 /*@
451   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
452   by succesively refining a coarse mesh, already defined in the DM object
453   provided by the user.
454 
455   Collective on MPI_Comm
456 
457   Input Parameter:
458 . dmb  - The DMMoab object
459 
460   Output Parameter:
461 . nlevels   - The number of levels of refinement needed to generate the hierarchy
462 + ldegrees  - The degree of refinement at each level in the hierarchy
463 
464   Level: beginner
465 
466 .keywords: DMMoab, create, refinement
467 @*/
468 PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx)
469 {
470   //DM_Moab        *dmmoab;
471 
472   PetscFunctionBegin;
473   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
474   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
475   //dmmoab = (DM_Moab*)(dm1)->data;
476 
477   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
478   PetscFunctionReturn(0);
479 }
480 
481 
482 #undef __FUNCT__
483 #define __FUNCT__ "DM_UMR_Moab_Private"
484 PetscErrorCode  DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref)
485 {
486   PetscErrorCode  ierr;
487   PetscInt        i,dim;
488   DM              dm2;
489   moab::ErrorCode merr;
490   DM_Moab        *dmb = (DM_Moab*)dm->data,*dd2;
491 
492   PetscFunctionBegin;
493   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
494   PetscValidPointer(dmref,4);
495 
496   if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) {
497     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);
498     if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1);
499     *dmref = PETSC_NULL;
500     PetscFunctionReturn(0);
501   }
502 
503   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
504   dd2 = (DM_Moab*)dm2->data;
505 
506   dd2->mbiface = dmb->mbiface;
507   dd2->pcomm = dmb->pcomm;
508   dd2->icreatedinstance = PETSC_FALSE;
509   dd2->nghostrings=dmb->nghostrings;
510 
511   /* set the new level based on refinement/coarsening */
512   if (refine) {
513     dd2->hlevel=dmb->hlevel+1;
514   }
515   else {
516     dd2->hlevel=dmb->hlevel-1;
517   }
518 
519   /* Copy the multilevel hierarchy pointers in MOAB */
520   dd2->hierarchy = dmb->hierarchy;
521   dd2->nhlevels = dmb->nhlevels;
522   ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr);
523   for (i=0; i<=dd2->nhlevels; i++) {
524     dd2->hsets[i]=dmb->hsets[i];
525   }
526   dd2->fileset = dd2->hsets[dd2->hlevel];
527 
528   /* do the remaining initializations for DMMoab */
529   dd2->bs = dmb->bs;
530   dd2->numFields = dmb->numFields;
531   dd2->rw_dbglevel = dmb->rw_dbglevel;
532   dd2->partition_by_rank = dmb->partition_by_rank;
533   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
534   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
535   dd2->read_mode = dmb->read_mode;
536   dd2->write_mode = dmb->write_mode;
537 
538   /* set global ID tag handle */
539   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
540 
541   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
542 
543   ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr);
544   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
545   ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr);
546 
547   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
548   dm2->ops->creatematrix = dm->ops->creatematrix;
549 
550   /* copy fill information if given */
551   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
552 
553   /* copy vector type information */
554   ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr);
555   ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr);
556   dd2->numFields = dmb->numFields;
557   if (dmb->numFields) {
558     ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr);
559   }
560 
561   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
562 
563   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
564   ierr = DMSetUp(dm2);CHKERRQ(ierr);
565 
566   *dmref = dm2;
567   PetscFunctionReturn(0);
568 }
569 
570 
571 #undef __FUNCT__
572 #define __FUNCT__ "DMRefine_Moab"
573 /*@
574   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
575   by succesively refining a coarse mesh, already defined in the DM object
576   provided by the user.
577 
578   Collective on DM
579 
580   Input Parameter:
581 + dm  - The DMMoab object
582 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
583 
584   Output Parameter:
585 . dmf - the refined DM, or NULL
586 
587   Note: If no refinement was done, the return value is NULL
588 
589   Level: developer
590 
591 .keywords: DMMoab, create, refinement
592 @*/
593 PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf)
594 {
595   PetscErrorCode  ierr;
596 
597   PetscFunctionBegin;
598   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
599 
600   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr);
601   PetscFunctionReturn(0);
602 }
603 
604 #undef __FUNCT__
605 #define __FUNCT__ "DMCoarsen_Moab"
606 /*@
607   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
608   by succesively refining a coarse mesh, already defined in the DM object
609   provided by the user.
610 
611   Collective on DM
612 
613   Input Parameter:
614 . dm  - The DMMoab object
615 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
616 
617   Output Parameter:
618 . dmf - the coarsened DM, or NULL
619 
620   Note: If no coarsening was done, the return value is NULL
621 
622   Level: developer
623 
624 .keywords: DMMoab, create, coarsening
625 @*/
626 PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc)
627 {
628   PetscErrorCode  ierr;
629 
630   PetscFunctionBegin;
631   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
632 
633   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr);
634   PetscFunctionReturn(0);
635 }
636