Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
mapping_q.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_mapping_q_h
17 #define dealii_mapping_q_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/fe/mapping_q_generic.h>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
28 
82 template <int dim, int spacedim = dim>
83 class MappingQ : public Mapping<dim, spacedim>
84 {
85 public:
103  MappingQ(const unsigned int polynomial_degree,
104  const bool use_mapping_q_on_all_cells = false);
105 
109  MappingQ(const MappingQ<dim, spacedim> &mapping);
110 
115  unsigned int
116  get_degree() const;
117 
122  virtual bool
123  preserves_vertex_locations() const override;
124 
129  virtual Point<spacedim>
131  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
132  const Point<dim> &p) const override;
133 
155  virtual Point<dim>
157  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
158  const Point<spacedim> &p) const override;
159 
160  // for documentation, see the Mapping base class
161  virtual void
162  transform(const ArrayView<const Tensor<1, dim>> & input,
163  const MappingType type,
164  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
165  const ArrayView<Tensor<1, spacedim>> &output) const override;
166 
167  // for documentation, see the Mapping base class
168  virtual void
170  const MappingType type,
171  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
172  const ArrayView<Tensor<2, spacedim>> &output) const override;
173 
174  // for documentation, see the Mapping base class
175  virtual void
176  transform(const ArrayView<const Tensor<2, dim>> & input,
177  const MappingType type,
178  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
179  const ArrayView<Tensor<2, spacedim>> &output) const override;
180 
181  // for documentation, see the Mapping base class
182  virtual void
184  const MappingType type,
185  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
186  const ArrayView<Tensor<3, spacedim>> &output) const override;
187 
188  // for documentation, see the Mapping base class
189  virtual void
190  transform(const ArrayView<const Tensor<3, dim>> & input,
191  const MappingType type,
192  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
193  const ArrayView<Tensor<3, spacedim>> &output) const override;
194 
200  virtual std::unique_ptr<Mapping<dim, spacedim>>
201  clone() const override;
202 
203 
209 protected:
230  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
231  {
232  public:
236  InternalData();
237 
238 
242  virtual std::size_t
243  memory_consumption() const override;
244 
252 
257  std::unique_ptr<typename MappingQGeneric<dim, spacedim>::InternalData>
259 
264  std::unique_ptr<typename MappingQGeneric<dim, spacedim>::InternalData>
266  };
267 
268 protected:
269  // documentation can be found in Mapping::requires_update_flags()
270  virtual UpdateFlags
271  requires_update_flags(const UpdateFlags update_flags) const override;
272 
273  // documentation can be found in Mapping::get_data()
274  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
275  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
276 
277  // documentation can be found in Mapping::get_face_data()
278  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
279  get_face_data(const UpdateFlags flags,
280  const Quadrature<dim - 1> &quadrature) const override;
281 
282  // documentation can be found in Mapping::get_subface_data()
283  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
284  get_subface_data(const UpdateFlags flags,
285  const Quadrature<dim - 1> &quadrature) const override;
286 
287  // documentation can be found in Mapping::fill_fe_values()
289  fill_fe_values(
290  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
291  const CellSimilarity::Similarity cell_similarity,
292  const Quadrature<dim> & quadrature,
293  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
295  &output_data) const override;
296 
297  // documentation can be found in Mapping::fill_fe_face_values()
298  virtual void
299  fill_fe_face_values(
300  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
301  const unsigned int face_no,
302  const Quadrature<dim - 1> & quadrature,
303  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
305  &output_data) const override;
306 
307  // documentation can be found in Mapping::fill_fe_subface_values()
308  virtual void
309  fill_fe_subface_values(
310  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
311  const unsigned int face_no,
312  const unsigned int subface_no,
313  const Quadrature<dim - 1> & quadrature,
314  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
316  &output_data) const override;
317 
322 protected:
327  const unsigned int polynomial_degree;
328 
334 
352  std::shared_ptr<const MappingQGeneric<dim, spacedim>> q1_mapping;
353 
368  std::shared_ptr<const MappingQGeneric<dim, spacedim>> qp_mapping;
369 };
370 
373 DEAL_II_NAMESPACE_CLOSE
374 
375 #endif
unsigned int get_degree() const
Definition: mapping_q.cc:123
virtual std::size_t memory_consumption() const override
Definition: mapping_q.cc:50
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:368
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Definition: mapping_q.cc:151
MappingType
Definition: mapping.h:61
bool use_mapping_q1_on_current_cell
Definition: mapping_q.h:251
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
Definition: mapping_q.cc:182
const bool use_mapping_q_on_all_cells
Definition: mapping_q.h:333
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:57
const unsigned int polynomial_degree
Definition: mapping_q.h:327
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_qp_data
Definition: mapping_q.h:265
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
Definition: mapping_q.cc:215
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
Definition: mapping_q.cc:504
virtual bool preserves_vertex_locations() const override
Definition: mapping_q.cc:132
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
Definition: mapping_q.cc:385
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_q.cc:536
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_q1_data
Definition: mapping_q.h:258
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Definition: mapping_q.cc:141
MappingQ(const unsigned int polynomial_degree, const bool use_mapping_q_on_all_cells=false)
Definition: mapping_q.cc:62
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
Definition: mapping_q.cc:521
std::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
Definition: mapping_q.h:352