Hi Alberto,

I have dealt with a similar problem where I needed to transmit solution values 
between two different problems that shared a common interface. If I remember 
correctly, the way that I did this was to precompute the positions at which I 
would need the solutions from problem 1 for problem 2. In a post processing 
step for problem 1 I then computed this data up front (cache it) and then fetch 
it from problem 2. This avoided the issue of going to look for which cell in 
problem 1 a point in problem 2 lay. I did this all manually, but I guess that 
this approach could be made simpler now that we have 
GridTools::compute_point_locations() 
<https://www.dealii.org/developer/doxygen/deal.II/namespaceGridTools.html#a8e8bb9211264d2106758ac4d7184117e>
 which allows a speedy lookup between a point and its containing cell. 

FEFieldFunction has a cache so you could also investigate using it (maybe it 
wouldn’t be too slow for your case). It also has the set_active_cell() 
<https://www.dealii.org/9.0.0/doxygen/deal.II/classFunctions_1_1FEFieldFunction.html#a0206a45c90d523792eea8bd725d14788>
 function so you could create a lookup between a point in problem 1 and a cell 
in problem 2 again using GridTools::Cache 
<https://www.dealii.org/developer/doxygen/deal.II/classGridTools_1_1Cache.html> 
(although I’m guessing that this is what it does internally).

I hope that this gives you some ideas to play with!

Best,
J-P

> On 07 Aug 2018, at 17:47, Alberto Salvadori <alberto.salvad...@unibs.it> 
> wrote:
> 
> Dear community,
> 
> I am asking some advices on the following issue. I am solving a simple 
> problem, say a Poisson problem in the unknown field u, and a more involved 
> problem separately. This second problem requires the values of u in the 
> Neumann boundary conditions.
> 
> Accordingly, I guess one could solve the Laplacian first and calculate the 
> numerical solution for u, say un. Afterwards one builds a solver for the more 
> complex operator and in the Neumann part of the code - that may look like 
> this for parallel::shared triangulations:
> 
>   for (unsigned int face_number=0; 
> face_number<GeometryInfo<dim>::faces_per_cell; ++face_number)
> 
>         if (
> 
>             cell->face(face_number)->at_boundary()
> 
>             &&
> 
>             cell->face(face_number)->boundary_id() == 2    // Neumann 
> boundaries
> 
>             )
> 
>         {
> 
>           
>           fe_face_values.reinit (cell, face_number);
> 
>           
>           // define points and normals
> 
>           
>           std::vector< Point<dim> >    points  = 
> fe_face_values.get_quadrature_points();
> 
>           std::vector< Tensor<1,dim> > normals = 
> fe_face_values.get_all_normal_vectors();
> 
>           
>           // calculate neumann values
> 
>           
>           for (unsigned int q_point=0; q_point<n_face_q_points; ++q_point)
> 
>           {
> 
>             
>             // values: mechanical
> 
>             Tensor<1,dim> mech_neumann_value;
> 
>             neumann_bc_for_mech.bc_value( points[q_point], normals[q_point], 
> mech_neumann_value );
> 
> 
> 
>              .....
> 
>      
> 
> one needs the value of the field un at points  points[q_point] to be passed 
> to neumann_bc_for_mech.bc_value.
> Which is an effective way to calculate this amount ?
> 
> Thank you.
> 
> Alberto
> 
> 
> 
> 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/ 
> <http://www.dealii.org/>
> For mailing list/forum options, see 
> https://groups.google.com/d/forum/dealii?hl=en 
> <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 
> <mailto:dealii+unsubscr...@googlegroups.com>.
> For more options, visit https://groups.google.com/d/optout 
> <https://groups.google.com/d/optout>.

-- 
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.
For more options, visit https://groups.google.com/d/optout.

Reply via email to