Reference documentation for deal.II version Git 5b8b897cb2 2021-04-22 22:24:19 -0400
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mapping_q.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2019 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 
23 
25 
28 
79 template <int dim, int spacedim = dim>
80 class MappingQ : public Mapping<dim, spacedim>
81 {
82 public:
100  MappingQ(const unsigned int polynomial_degree,
101  const bool use_mapping_q_on_all_cells = false);
102 
106  MappingQ(const MappingQ<dim, spacedim> &mapping);
107 
112  unsigned int
113  get_degree() const;
114 
119  virtual bool
120  preserves_vertex_locations() const override;
121 
122  // for documentation, see the Mapping base class
123  virtual BoundingBox<spacedim>
125  &cell) const override;
126 
127  virtual bool
128  is_compatible_with(const ReferenceCell &reference_cell) const override;
129 
134  virtual Point<spacedim>
136  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
137  const Point<dim> &p) const override;
138 
160  virtual Point<dim>
162  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
163  const Point<spacedim> &p) const override;
164 
165  // for documentation, see the Mapping base class
166  virtual void
167  transform(const ArrayView<const Tensor<1, dim>> & input,
168  const MappingKind kind,
169  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
170  const ArrayView<Tensor<1, spacedim>> &output) const override;
171 
172  // for documentation, see the Mapping base class
173  virtual void
175  const MappingKind kind,
176  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
177  const ArrayView<Tensor<2, spacedim>> &output) const override;
178 
179  // for documentation, see the Mapping base class
180  virtual void
181  transform(const ArrayView<const Tensor<2, dim>> & input,
182  const MappingKind kind,
183  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
184  const ArrayView<Tensor<2, spacedim>> &output) const override;
185 
186  // for documentation, see the Mapping base class
187  virtual void
189  const MappingKind kind,
190  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
191  const ArrayView<Tensor<3, spacedim>> &output) const override;
192 
193  // for documentation, see the Mapping base class
194  virtual void
195  transform(const ArrayView<const Tensor<3, dim>> & input,
196  const MappingKind kind,
197  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
198  const ArrayView<Tensor<3, spacedim>> &output) const override;
199 
205  virtual std::unique_ptr<Mapping<dim, spacedim>>
206  clone() const override;
207 
208 
214 protected:
235  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
236  {
237  public:
241  InternalData();
242 
243 
247  virtual std::size_t
248  memory_consumption() const override;
249 
257 
262  std::unique_ptr<typename MappingQGeneric<dim, spacedim>::InternalData>
264 
269  std::unique_ptr<typename MappingQGeneric<dim, spacedim>::InternalData>
271  };
272 
273 protected:
274  // documentation can be found in Mapping::requires_update_flags()
275  virtual UpdateFlags
276  requires_update_flags(const UpdateFlags update_flags) const override;
277 
278  // documentation can be found in Mapping::get_data()
279  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
280  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
281 
283 
284  // documentation can be found in Mapping::get_face_data()
285  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
286  get_face_data(const UpdateFlags flags,
287  const hp::QCollection<dim - 1> &quadrature) const override;
288 
289  // documentation can be found in Mapping::get_subface_data()
290  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
291  get_subface_data(const UpdateFlags flags,
292  const Quadrature<dim - 1> &quadrature) const override;
293 
294  // documentation can be found in Mapping::fill_fe_values()
297  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
298  const CellSimilarity::Similarity cell_similarity,
299  const Quadrature<dim> & quadrature,
300  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
302  &output_data) const override;
303 
305 
306  // documentation can be found in Mapping::fill_fe_face_values()
307  virtual void
309  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
310  const unsigned int face_no,
311  const hp::QCollection<dim - 1> & quadrature,
312  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
314  &output_data) const override;
315 
316  // documentation can be found in Mapping::fill_fe_subface_values()
317  virtual void
319  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
320  const unsigned int face_no,
321  const unsigned int subface_no,
322  const Quadrature<dim - 1> & quadrature,
323  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
325  &output_data) const override;
326 
331 protected:
336  const unsigned int polynomial_degree;
337 
343 
361  std::shared_ptr<const MappingQGeneric<dim, spacedim>> q1_mapping;
362 
377  std::shared_ptr<const MappingQGeneric<dim, spacedim>> qp_mapping;
378 };
379 
383 
384 #endif
unsigned int get_degree() const
Definition: mapping_q.cc:122
virtual std::size_t memory_consumption() const override
Definition: mapping_q.cc:49
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:377
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Definition: mapping_q.cc:150
bool use_mapping_q1_on_current_cell
Definition: mapping_q.h:256
virtual void transform(const ArrayView< const Tensor< 1, dim >> &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim >> &output) const override
Definition: mapping_q.cc:388
const bool use_mapping_q_on_all_cells
Definition: mapping_q.h:342
MappingKind
Definition: mapping.h:64
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
Definition: mapping_q.cc:181
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
Definition: mapping_q.cc:342
UpdateFlags
void reference_cell(const ReferenceCell &reference_cell, Triangulation< dim, spacedim > &tria)
Abstract base class for mapping classes.
Definition: mapping.h:303
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:396
const unsigned int polynomial_degree
Definition: mapping_q.h:336
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
Definition: mapping_q.cc:549
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_qp_data
Definition: mapping_q.h:270
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:218
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
Definition: mapping_q.cc:564
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:507
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
Definition: mapping_q.cc:306
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:395
virtual bool preserves_vertex_locations() const override
Definition: mapping_q.cc:131
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_q.cc:539
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
Definition: mapping_q.cc:252
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_q1_data
Definition: mapping_q.h:263
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Definition: mapping_q.cc:140
MappingQ(const unsigned int polynomial_degree, const bool use_mapping_q_on_all_cells=false)
Definition: mapping_q.cc:61
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:524
std::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
Definition: mapping_q.h:361