#include <deal.II/fe/fe.h>
#include <deal.II/fe/fe_tools.h>
#include <deal.II/fe/mapping_q_eulerian.h>
+#include <deal.II/fe/mapping_q1_eulerian.h>
DEAL_II_NAMESPACE_OPEN
fe_values(euler_dof_handler.get_fe(),
support_quadrature,
update_values | update_q_points)
-{ }
+{
+ // reset the q1 mapping we use for interior cells (and previously
+ // set by the MappingQ constructor) to a MappingQ1Eulerian with the
+ // current vector
+ this->q1_mapping.reset (new MappingQ1Eulerian<dim,EulerVectorType,spacedim>(euler_vector,
+ euler_dof_handler));
+}
fe_values(euler_dof_handler.get_fe(),
support_quadrature,
update_values | update_q_points)
-{ }
+{
+ // reset the q1 mapping we use for interior cells (and previously
+ // set by the MappingQ constructor) to a MappingQ1Eulerian with the
+ // current vector
+ this->q1_mapping.reset (new MappingQ1Eulerian<dim,EulerVectorType,spacedim>(euler_vector,
+ euler_dof_handler));
+}
for (unsigned int i=1; i<dpo.size(); ++i)
dpo[i]=dpo[i-1]*(map_degree-1);
- FETools::lexicographic_to_hierarchic_numbering (
- FiniteElementData<dim> (dpo, 1, map_degree), renumber);
+ FETools::lexicographic_to_hierarchic_numbering (FiniteElementData<dim> (dpo, 1, map_degree),
+ renumber);
// finally we assign the quadrature points in the required order.
for (unsigned int q=0; q<n_q_points; ++q)
template <int dim, class EulerVectorType, int spacedim>
void
MappingQEulerian<dim, EulerVectorType, spacedim>::
-compute_mapping_support_points
-(const typename Triangulation<dim,spacedim>::cell_iterator &cell,
- std::vector<Point<spacedim> > &a) const
+compute_mapping_support_points (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
+ std::vector<Point<spacedim> > &a) const
{
// first, basic assertion with respect to vector size,
Assert (n_components >= spacedim, ExcDimensionMismatch(n_components, spacedim) );
- std::vector<Vector<double> > shift_vector(n_support_pts,Vector<double>(n_components));
+ std::vector<Vector<typename EulerVectorType::value_type> >
+ shift_vector(n_support_pts,
+ Vector<typename EulerVectorType::value_type>(n_components));
// fill shift vector for each support point using an fe_values object. make
// sure that the fe_values variable isn't used simultaneously from different
// ---------------------------------------------------------------------
//
-// Copyright (C) 1998 - 2014 by the deal.II authors
+// Copyright (C) 1998 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
{
#if deal_II_dimension <= deal_II_space_dimension
template class MappingQEulerian<deal_II_dimension, Vector<double>, deal_II_space_dimension>;
+ template class MappingQEulerian<deal_II_dimension, Vector<float>, deal_II_space_dimension>;
# ifdef DEAL_II_WITH_PETSC
template class MappingQEulerian<deal_II_dimension,
PETScWrappers::Vector, deal_II_space_dimension>;