16 #include <deal.II/base/array_view.h> 17 #include <deal.II/base/polynomial.h> 18 #include <deal.II/base/quadrature.h> 19 #include <deal.II/base/quadrature_lib.h> 20 #include <deal.II/base/memory_consumption.h> 21 #include <deal.II/base/std_cxx14/memory.h> 22 #include <deal.II/base/tensor_product_polynomials.h> 23 #include <deal.II/base/utilities.h> 24 #include <deal.II/lac/full_matrix.h> 25 #include <deal.II/grid/tria_iterator.h> 26 #include <deal.II/dofs/dof_accessor.h> 27 #include <deal.II/fe/mapping_q.h> 28 #include <deal.II/fe/fe_q.h> 29 #include <deal.II/fe/fe_values.h> 34 DEAL_II_NAMESPACE_OPEN
37 template <
int dim,
int spacedim>
40 use_mapping_q1_on_current_cell(false)
45 template <
int dim,
int spacedim>
57 template <
int dim,
int spacedim>
88 template <
int dim,
int spacedim>
91 polynomial_degree (mapping.polynomial_degree),
92 use_mapping_q_on_all_cells (mapping.use_mapping_q_on_all_cells)
96 std::shared_ptr<const Mapping<dim,spacedim>> other_q1_map = mapping.
q1_mapping->clone();
108 std::shared_ptr<const Mapping<dim,spacedim>> other_qp_map = mapping.
qp_mapping->clone();
116 template <
int dim,
int spacedim>
120 return polynomial_degree;
125 template <
int dim,
int spacedim>
135 template <
int dim,
int spacedim>
139 return (q1_mapping->requires_update_flags(in)
141 qp_mapping->requires_update_flags(in));
146 template <
int dim,
int spacedim>
147 std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
151 auto data = std_cxx14::make_unique<InternalData>();
160 if (!use_mapping_q_on_all_cells)
162 (std::move(q1_mapping->get_data (update_flags, quadrature)));
166 (std::move(do_get_data.return_value()));
167 return std::move(data);
172 template <
int dim,
int spacedim>
173 std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
177 auto data = std_cxx14::make_unique<InternalData>();
186 if (!use_mapping_q_on_all_cells)
188 (std::move(q1_mapping->get_face_data (update_flags, quadrature)));
192 (std::move(do_get_data.return_value()));
193 return std::move(data);
198 template <
int dim,
int spacedim>
199 std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
203 auto data = std_cxx14::make_unique<InternalData>();
212 if (!use_mapping_q_on_all_cells)
214 (std::move(q1_mapping->get_subface_data (update_flags, quadrature)));
218 (std::move(do_get_data.return_value()));
219 return std::move(data);
225 template <
int dim,
int spacedim>
237 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
241 data.use_mapping_q1_on_current_cell = !(use_mapping_q_on_all_cells
242 || cell->has_boundary_lines());
252 = ((data.use_mapping_q1_on_current_cell ==
false)
254 (this->polynomial_degree > 1)
262 if (data.use_mapping_q1_on_current_cell)
263 q1_mapping->fill_fe_values (cell,
264 updated_cell_similarity,
266 *data.mapping_q1_data,
269 qp_mapping->fill_fe_values(cell,
270 updated_cell_similarity,
272 *data.mapping_qp_data,
275 return updated_cell_similarity;
280 template <
int dim,
int spacedim>
284 const unsigned int face_no,
291 Assert (dynamic_cast<const InternalData *> (&internal_data) !=
nullptr,
293 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
301 data.use_mapping_q1_on_current_cell = !(use_mapping_q_on_all_cells
302 || cell->has_boundary_lines());
306 if (data.use_mapping_q1_on_current_cell)
307 q1_mapping->fill_fe_face_values (cell,
310 *data.mapping_q1_data,
313 qp_mapping->fill_fe_face_values(cell,
316 *data.mapping_qp_data,
321 template <
int dim,
int spacedim>
325 const unsigned int face_no,
326 const unsigned int subface_no,
333 Assert (dynamic_cast<const InternalData *> (&internal_data) !=
nullptr,
335 const InternalData &data =
static_cast<const InternalData &
> (internal_data);
343 data.use_mapping_q1_on_current_cell = !(use_mapping_q_on_all_cells
344 || cell->has_boundary_lines());
348 if (data.use_mapping_q1_on_current_cell)
349 q1_mapping->fill_fe_subface_values (cell,
353 *data.mapping_q1_data,
356 qp_mapping->fill_fe_subface_values(cell,
360 *data.mapping_qp_data,
366 template <
int dim,
int spacedim>
380 if (data->use_mapping_q1_on_current_cell)
381 q1_mapping->transform (input, mapping_type, *data->mapping_q1_data, output);
384 qp_mapping->transform(input, mapping_type, *data->mapping_qp_data, output);
389 template <
int dim,
int spacedim>
405 q1_mapping->transform (input, mapping_type, *data->
mapping_q1_data, output);
408 qp_mapping->transform(input, mapping_type, *data->
mapping_qp_data, output);
413 template <
int dim,
int spacedim>
429 q1_mapping->transform (input, mapping_type, *data->
mapping_q1_data, output);
432 qp_mapping->transform(input, mapping_type, *data->
mapping_qp_data, output);
437 template <
int dim,
int spacedim>
453 q1_mapping->transform (input, mapping_type, *data->
mapping_q1_data, output);
456 qp_mapping->transform(input, mapping_type, *data->
mapping_qp_data, output);
461 template <
int dim,
int spacedim>
476 q1_mapping->transform (input, mapping_type, *data->
mapping_q1_data, output);
479 qp_mapping->transform(input, mapping_type, *data->
mapping_qp_data, output);
484 template <
int dim,
int spacedim>
493 if (use_mapping_q_on_all_cells || cell->has_boundary_lines())
494 return qp_mapping->transform_unit_to_real_cell (cell, p);
496 return q1_mapping->transform_unit_to_real_cell (cell, p);
501 template <
int dim,
int spacedim>
507 if (cell->has_boundary_lines()
509 use_mapping_q_on_all_cells
512 return qp_mapping->transform_real_to_unit_cell(cell, p);
514 return q1_mapping->transform_real_to_unit_cell(cell, p);
519 template <
int dim,
int spacedim>
520 std::unique_ptr<Mapping<dim,spacedim> >
523 return std_cxx14::make_unique<MappingQ<dim,spacedim>>(this->polynomial_degree,
524 this->use_mapping_q_on_all_cells);
530 #include "mapping_q.inst" 533 DEAL_II_NAMESPACE_CLOSE
unsigned int get_degree() const
virtual std::size_t memory_consumption() const override
#define AssertDimension(dim1, dim2)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Task< RT > new_task(const std::function< RT()> &function)
bool use_mapping_q1_on_current_cell
const bool use_mapping_q_on_all_cells
std::unique_ptr< To > dynamic_unique_cast(std::unique_ptr< From > &&p)
std::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
#define Assert(cond, exc)
Abstract base class for mapping classes.
const unsigned int polynomial_degree
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_q1_data
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_qp_data
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
virtual bool preserves_vertex_locations() const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
MappingQ(const unsigned int polynomial_degree, const bool use_mapping_q_on_all_cells=false)
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()