Alberto, 

be aware of the fact that, when you extract the boundary mesh, the orientation 
of the cells may be different w.r.t. to the corresponding face on the 
volumetric mesh. In other words, using 

FEValues<dim, spacedim>  codim_one(…);
FEFaceValues<dim>        codim_zero(…);

codim_one.reinit(codim_one_cell);
codim_zero.reinit(codmi_zero_cell, corresponding_face);

in general, you will *NOT* get the same ordering of the quadrature points, and 
of the dof indices.

In the past, we proceeded by constructing a dof_map between codim zero and 
codim one (instead of trying to assemble together fe values coming from grids 
with different dimensions), or by constructing (by hand) the permutation of the 
dof indices.

As far as your question is concerned:

  c1_to_c0 = GridGenerator::extract_boundary_mesh(triangulation, c1_tria);

  // map volume mesh face -> codimension 1 mesh cell
  for (auto c1_cell : c1_tria.active_cell_iterators())
    c0_to_c1[c1_to_c0[c1_cell]] = c1_cell;

  // generate a mapping that maps codimension-1 cells
  // to codimension-0 cells and faces
  for (auto cell :
       triangulation
         .active_cell_iterators()) // disp_dof.active_cell_iterators())
    for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
      if (cell->face(f)->at_boundary())
        {
          auto c1_cell                      = c0_to_c1[cell->face(f)];
          c1_to_c0_cells_and_faces[c1_cell] = {cell, f};
        }


L.


> On 10 Sep 2019, at 12:15, Alberto Salvadori <alberto.salvad...@unibs.it> 
> wrote:
> 
> Thank you very much, Wolfgang, crystal clear as always.
> 
> It is still unclear to me, though, how to reinit() the FEFaceValues 
> for the cell making use of the 
> map returned by the extract_boundary_mesh() function. I learned to reinit 
> FEFaceValues as:
> 
> fe_face_values.reinit (cell, face_number);
> 
> Is there a "simple way" to get the (volume) cell iterator and the face_number 
> from the 
> parallel::shared::Triangulation<dim>::face_iterator in the map returned by 
> the extract_boundary_mesh() function? 
> 
> Or shall I use a different way to reinit FEFaceValues?
> 
> Thank you
> 
> Alberto
> 

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/1DD9B280-CD2F-4C83-81FE-3FC00A522B91%40gmail.com.

Reply via email to