Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/fe/mapping_q_generic.h>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
27 
81 template <int dim, int spacedim=dim>
82 class MappingQ : public Mapping<dim,spacedim>
83 {
84 public:
102  MappingQ (const unsigned int polynomial_degree,
103  const bool use_mapping_q_on_all_cells = false);
104 
108  MappingQ (const MappingQ<dim,spacedim> &mapping);
109 
114  unsigned int get_degree () const;
115 
120  virtual
121  bool preserves_vertex_locations () const override;
122 
127  virtual
130  const typename Triangulation<dim,spacedim>::cell_iterator &cell,
131  const Point<dim> &p) const override;
132 
154  virtual
155  Point<dim>
157  const Point<spacedim> &p) const override;
158 
159  // for documentation, see the Mapping base class
160  virtual
161  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
169  void
171  const MappingType type,
172  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
173  const ArrayView<Tensor<2,spacedim> > &output) const override;
174 
175  // for documentation, see the Mapping base class
176  virtual
177  void
178  transform (const ArrayView<const Tensor<2, dim> > &input,
179  const MappingType type,
180  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
181  const ArrayView<Tensor<2,spacedim> > &output) const override;
182 
183  // for documentation, see the Mapping base class
184  virtual
185  void
187  const MappingType type,
188  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
189  const ArrayView<Tensor<3,spacedim> > &output) const override;
190 
191  // for documentation, see the Mapping base class
192  virtual
193  void
194  transform (const ArrayView<const Tensor<3, dim> > &input,
195  const MappingType type,
196  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
197  const ArrayView<Tensor<3,spacedim> > &output) const override;
198 
204  virtual
205  std::unique_ptr<Mapping<dim,spacedim>> clone () const override;
206 
207 
213 protected:
214 
235  class InternalData : public Mapping<dim,spacedim>::InternalDataBase
236  {
237  public:
241  InternalData ();
242 
243 
247  virtual
248  std::size_t memory_consumption () const override;
249 
257 
262  std::unique_ptr<typename MappingQGeneric<dim,spacedim>::InternalData> mapping_q1_data;
263 
268  std::unique_ptr<typename MappingQGeneric<dim,spacedim>::InternalData> mapping_qp_data;
269  };
270 
271 protected:
272 
273  // documentation can be found in Mapping::requires_update_flags()
274  virtual
276  requires_update_flags (const UpdateFlags update_flags) const override;
277 
278  // documentation can be found in Mapping::get_data()
279  virtual
280  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
281  get_data (const UpdateFlags,
282  const Quadrature<dim> &quadrature) const override;
283 
284  // documentation can be found in Mapping::get_face_data()
285  virtual
286  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
287  get_face_data (const UpdateFlags flags,
288  const Quadrature<dim-1>& quadrature) const override;
289 
290  // documentation can be found in Mapping::get_subface_data()
291  virtual
292  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
293  get_subface_data (const UpdateFlags flags,
294  const Quadrature<dim-1>& quadrature) const override;
295 
296  // documentation can be found in Mapping::fill_fe_values()
297  virtual
299  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
300  const CellSimilarity::Similarity cell_similarity,
301  const Quadrature<dim> &quadrature,
302  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
304 
305  // documentation can be found in Mapping::fill_fe_face_values()
306  virtual
307  void
308  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
309  const unsigned int face_no,
310  const Quadrature<dim-1> &quadrature,
311  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
313 
314  // documentation can be found in Mapping::fill_fe_subface_values()
315  virtual
316  void
317  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
318  const unsigned int face_no,
319  const unsigned int subface_no,
320  const Quadrature<dim-1> &quadrature,
321  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
323 
328 protected:
329 
334  const unsigned int polynomial_degree;
335 
341 
359  std::shared_ptr<const MappingQGeneric<dim,spacedim> > q1_mapping;
360 
375  std::shared_ptr<const MappingQGeneric<dim,spacedim> > qp_mapping;
376 };
377 
380 DEAL_II_NAMESPACE_CLOSE
381 
382 #endif
unsigned int get_degree() const
Definition: mapping_q.cc:118
virtual std::size_t memory_consumption() const override
Definition: mapping_q.cc:47
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Definition: mapping_q.cc:148
MappingType
Definition: mapping.h:51
bool use_mapping_q1_on_current_cell
Definition: mapping_q.h:256
const bool use_mapping_q_on_all_cells
Definition: mapping_q.h:340
std::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
Definition: mapping_q.h:359
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
const unsigned int polynomial_degree
Definition: mapping_q.h:334
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_q1_data
Definition: mapping_q.h:262
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:487
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_qp_data
Definition: mapping_q.h:268
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:375
virtual bool preserves_vertex_locations() const override
Definition: mapping_q.cc:128
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_q.cc:521
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:369
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:174
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:200
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Definition: mapping_q.cc:137
MappingQ(const unsigned int polynomial_degree, const bool use_mapping_q_on_all_cells=false)
Definition: mapping_q.cc:58
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:504