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