16 #include <deal.II/base/utilities.h> 17 #include <deal.II/base/quadrature_lib.h> 18 #include <deal.II/base/std_cxx14/memory.h> 19 #include <deal.II/lac/vector.h> 20 #include <deal.II/lac/block_vector.h> 21 #include <deal.II/lac/la_vector.h> 22 #include <deal.II/lac/parallel_vector.h> 23 #include <deal.II/lac/parallel_block_vector.h> 24 #include <deal.II/lac/petsc_parallel_vector.h> 25 #include <deal.II/lac/petsc_parallel_block_vector.h> 26 #include <deal.II/lac/trilinos_vector.h> 27 #include <deal.II/lac/trilinos_parallel_block_vector.h> 28 #include <deal.II/lac/trilinos_parallel_block_vector.h> 29 #include <deal.II/grid/tria_iterator.h> 30 #include <deal.II/dofs/dof_handler.h> 31 #include <deal.II/dofs/dof_accessor.h> 32 #include <deal.II/fe/fe.h> 33 #include <deal.II/fe/fe_tools.h> 34 #include <deal.II/fe/mapping_q_eulerian.h> 35 #include <deal.II/fe/mapping_q1_eulerian.h> 37 DEAL_II_NAMESPACE_OPEN
43 template <
int dim,
class VectorType,
int spacedim>
49 mapping_q_eulerian (mapping_q_eulerian),
50 support_quadrature(degree),
51 fe_values(mapping_q_eulerian.euler_dof_handler->get_fe(),
58 template <
int dim,
class VectorType,
int spacedim>
63 const unsigned int level)
65 MappingQ<dim,spacedim>(degree, true),
73 this->
q1_mapping = std::make_shared<MappingQ1Eulerian<dim,VectorType,spacedim>>
77 this->
qp_mapping = std::make_shared<MappingQEulerianGeneric>(degree,*
this);
82 template <
int dim,
class VectorType,
int spacedim>
83 std::unique_ptr<Mapping<dim,spacedim> >
86 return std_cxx14::make_unique<MappingQEulerian<dim,VectorType,spacedim>>
87 (this->get_degree(), *euler_dof_handler, *euler_vector);
94 template <
int dim,
class VectorType,
int spacedim>
105 const unsigned int n_q_points = q_iterated.
size();
110 std::vector<unsigned int> renumber(n_q_points);
111 std::vector<unsigned int> dpo(dim+1, 1U);
112 for (
unsigned int i=1; i<dpo.size(); ++i)
113 dpo[i]=dpo[i-1]*(map_degree-1);
119 for (
unsigned int q=0; q<n_q_points; ++q)
127 template <
int dim,
class VectorType,
int spacedim>
134 const std::vector<Point<spacedim> > a
138 std::copy (a.begin(),
140 vertex_locations.begin());
142 return vertex_locations;
147 template <
int dim,
class VectorType,
int spacedim>
158 template <
int dim,
class VectorType,
int spacedim>
168 template <
int dim,
class VectorType,
int spacedim>
169 std::vector<Point<spacedim> >
204 const unsigned int n_components =
mapping_q_eulerian.euler_dof_handler->get_fe(0).n_components();
208 std::vector<Vector<typename VectorType::value_type> >
209 shift_vector(n_support_pts,
212 std::vector<types::global_dof_index> dof_indices(
mapping_q_eulerian.euler_dof_handler->get_fe(0).dofs_per_cell);
220 dof_cell->get_mg_dof_indices(dof_indices);
230 std::vector<Point<spacedim> > a(n_support_pts);
231 for (
unsigned int q=0; q<n_support_pts; ++q)
234 for (
unsigned int d=0; d<spacedim; ++d)
235 a[q](d) += shift_vector[q](d);
243 template <
int dim,
class VectorType,
int spacedim>
269 #include "mapping_q_eulerian.inst" 272 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
#define AssertDimension(dim1, dim2)
FEValues< dim, spacedim > fe_values
Threads::Mutex fe_values_mutex
static ::ExceptionBase & ExcInactiveCell()
const Point< dim > & point(const unsigned int i) const
Transformed quadrature points.
const SupportQuadrature support_quadrature
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
std::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
ActiveSelector::cell_iterator cell_iterator
MappingQEulerianGeneric(const unsigned int degree, const MappingQEulerian< dim, VectorType, spacedim > &mapping_q_eulerian)
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
unsigned int global_dof_index
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
MappingQEulerian(const unsigned int degree, const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector, const unsigned int level=numbers::invalid_unsigned_int)
SmartPointer< const VectorType, MappingQEulerian< dim, VectorType, spacedim > > euler_vector
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
std::vector< Point< dim > > quadrature_points
unsigned int size() const
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim > > euler_dof_handler
virtual bool preserves_vertex_locations() const override
SupportQuadrature(const unsigned int map_degree)
const MappingQEulerian< dim, VectorType, spacedim > & mapping_q_eulerian