Jane,

> I am solving some equations on a moving domain which essentially compacts. I 
> am using something similar to that in step-18 to move my mesh around. So i 
> have something like
> pcout << " Moving mesh..." << std::endl;
> std::vector<bool> vertex_touched (triangulation.n_vertices 
> <https://www.dealii.org/8.5.0/doxygen/deal.II/classTriangulation.html#ad49cc6eadbda275e0c2ef0155469d656>(),
> false);
> for (typename DoFHandler<dim>::active_cell_iterator 
> <https://www.dealii.org/8.5.0/doxygen/deal.II/classDoFHandler.html>
> cell = dof_handler.begin_active 
> <https://www.dealii.org/8.5.0/doxygen/deal.II/classDoFHandler.html#ad2df51a32781906218910731e063ac62>
>  
> ();
> cell != dof_handler.end 
> <https://www.dealii.org/8.5.0/doxygen/deal.II/classDoFHandler.html#a05c70e1862a2ad145f91f9da1f44cc28>();
>  
> ++cell)
> for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
> if (vertex_touched[cell->vertex_index(v)] == false)
> {
> vertex_touched[cell->vertex_index(v)] = true;
> Point<dim> <https://www.dealii.org/8.5.0/doxygen/deal.II/classPoint.html> 
> vertex_displacement;
> for (unsigned int d=0; d<dim; ++d)
> vertex_displacement[d]
> = solution.block(0)(cell->vertex_dof_index(v,d));
> cell->vertex(v) += vertex_displacement*timestep_size;
> }
> 
> However, I realised that this moves each vertex by the velocity at that 
> point, 
> but I need, for example, a vertex on the top boundary to have moved as much 
> as 
> the ones below it has moved as well as the velocity at the points itself.

This is actually an interesting question. I had a visitor -- Sean McGovern in 
my lab earlier this year -- who was looking at exactly this kind of problem.

The issue is that the displacement you need to compute is given by

   d(x,y,z) = \int_z^0 v(x,y,zeta) dzeta

where v(.,.,.) is the velocity field and you need to integrate this on a line 
from the point where you need the displacement to a point on the surface (or, 
at the bottom, depending on your point of view). This is, in general, a 
quantity that is hard to compute in the finite element method, in particular 
if you happen to have an unstructured mesh.

But it also turns out to satisfy the differential equation

   d/dz d(x,y,z) = v(x,y,z)

(possibly with a negative sign) or, if you're into partial differential 
equations,

   e_z . grad d = v

where e_z = (0,0,1). This now is something we know how to solve with the 
finite element method: It's just a regular advection equation.

In other words, one way of solving what you want to solve for is to use the 
power of reformulating the problem in a way that makes it more convenient to 
solve the problem -- namely, by adding one more variable to your system of PDEs.

Best
  W.


-- 
------------------------------------------------------------------------
Wolfgang Bangerth          email:                 [email protected]
                            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 [email protected].
For more options, visit https://groups.google.com/d/optout.

Reply via email to