Dear all,
I'm modifying MappingQ<dim, spacedim> to make it work for codimension one
meshes, and I'm trying to apply it to step-34 example.
Skimming through the implementation of MappingQ, I can see that the
construction of the laplace support points for dim = 2 is limited to the lines
of the quadrilaterals, even in the case where we want to deform also the inner
cells and not only the boundary cells.
A consequence of this is that, by instantiating MappingQ<2,3> (and changing a
few dim into spacedim), everything works perfectly fine on the edges of quads,
but no support points are computed in the interior of the cell. While this
might be enough in some cases, when we are talking about codimension one
meshes, the central support point of a quadrilateral element deformed with a Q2
mapping remains "flat" (in the sense that it is the average of the 4 *vertices*
of the cell).
The first question I have is: was this behavior intended? If so, why? The
current implementation is limited to an incomplete tensor product space, where
all the "bubble elements" are ignored.
Second question: people using MappingQEulerian should have problems with this
implementation, since the finite element space describing the displacement of
the support points is complete, but the mapping is not, leading to a
quasi-isoparametric implementation, where the interior of the cell is not
deformed according to the eulerian mapping.
I am not overly familiar with the internals of MappingQ, and I don't want to
start breaking things apart just to discover that what I want to accomplish is
more complicated than I thought it was, and I end up throwing weeks away...
This is what I'd do to fix this. Can anyone who has a better knowledge of the
MappingQ class point any flaw in my thoughts?
I'd basically just add the part of the code that computes the support points on
the interior of the cell, by adapting the one that is used on the boundary of
hexes, to the function
compute_support_points_laplace(const typename
Triangulation<dim,spacedim>::cell_iterator &cell,
std::vector<Point<spacedim> > &a) const
that is:
add_line_support_points (cell, a);
// -> add this: add_quad_support_points (cell, a);
apply_laplace_vector (laplace_on_quad_vector,a);
break;
which would require some resizing and rechecking here and there.
Is this all there is to it, or do I need to do something else?
Thanks,
Luca.
--
Luca Heltai <[email protected]>
http://people.sissa.it/~heltai/
Scuola Internazionale Superiore di Studi Avanzati
Phone: +39 040 3787 449, Office: 732
--
There are no answers, only cross references.
_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii