From: Luca Heltai Date: Sat, 25 Jul 2015 16:55:21 +0000 (+0200) Subject: Fixed documentation and a few const. X-Git-Tag: v8.3.0-rc2~1^2~6 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=9953f83c7d3165a0b51d0664d3a479a268fb2a6d;p=dealii.git Fixed documentation and a few const. Conflicts: doc/news/changes.h --- diff --git a/doc/news/8.2.1-vs-8.3.0.h b/doc/news/8.2.1-vs-8.3.0.h index 5334f601f3..bc8898b028 100644 --- a/doc/news/8.2.1-vs-8.3.0.h +++ b/doc/news/8.2.1-vs-8.3.0.h @@ -548,6 +548,14 @@ inconvenience this causes.
    +
  1. New: VectorTools::get_position_vector now works with arbitrary + FESystems, provided that the geometrical components are primitive, + and that you provide a component mask to select what components of + the finite element to use for the geometrical interpolation. +
    + (Luca Heltai, 2015/07/25) +
  2. +
  3. Fixed: parallel::distributed::refine_and_coarsen_fixed_fraction() in rare circumstances decided to not refine any cells at all, even if the refinement threshold was nonzero. This is now fixed. diff --git a/include/deal.II/numerics/vector_tools.h b/include/deal.II/numerics/vector_tools.h index 142950baae..7d0608694a 100644 --- a/include/deal.II/numerics/vector_tools.h +++ b/include/deal.II/numerics/vector_tools.h @@ -2345,17 +2345,15 @@ namespace VectorTools */ //@{ /** - * Given a DoFHandler containing at least a spacedim vector field, this - * function interpolates the Triangulation at the support points of a FE_Q() - * finite element of the same degree as dh->get_fe().degree. + * Given a DoFHandler containing at least a spacedim vector field, + * this function interpolates the Triangulation at the support + * points of a FE_Q() finite element of the same degree as the + * degree of the required components. * * Curved manifold are respected, and the resulting VECTOR will be - * geometrically consistent. - * - * The resulting map is guaranteed to be interpolatory at the support points - * of a FE_Q() finite element of the same degree as dh->get_fe().degree. - * Notice that this may or may not be meaningful, depending on the - * FiniteElement you have distributed in dh. + * geometrically consistent. The resulting map is guaranteed to be + * interpolatory at the support points of a FE_Q() finite element of + * the same degree as the degree of the required components. * * If the underlying finite element is an FE_Q(1)^spacedim, then the * resulting VECTOR is a finite element field representation of the vertices @@ -2367,6 +2365,9 @@ namespace VectorTools * the first spacedim components of the FiniteElement are assumed to * represent the geometry of the problem. * + * This function is only implemented for FiniteElements where the + * specified components are primitive. + * * @author Luca Heltai, 2015 */ template diff --git a/include/deal.II/numerics/vector_tools.templates.h b/include/deal.II/numerics/vector_tools.templates.h index 714cf9a6e7..452424eff1 100644 --- a/include/deal.II/numerics/vector_tools.templates.h +++ b/include/deal.II/numerics/vector_tools.templates.h @@ -7183,13 +7183,19 @@ namespace VectorTools for (unsigned int i=0; icomponent_to_base_index(i).first; + const unsigned int base_i = fe_system->component_to_base_index(i).first; Assert(degree == numbers::invalid_unsigned_int || degree == fe_system->base_element(base_i).degree, ExcNotImplemented()); + Assert(fe_system->base_element(base_i).is_primitive(), + ExcNotImplemented()); degree = fe_system->base_element(base_i).degree; } + // We create an intermediate FE_Q vector space, and then + // interpolate from that vector space to this one, by + // carefully selecting the right components. + FESystem feq(FE_Q(degree), spacedim); DH dhq(dh.get_tria()); dhq.distribute_dofs(feq); @@ -7211,11 +7217,17 @@ namespace VectorTools // v_j phi_j(q_i) = w_k psi_k(q_i) = w_k delta_ki = w_i // // where psi_k are the basis functions for fe_q, and phi_i are - // the basis functions of the target space. + // the basis functions of the target space while q_i are the + // support points for the fe_q space. With this choice of + // interpolation points, on the matrices is the identity + // matrix, and we have to invert only one matrix. The + // resulting vector will be interpolatory at the support + // points of fe_q, even if the finite element does not have + // support points. // // Morally, we should invert the matrix T_ij = phi_i(q_j), // however in general this matrix is not invertible, since - // there may by components which do not contribute to the + // there may be components which do not contribute to the // displacement vector. Since we are not interested in those // components, we construct a square matrix with the same // number of components of the FE_Q system. The FE_Q system @@ -7239,7 +7251,7 @@ namespace VectorTools for (unsigned int j=0; j