Reference documentation for deal.II version 8.5.1
mapping_q.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2017 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;
122 
127  virtual
130  const typename Triangulation<dim,spacedim>::cell_iterator &cell,
131  const Point<dim> &p) const;
132 
154  virtual
155  Point<dim>
157  const Point<spacedim> &p) const;
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;
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;
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;
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;
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;
198 
203  virtual
204  Mapping<dim,spacedim> *clone () const;
205 
206 
212 protected:
213 
234  class InternalData : public Mapping<dim,spacedim>::InternalDataBase
235  {
236  public:
240  InternalData ();
241 
242 
246  virtual std::size_t memory_consumption () const;
247 
255 
260  std_cxx11::unique_ptr<typename MappingQGeneric<dim,spacedim>::InternalData> mapping_q1_data;
261 
266  std_cxx11::unique_ptr<typename MappingQGeneric<dim,spacedim>::InternalData> mapping_qp_data;
267  };
268 
269 protected:
270 
271  // documentation can be found in Mapping::requires_update_flags()
272  virtual
274  requires_update_flags (const UpdateFlags update_flags) const;
275 
276  // documentation can be found in Mapping::get_data()
277  virtual
278  InternalData *
279  get_data (const UpdateFlags,
280  const Quadrature<dim> &quadrature) const;
281 
282  // documentation can be found in Mapping::get_face_data()
283  virtual
284  InternalData *
285  get_face_data (const UpdateFlags flags,
286  const Quadrature<dim-1>& quadrature) const;
287 
288  // documentation can be found in Mapping::get_subface_data()
289  virtual
290  InternalData *
291  get_subface_data (const UpdateFlags flags,
292  const Quadrature<dim-1>& quadrature) const;
293 
294  // documentation can be found in Mapping::fill_fe_values()
295  virtual
297  fill_fe_values (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 
303  // documentation can be found in Mapping::fill_fe_face_values()
304  virtual void
305  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
306  const unsigned int face_no,
307  const Quadrature<dim-1> &quadrature,
308  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
310 
311  // documentation can be found in Mapping::fill_fe_subface_values()
312  virtual void
313  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
314  const unsigned int face_no,
315  const unsigned int subface_no,
316  const Quadrature<dim-1> &quadrature,
317  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
319 
324 protected:
325 
330  const unsigned int polynomial_degree;
331 
337 
355  std_cxx11::shared_ptr<const MappingQGeneric<dim,spacedim> > q1_mapping;
356 
371  std_cxx11::shared_ptr<const MappingQGeneric<dim,spacedim> > qp_mapping;
372 };
373 
376 DEAL_II_NAMESPACE_CLOSE
377 
378 #endif
unsigned int get_degree() const
Definition: mapping_q.cc:107
virtual InternalData * get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
Definition: mapping_q.cc:161
std_cxx11::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
Definition: mapping_q.h:355
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const
Definition: mapping_q.cc:488
MappingType
Definition: mapping.h:50
virtual InternalData * get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
Definition: mapping_q.cc:185
bool use_mapping_q1_on_current_cell
Definition: mapping_q.h:254
std_cxx11::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_q1_data
Definition: mapping_q.h:260
const bool use_mapping_q_on_all_cells
Definition: mapping_q.h:336
virtual std::size_t memory_consumption() const
Definition: mapping_q.cc:45
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
Definition: mapping_q.cc:126
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
Definition: mapping_q.cc:352
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
const unsigned int polynomial_degree
Definition: mapping_q.h:330
virtual Mapping< dim, spacedim > * clone() const
Definition: mapping_q.cc:505
std_cxx11::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_qp_data
Definition: mapping_q.h:266
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const
Definition: mapping_q.cc:471
virtual InternalData * get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const
Definition: mapping_q.cc:137
std_cxx11::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:371
MappingQ(const unsigned int polynomial_degree, const bool use_mapping_q_on_all_cells=false)
Definition: mapping_q.cc:56
virtual bool preserves_vertex_locations() const
Definition: mapping_q.cc:117