47template <
int dim,
class VectorType,
int spacedim>
50 const unsigned int degree,
53 , mapping_q_eulerian(mapping_q_eulerian)
54 , support_quadrature(degree)
55 , fe_values(mapping_q_eulerian.euler_dof_handler->get_fe(),
62template <
int dim,
class VectorType,
int spacedim>
64 const unsigned int degree,
66 const VectorType & euler_vector,
67 const unsigned int level)
68 :
MappingQ<dim, spacedim>(degree, true)
69 , euler_vector(&euler_vector)
70 , euler_dof_handler(&euler_dof_handler)
77 std::make_shared<MappingQ1Eulerian<dim, VectorType, spacedim>>(
78 euler_dof_handler, euler_vector);
81 this->qp_mapping = std::make_shared<MappingQEulerianGeneric>(degree, *
this);
86template <
int dim,
class VectorType,
int spacedim>
87std::unique_ptr<Mapping<dim, spacedim>>
90 return std::make_unique<MappingQEulerian<dim, VectorType, spacedim>>(
98template <
int dim,
class VectorType,
int spacedim>
107 const unsigned int n_q_points = q_iterated.size();
111 const std::vector<unsigned int> renumber =
112 FETools::lexicographic_to_hierarchic_numbering<dim>(map_degree);
115 for (
unsigned int q = 0; q < n_q_points; ++q)
123template <
int dim,
class VectorType,
int spacedim>
124boost::container::small_vector<Point<spacedim>,
130 const std::vector<Point<spacedim>> a =
131 dynamic_cast<const MappingQEulerianGeneric &
>(*this->
qp_mapping)
132 .compute_mapping_support_points(cell);
134 boost::container::small_vector<Point<spacedim>,
136 vertex_locations(a.begin(),
139 return vertex_locations;
144template <
int dim,
class VectorType,
int spacedim>
145boost::container::small_vector<Point<spacedim>,
151 return mapping_q_eulerian.get_vertices(cell);
156template <
int dim,
class VectorType,
int spacedim>
166template <
int dim,
class VectorType,
int spacedim>
167std::vector<Point<spacedim>>
172 const bool mg_vector =
177 mapping_q_eulerian.euler_dof_handler->n_dofs(mapping_q_eulerian.level) :
178 mapping_q_eulerian.euler_dof_handler->n_dofs();
180 mapping_q_eulerian.euler_vector->size();
190 *cell, mapping_q_eulerian.euler_dof_handler);
192 Assert(mg_vector || dof_cell->is_active() ==
true, ExcInactiveCell());
205 const unsigned int n_support_pts = support_quadrature.size();
206 const unsigned int n_components =
207 mapping_q_eulerian.euler_dof_handler->get_fe(0).n_components();
209 Assert(n_components >= spacedim,
212 std::vector<Vector<typename VectorType::value_type>> shift_vector(
215 std::vector<types::global_dof_index> dof_indices(
216 mapping_q_eulerian.euler_dof_handler->get_fe(0).n_dofs_per_cell());
220 std::lock_guard<std::mutex> lock(fe_values_mutex);
221 fe_values.reinit(dof_cell);
224 dof_cell->get_mg_dof_indices(dof_indices);
225 fe_values.get_function_values(*mapping_q_eulerian.euler_vector,
230 fe_values.get_function_values(*mapping_q_eulerian.euler_vector,
235 std::vector<Point<spacedim>> a(n_support_pts);
236 for (
unsigned int q = 0; q < n_support_pts; ++q)
238 a[q] = fe_values.quadrature_point(q);
239 for (
unsigned int d = 0;
d < spacedim; ++
d)
240 a[q](
d) += shift_vector[q](
d);
248template <
int dim,
class VectorType,
int spacedim>
274template <
int dim,
class VectorType,
int spacedim>
280 dynamic_cast<const MappingQEulerianGeneric &
>(*this->
qp_mapping)
281 .compute_mapping_support_points(cell));
287#include "mapping_q_eulerian.inst"
SupportQuadrature(const unsigned int map_degree)
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
virtual bool preserves_vertex_locations() const override
MappingQEulerianGeneric(const unsigned int degree, const MappingQEulerian< dim, VectorType, spacedim > &mapping_q_eulerian)
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim > > euler_dof_handler
virtual BoundingBox< spacedim > get_bounding_box(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::unique_ptr< Mapping< dim, spacedim > > clone() const override
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
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
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 > > qp_mapping
unsigned int get_degree() const
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ update_values
Shape function values.
@ update_quadrature_points
Transformed quadrature points.
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const std::vector< std::vector< double > > &properties={})
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static const unsigned int invalid_unsigned_int