Timo,
thank you. These two files are extremely helpful and I guess I understood a
lot from them. Extremely well written, by the way.
I really appreciate it.

 I would like to associate to each gauss point on the interface, identified
by the command
      const auto &q_points = fiv.get_quadrature_points();
 the values of the unknown fields at the two faces of the two cells that
form the interface itself.
 It should not be too difficult I think.
 Perhaps I should use

      FEFaceValues<dim> fe_face0_values (fe, face_quadrature_formula,
                                        update_values         |
update_quadrature_points  |
                                        update_normal_vectors |
update_JxW_values);
      FEFaceValues<dim> fe_face1_values (fe, face_quadrature_formula,
                                        update_values         |
update_quadrature_points  |
                                        update_normal_vectors |
update_JxW_values);

      fe_face0_values.reinit (cell0, face_number_on_cell0);
      fe_face1_values.reinit (cell0, face_number_on_cell1);

      std::vector< Vector< double > > old_solution_values_at_face_0 (
face_quadrature_formula.size(), Vector< double >(dofs.GPDofs) );
      fe_face0_values.get_function_values ( accumulated_displacement,
old_solution_values_at_face_0 );
      std::vector< Vector< double > > old_solution_values_at_face_1 (
face_quadrature_formula.size(), Vector< double >(dofs.GPDofs) );
      fe_face1_values.get_function_values ( accumulated_displacement,
old_solution_values_at_face_0 );

I am concerned about the correspondence of quadrature points on the
interface and on the cell faces. Is there a simple and effective way to
find their association?

Thanks,
Alberto


*Alberto Salvadori* Dipartimento di Ingegneria Meccanica e Industriale
(DIMI)
 Universita` di Brescia, via Branze 43, 25123 Brescia
 Italy
 tel 030 3715426

e-mail:
 alberto.salvad...@unibs.it
web-page:
 http://m4lab.unibs.it/faculty.html



On Wed, Dec 16, 2020 at 11:25 PM Timo Heister <heis...@clemson.edu> wrote:

> Alberto,
>
> > I would like to have a few more details on the class
> FEInterfaceValues<dim>.
> > since it is my feeling that  FEInterfaceValues<dim>
> > might be conveniently used for interfaces, even out of DG methods.
>
> Yes, I would think so.
>
> > Specifically I'd like to figure out if this notion can be used without
> recurring to the MeshWorker::mesh_loop,
> > which I am not quite confident at present.
>
> Note that we are not using the MeshWorker machinery but a single
> function that contains the loops over the cells. The code of the
> function is "relatively" easy to understand and you will appreciate
> the logic for the following things you would have to deal with when
> looping manually:
> - assembling from the finer cell for interfaces with adaptive refinement
> - correctly handling each facet once (or twice if desired) even in parallel
> - handling of periodic boundaries
>
> > I'd like to play autonomously with the class FEInterfaceValues<dim>
> > but unfortunately I cannot grasp the notion of subface, which is
> apparently
> > required by the reinit as in step 12:
> >  fe_iv.reinit(cell, f, sf, ncell, nf, nsf);
>
> With adaptivity you will need to pass the subface as returned by
> neighbor_of_coarser_neighbor, otherwise it is
> numbers::invalid_unsigned_int. See
>
> https://github.com/dealii/dealii/blob/691ff4907964605dd21e94f06c880786af2cd612/tests/feinterface/fe_interface_values_03.cc#L94-L100
> and
>
> https://github.com/dealii/dealii/blob/691ff4907964605dd21e94f06c880786af2cd612/tests/feinterface/fe_interface_values_02.cc#L88-L93
>
> > Can anyone address me some reading on this notion and how it is dealt
> with
> > by the MeshWorker::mesh_loop ?
>
> Take a look at the tests in tests/feinterface/:
> https://github.com/dealii/dealii/tree/master/tests/feinterface
>
> Many of them use FEInterfaceValues directly on a small test mesh.
>
> --
> Timo Heister
> http://www.math.clemson.edu/~heister/
>
> --
> 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/CAMRj59Fq2wUezf%3DikmWF5_BPB-Q1zzYgjOgd4XESUYdwX_pJ_A%40mail.gmail.com
> .
>

-- 


Informativa sulla Privacy: http://www.unibs.it/node/8155 
<http://www.unibs.it/node/8155>

-- 
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/CABcATpe%2BVBb%3DPt9WxDjnDGRpkfM9p-WGsf%2BZn_KSRYeLkZYi6w%40mail.gmail.com.

Reply via email to