Reference documentation for deal.II version 8.5.1
mapping_fe_field.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2016 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_fe_h
17 #define dealii__mapping_fe_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/table.h>
22 #include <deal.II/base/thread_management.h>
23 #include <deal.II/lac/vector.h>
24 #include <deal.II/dofs/dof_handler.h>
25 #include <deal.II/fe/mapping.h>
26 #include <deal.II/fe/fe.h>
27 
28 #include <deal.II/base/std_cxx11/array.h>
29 
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 
36 
74 template <int dim, int spacedim=dim,
75  typename VectorType=Vector<double>,
76  typename DoFHandlerType=DoFHandler<dim,spacedim> >
77 class MappingFEField : public Mapping<dim,spacedim>
78 {
79 public:
112  MappingFEField (const DoFHandlerType &euler_dof_handler,
113  const VectorType &euler_vector,
114  const ComponentMask mask = ComponentMask());
115 
120 
125  virtual
126  Mapping<dim,spacedim> *clone () const;
127 
128 
129 
133  virtual
134  bool preserves_vertex_locations () const;
135 
136 
142  // for documentation, see the Mapping base class
143  virtual
146  const Point<dim> &p) const;
147 
148  // for documentation, see the Mapping base class
149  virtual
150  Point<dim>
152  const Point<spacedim> &p) const;
153 
163  // for documentation, see the Mapping base class
164  virtual
165  void
166  transform (const ArrayView<const Tensor<1,dim> > &input,
167  const MappingType type,
168  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
169  const ArrayView<Tensor<1,spacedim> > &output) const;
170 
171  // for documentation, see the Mapping base class
172  virtual
173  void
175  const MappingType type,
176  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
177  const ArrayView<Tensor<2,spacedim> > &output) const;
178 
179  // for documentation, see the Mapping base class
180  virtual
181  void
182  transform (const ArrayView<const Tensor<2, dim> > &input,
183  const MappingType type,
184  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
185  const ArrayView<Tensor<2,spacedim> > &output) const;
186 
187  // for documentation, see the Mapping base class
188  virtual
189  void
191  const MappingType type,
192  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
193  const ArrayView<Tensor<3,spacedim> > &output) const;
194 
195  // for documentation, see the Mapping base class
196  virtual
197  void
198  transform (const ArrayView<const Tensor<3, dim> > &input,
199  const MappingType type,
200  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
201  const ArrayView<Tensor<3,spacedim> > &output) const;
202 
212  unsigned int get_degree () const;
213 
219 
224 
225 private:
226 
232  // documentation can be found in Mapping::requires_update_flags()
233  virtual
235  requires_update_flags (const UpdateFlags update_flags) const;
236 
237 public:
249  class InternalData : public Mapping<dim,spacedim>::InternalDataBase
250  {
251  public:
256  const ComponentMask mask);
257 
262  const double &shape (const unsigned int qpoint,
263  const unsigned int shape_nr) const;
264 
268  double &shape (const unsigned int qpoint,
269  const unsigned int shape_nr);
270 
274  const Tensor<1,dim> &derivative (const unsigned int qpoint,
275  const unsigned int shape_nr) const;
276 
280  Tensor<1,dim> &derivative (const unsigned int qpoint,
281  const unsigned int shape_nr);
282 
286  const Tensor<2,dim> &second_derivative (const unsigned int qpoint,
287  const unsigned int shape_nr) const;
288 
292  Tensor<2,dim> &second_derivative (const unsigned int qpoint,
293  const unsigned int shape_nr);
294 
298  const Tensor<3,dim> &third_derivative (const unsigned int qpoint,
299  const unsigned int shape_nr) const;
300 
304  Tensor<3,dim> &third_derivative (const unsigned int qpoint,
305  const unsigned int shape_nr);
306 
310  const Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
311  const unsigned int shape_nr) const;
312 
316  Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
317  const unsigned int shape_nr);
318 
322  virtual std::size_t memory_consumption () const;
323 
329  std::vector<double> shape_values;
330 
336  std::vector<Tensor<1,dim> > shape_derivatives;
337 
344  std::vector<Tensor<2,dim> > shape_second_derivatives;
345 
352  std::vector<Tensor<3,dim> > shape_third_derivatives;
353 
360  std::vector<Tensor<4,dim> > shape_fourth_derivatives;
361 
375  std_cxx11::array<std::vector<Tensor<1,dim> >, GeometryInfo<dim>::faces_per_cell *(dim-1)> unit_tangentials;
376 
383  unsigned int n_shape_functions;
384 
398 
408  mutable std::vector<DerivativeForm<1,dim, spacedim > > covariant;
409 
417  mutable std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
418 
423  mutable std::vector<double> volume_elements;
424 
428  mutable std::vector<std::vector<Tensor<1,spacedim> > > aux;
429 
433  mutable std::vector<types::global_dof_index> local_dof_indices;
434 
438  mutable std::vector<double> local_dof_values;
439  };
440 
441 private:
442 
443  // documentation can be found in Mapping::get_data()
444  virtual
445  InternalData *
446  get_data (const UpdateFlags,
447  const Quadrature<dim> &quadrature) const;
448 
449  // documentation can be found in Mapping::get_face_data()
450  virtual
452  get_face_data (const UpdateFlags flags,
453  const Quadrature<dim-1>& quadrature) const;
454 
455  // documentation can be found in Mapping::get_subface_data()
456  virtual
458  get_subface_data (const UpdateFlags flags,
459  const Quadrature<dim-1>& quadrature) const;
460 
461  // documentation can be found in Mapping::fill_fe_values()
462  virtual
464  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
465  const CellSimilarity::Similarity cell_similarity,
466  const Quadrature<dim> &quadrature,
467  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
469 
470  // documentation can be found in Mapping::fill_fe_face_values()
471  virtual void
472  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
473  const unsigned int face_no,
474  const Quadrature<dim-1> &quadrature,
475  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
477 
478  // documentation can be found in Mapping::fill_fe_subface_values()
479  virtual void
480  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
481  const unsigned int face_no,
482  const unsigned int subface_no,
483  const Quadrature<dim-1> &quadrature,
484  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
486 
496 
505 
506 
511 
512 private:
529  do_transform_unit_to_real_cell (const InternalData &mdata) const;
530 
531 
546  Point<dim>
548  const Point<spacedim> &p,
549  const Point<dim> &initial_p_unit,
550  InternalData &mdata) const;
551 
557 
561  virtual void
562  compute_shapes_virtual (const std::vector<Point<dim> > &unit_points,
564 
565  /*
566  * Which components to use for the mapping.
567  */
568  const ComponentMask fe_mask;
569 
579  std::vector<unsigned int> fe_to_real;
580 
581  void
582  compute_data (const UpdateFlags update_flags,
583  const Quadrature<dim> &q,
584  const unsigned int n_original_q_points,
585  InternalData &data) const;
586 
587  void
588  compute_face_data (const UpdateFlags update_flags,
589  const Quadrature<dim> &q,
590  const unsigned int n_original_q_points,
591  InternalData &data) const;
592 
593 
597  template <int,int,class,class> friend class MappingFEField;
598 };
599 
602 /* -------------- declaration of explicit specializations ------------- */
603 
604 #ifndef DOXYGEN
605 
606 
607 
608 
609 #endif // DOXYGEN
610 
611 DEAL_II_NAMESPACE_CLOSE
612 
613 #endif
virtual Mapping< dim, spacedim >::InternalDataBase * get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
virtual std::size_t memory_consumption() const
MappingType
Definition: mapping.h:50
std::vector< unsigned int > fe_to_real
SmartPointer< const FiniteElement< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > fe
Point< spacedim > do_transform_unit_to_real_cell(const InternalData &mdata) const
unsigned int get_degree() const
virtual Mapping< dim, spacedim > * clone() const
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const
InternalData(const FiniteElement< dim, spacedim > &fe, const ComponentMask mask)
std::vector< Tensor< 3, dim > > shape_third_derivatives
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const
#define DeclException0(Exception0)
Definition: exceptions.h:541
std::vector< double > volume_elements
std_cxx11::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim-1)> unit_tangentials
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
virtual InternalData * get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const
virtual Mapping< dim, spacedim >::InternalDataBase * get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
SmartPointer< const DoFHandlerType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_dof_handler
virtual bool preserves_vertex_locations() const
std::vector< double > local_dof_values
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< double > shape_values
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
friend class MappingFEField
std::vector< types::global_dof_index > local_dof_indices
virtual void compute_shapes_virtual(const std::vector< Point< dim > > &unit_points, typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
SmartPointer< const VectorType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_vector
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
ComponentMask get_component_mask() const
std::vector< Tensor< 2, dim > > shape_second_derivatives
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
void update_internal_dofs(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
static ::ExceptionBase & ExcInactiveCell()
std::vector< Tensor< 1, dim > > shape_derivatives
Point< dim > do_transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit, InternalData &mdata) const
std::vector< std::vector< Tensor< 1, spacedim > > > aux