On 08/08/2018 02:03 AM, Jean-Paul Pelteret wrote:
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).
This is indeed what you would do if you had separate meshes. The better
approach, of course, is to use the same mesh for both problems. In that
case, you can use this to obtain the gradients of the Poisson solution
un at the quadrature points on the faces where you need these values for
the second problem:
FEFaceValues poisson_fe_face_values (...);
std::vector<Tensor<1,dim>> un_gradients (...);
std::vector<double> un_dot_n_values (...);
for (cell=...)
for (face=...)
if (face is interesting)
{
poisson_fe_face_values.reinit (cell, face);
poisson_fe_face_values.get_function_gradients (poisson_solution,
un_gradients);
for (q=0...)
un_dot_n_values[q] = un_gradients[q] *
poisson_fe_face_values.normal_vectors(q);
...do assembly for the second problem using the values of
(grad un).n just computed...
}
Best
W.
--
------------------------------------------------------------------------
Wolfgang Bangerth email: bange...@colostate.edu
www: http://www.math.colostate.edu/~bangerth/
--
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.