xref: /petsc/src/dm/impls/moab/dmmbmg.cxx (revision e882eb38466228dc46cfd515783348bd1ccf5a7f)
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->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   PetscBool       eonbnd,dbdry[27];
159   std::vector<int> bndrows;
160 
161   PetscFunctionBegin;
162   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
163   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
164   dmb1 = (DM_Moab*)(dm1)->data;
165   dmb2 = (DM_Moab*)(dm2)->data;
166 
167   ierr = MatCreate(PetscObjectComm((PetscObject)dm1), interpl);CHKERRQ(ierr);
168   ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr);
169   ierr = MatSetSizes(*interpl, dmb2->nloc, dmb1->nloc, dmb2->n, dmb1->n);CHKERRQ(ierr);
170 
171   /* TODO: This is a hack for the rectangular system - decipher NNZ pattern better */
172   ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr);
173 
174   ierr = MatSeqAIJSetPreallocation(*interpl,dmb2->nloc,0);CHKERRQ(ierr);
175   ierr = MatMPIAIJSetPreallocation(*interpl,dmb2->nloc,0,dmb2->nghost,0);CHKERRQ(ierr);
176 
177   /* set up internal matrix data-structures */
178   ierr = MatSetUp(*interpl);CHKERRQ(ierr);
179 
180   ierr = MatZeroEntries(*interpl);CHKERRQ(ierr);
181 
182   ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr);
183 
184   PetscPrintf(PETSC_COMM_WORLD, "Creating a %D DIM matrix that is %D X %D and setting diagonal to 1.0\n", dim, dmb1->nloc, dmb2->nloc);
185 
186   double factor = std::pow(2.0,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0);
187 
188   /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */
189   for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) {
190 
191     const moab::EntityHandle ehandle = *iter;
192     std::vector<moab::EntityHandle> children;
193     std::vector<moab::EntityHandle> connp, connc;
194 
195     /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */
196     merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr);
197 
198     /* Get connectivity and coordinates of the parent vertices */
199     merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr);
200     for (unsigned ic=0; ic < children.size(); ic++) {
201       std::vector<moab::EntityHandle> tconnc;
202       /* Get coordinates of the parent vertices in canonical order */
203       merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr);
204       for (unsigned tc=0; tc<tconnc.size(); tc++) {
205         connc.push_back(tconnc[tc]);
206       }
207     }
208 
209     std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size());
210     /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */
211     merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr);
212     merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr);
213 
214     std::vector<int> dofsp(connp.size()), dofsc(connc.size());
215     /* TODO: specific to scalar system - use GetDofs */
216     ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr);
217     ierr = DMMoabGetFieldDofs(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr);
218 
219     /* Compute the interpolation weights by determining distance of 1-ring
220        neighbor vertices from current vertex */
221     for (unsigned i=0;i<connp.size(); i++) {
222       double normsum=0.0;
223       for (unsigned j=0;j<connc.size(); j++) {
224         unsigned offset=j;
225         values_phi[offset] = 0.0;
226         for (unsigned k=0;k<3; k++)
227           values_phi[offset] += (pcoords[i*3+k]-ccoords[k+j*3])*(pcoords[i*3+k]-ccoords[k+j*3]);
228         if (values_phi[offset] < 1e-12) {
229           values_phi[offset] = 1e12;
230         }
231         else {
232           values_phi[offset] = 1.0/(values_phi[offset]);
233           normsum += values_phi[offset];
234         }
235       }
236       for (unsigned j=0;j<connc.size(); j++) {
237         unsigned offset=j;
238         if (values_phi[offset] > 1e11)
239           values_phi[offset] = factor*0.5/connc.size();
240         else
241           values_phi[offset] = factor*values_phi[offset]*0.5/(connc.size()*normsum);
242       }
243       ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr);
244     }
245 
246     /* check if element is on the boundary */
247     //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr);
248     ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr);
249     eonbnd=PETSC_FALSE;
250     for (unsigned i=0; i< connc.size(); ++i)
251       if (dbdry[i]) eonbnd=PETSC_TRUE;
252 
253     values_phi.clear();
254     values_phi.resize(connp.size());
255     /* apply dirichlet boundary conditions */
256     if (eonbnd) {
257 
258       ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
259       ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
260       /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */
261       //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr);
262       for (unsigned i=0; i < connc.size(); i++) {
263         if (dmb2->hierarchy->is_boundary_vertex(connc[i])) {  /* dirichlet node */
264           /* think about strongly imposing dirichlet */
265           //bndrows.push_back(dofsc[i]);
266 
267           ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
268           //values_phi[0]=1.0;
269           //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
270         }
271       }
272 
273       ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
274       ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
275     }
276 
277     //get interpolation weights
278     //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr);
279     // for (int j=0;j<dofs_per_element; j++)
280     //  std::cout<<"values "<<values_phi[j]<<std::endl;
281 
282     //get row and column indices, zero weights are ignored
283     /*
284     int nz_ind = 0;
285     idx = dmb2->vowned->index(vhandle);
286     for (int j=0;j<dofs_per_element; j++){
287       idy[nz_ind] = dmb1->vowned->index(connectivity[j]);
288       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]);
289       //values_phi[nz_ind] = values_phi[j];
290       nz_ind = nz_ind+1;
291     }
292     */
293 
294     //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr);
295     //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr);
296   }
297 
298   //PetscPrintf(PETSC_COMM_WORLD, "[Boundary vertices = %D] :: A few: %D %D %D %D \n", bndrows.size(), bndrows[0], bndrows[1], bndrows[2], bndrows[3]);
299   //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr);
300 
301   ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
302   ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
303   PetscFunctionReturn(0);
304 }
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "DMCreateInjection_Moab"
308 /*@
309   DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy
310   by succesively refining a coarse mesh, already defined in the DM object
311   provided by the user.
312 
313   Collective on MPI_Comm
314 
315   Input Parameter:
316 . dmb  - The DMMoab object
317 
318   Output Parameter:
319 . nlevels   - The number of levels of refinement needed to generate the hierarchy
320 + ldegrees  - The degree of refinement at each level in the hierarchy
321 
322   Level: beginner
323 
324 .keywords: DMMoab, create, refinement
325 @*/
326 PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx)
327 {
328   //DM_Moab        *dmmoab;
329 
330   PetscFunctionBegin;
331   PetscValidHeaderSpecific(dm1,DM_CLASSID,1);
332   PetscValidHeaderSpecific(dm2,DM_CLASSID,2);
333   //dmmoab = (DM_Moab*)(dm1)->data;
334 
335   PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n");
336   PetscFunctionReturn(0);
337 }
338 
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "DM_UMR_Moab_Private"
342 PetscErrorCode  DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref)
343 {
344   PetscErrorCode  ierr;
345   PetscInt        i,dim;
346   DM              dm2;
347   moab::ErrorCode merr;
348   DM_Moab        *dmb = (DM_Moab*)dm->data,*dd2;
349 
350   PetscFunctionBegin;
351   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
352   PetscValidPointer(dmref,4);
353 
354   if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) {
355     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);
356     if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1);
357     *dmref = PETSC_NULL;
358     PetscFunctionReturn(0);
359   }
360 
361   ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr);
362   dd2 = (DM_Moab*)dm2->data;
363 
364   dd2->mbiface = dmb->mbiface;
365   dd2->pcomm = dmb->pcomm;
366   dd2->icreatedinstance = PETSC_FALSE;
367 
368   /* set the new level based on refinement/coarsening */
369   if (refine) {
370     dd2->hlevel=dmb->hlevel+1;
371   }
372   else {
373     dd2->hlevel=dmb->hlevel-1;
374   }
375 
376   /* Copy the multilevel hierarchy pointers in MOAB */
377   dd2->hierarchy = dmb->hierarchy;
378   dd2->nhlevels = dmb->nhlevels;
379   ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr);
380   for (i=0; i<=dd2->nhlevels; i++) {
381     dd2->hsets[i]=dmb->hsets[i];
382   }
383   dd2->fileset = dd2->hsets[dd2->hlevel];
384 
385   /* do the remaining initializations for DMMoab */
386   dd2->bs = dmb->bs;
387   dd2->numFields = dmb->numFields;
388   dd2->rw_dbglevel = dmb->rw_dbglevel;
389   dd2->partition_by_rank = dmb->partition_by_rank;
390   ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr);
391   ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr);
392   dd2->read_mode = dmb->read_mode;
393   dd2->write_mode = dmb->write_mode;
394 
395   /* set global ID tag handle */
396   ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr);
397 
398   merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr);
399 
400   ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr);
401   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
402   ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr);
403 
404   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
405   dm2->ops->creatematrix = dm->ops->creatematrix;
406 
407   /* copy fill information if given */
408   ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr);
409 
410   /* copy vector type information */
411   ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr);
412   ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr);
413   ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr);
414 
415   ierr = DMSetFromOptions(dm2);CHKERRQ(ierr);
416 
417   /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */
418   ierr = DMSetUp(dm2);CHKERRQ(ierr);
419 
420   *dmref = dm2;
421   PetscFunctionReturn(0);
422 }
423 
424 
425 #undef __FUNCT__
426 #define __FUNCT__ "DMRefine_Moab"
427 /*@
428   DMRefine_Moab - Generate a multi-level uniform refinement hierarchy
429   by succesively refining a coarse mesh, already defined in the DM object
430   provided by the user.
431 
432   Collective on DM
433 
434   Input Parameter:
435 + dm  - The DMMoab object
436 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
437 
438   Output Parameter:
439 . dmf - the refined DM, or NULL
440 
441   Note: If no refinement was done, the return value is NULL
442 
443   Level: developer
444 
445 .keywords: DMMoab, create, refinement
446 @*/
447 PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf)
448 {
449   PetscErrorCode  ierr;
450 
451   PetscFunctionBegin;
452   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
453 
454   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr);
455   PetscFunctionReturn(0);
456 }
457 
458 #undef __FUNCT__
459 #define __FUNCT__ "DMCoarsen_Moab"
460 /*@
461   DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy
462   by succesively refining a coarse mesh, already defined in the DM object
463   provided by the user.
464 
465   Collective on DM
466 
467   Input Parameter:
468 . dm  - The DMMoab object
469 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL)
470 
471   Output Parameter:
472 . dmf - the coarsened DM, or NULL
473 
474   Note: If no coarsening was done, the return value is NULL
475 
476   Level: developer
477 
478 .keywords: DMMoab, create, coarsening
479 @*/
480 PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc)
481 {
482   PetscErrorCode  ierr;
483 
484   PetscFunctionBegin;
485   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
486 
487   ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr);
488   PetscFunctionReturn(0);
489 }
490