Reference documentation for deal.II version 8.4.2
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/fe/mapping.h>
23 #include <deal.II/fe/fe.h>
24 #include <deal.II/base/qprojector.h>
25 #include <deal.II/base/thread_management.h>
26 
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 
33 
71 template <int dim, int spacedim=dim,
72  typename VectorType=Vector<double>,
73  typename DoFHandlerType=DoFHandler<dim,spacedim> >
74 class MappingFEField : public Mapping<dim,spacedim>
75 {
76 public:
109  MappingFEField (const DoFHandlerType &euler_dof_handler,
110  const VectorType &euler_vector,
111  const ComponentMask mask = ComponentMask());
112 
117 
122  virtual
123  Mapping<dim,spacedim> *clone () const;
124 
125 
126 
130  virtual
131  bool preserves_vertex_locations () const;
132 
133 
139  // for documentation, see the Mapping base class
140  virtual
143  const Point<dim> &p) const;
144 
145  // for documentation, see the Mapping base class
146  virtual
147  Point<dim>
149  const Point<spacedim> &p) const;
150 
160  // for documentation, see the Mapping base class
161  virtual
162  void
163  transform (const ArrayView<const Tensor<1,dim> > &input,
164  const MappingType type,
165  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
166  const ArrayView<Tensor<1,spacedim> > &output) const;
167 
168  // for documentation, see the Mapping base class
169  virtual
170  void
172  const MappingType type,
173  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
174  const ArrayView<Tensor<2,spacedim> > &output) const;
175 
176  // for documentation, see the Mapping base class
177  virtual
178  void
179  transform (const ArrayView<const Tensor<2, dim> > &input,
180  const MappingType type,
181  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
182  const ArrayView<Tensor<2,spacedim> > &output) const;
183 
184  // for documentation, see the Mapping base class
185  virtual
186  void
188  const MappingType type,
189  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
190  const ArrayView<Tensor<3,spacedim> > &output) const;
191 
192  // for documentation, see the Mapping base class
193  virtual
194  void
195  transform (const ArrayView<const Tensor<3, dim> > &input,
196  const MappingType type,
197  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
198  const ArrayView<Tensor<3,spacedim> > &output) const;
199 
209  unsigned int get_degree () const;
210 
216 
220  DeclException0(ExcInactiveCell);
221 
222 private:
223 
229  // documentation can be found in Mapping::requires_update_flags()
230  virtual
232  requires_update_flags (const UpdateFlags update_flags) const;
233 
234 public:
246  class InternalData : public Mapping<dim,spacedim>::InternalDataBase
247  {
248  public:
253  const ComponentMask mask);
254 
259  const double &shape (const unsigned int qpoint,
260  const unsigned int shape_nr) const;
261 
265  double &shape (const unsigned int qpoint,
266  const unsigned int shape_nr);
267 
271  const Tensor<1,dim> &derivative (const unsigned int qpoint,
272  const unsigned int shape_nr) const;
273 
277  Tensor<1,dim> &derivative (const unsigned int qpoint,
278  const unsigned int shape_nr);
279 
283  const Tensor<2,dim> &second_derivative (const unsigned int qpoint,
284  const unsigned int shape_nr) const;
285 
289  Tensor<2,dim> &second_derivative (const unsigned int qpoint,
290  const unsigned int shape_nr);
291 
295  const Tensor<3,dim> &third_derivative (const unsigned int qpoint,
296  const unsigned int shape_nr) const;
297 
301  Tensor<3,dim> &third_derivative (const unsigned int qpoint,
302  const unsigned int shape_nr);
303 
307  const Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
308  const unsigned int shape_nr) const;
309 
313  Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
314  const unsigned int shape_nr);
315 
319  virtual std::size_t memory_consumption () const;
320 
326  std::vector<double> shape_values;
327 
333  std::vector<Tensor<1,dim> > shape_derivatives;
334 
341  std::vector<Tensor<2,dim> > shape_second_derivatives;
342 
349  std::vector<Tensor<3,dim> > shape_third_derivatives;
350 
357  std::vector<Tensor<4,dim> > shape_fourth_derivatives;
358 
372  std::vector<std::vector<Tensor<1,dim> > > unit_tangentials;
373 
380  unsigned int n_shape_functions;
381 
395 
405  mutable std::vector<DerivativeForm<1,dim, spacedim > > covariant;
406 
414  mutable std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
415 
420  mutable std::vector<double> volume_elements;
421 
425  mutable std::vector<std::vector<Tensor<1,spacedim> > > aux;
426 
430  mutable std::vector<types::global_dof_index> local_dof_indices;
431 
435  mutable std::vector<double> local_dof_values;
436  };
437 
438 private:
439 
440  // documentation can be found in Mapping::get_data()
441  virtual
442  InternalData *
443  get_data (const UpdateFlags,
444  const Quadrature<dim> &quadrature) const;
445 
446  // documentation can be found in Mapping::get_face_data()
447  virtual
449  get_face_data (const UpdateFlags flags,
450  const Quadrature<dim-1>& quadrature) const;
451 
452  // documentation can be found in Mapping::get_subface_data()
453  virtual
455  get_subface_data (const UpdateFlags flags,
456  const Quadrature<dim-1>& quadrature) const;
457 
458  // documentation can be found in Mapping::fill_fe_values()
459  virtual
460  CellSimilarity::Similarity
461  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
462  const CellSimilarity::Similarity cell_similarity,
463  const Quadrature<dim> &quadrature,
464  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
466 
467  // documentation can be found in Mapping::fill_fe_face_values()
468  virtual void
469  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
470  const unsigned int face_no,
471  const Quadrature<dim-1> &quadrature,
472  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
474 
475  // documentation can be found in Mapping::fill_fe_subface_values()
476  virtual void
477  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
478  const unsigned int face_no,
479  const unsigned int subface_no,
480  const Quadrature<dim-1> &quadrature,
481  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
483 
493 
502 
503 
508 
509 private:
526  do_transform_unit_to_real_cell (const InternalData &mdata) const;
527 
528 
543  Point<dim>
545  const Point<spacedim> &p,
546  const Point<dim> &initial_p_unit,
547  InternalData &mdata) const;
548 
553  const typename MappingFEField<dim, spacedim>::InternalData &data) const;
554 
558  virtual void
559  compute_shapes_virtual (const std::vector<Point<dim> > &unit_points,
560  typename MappingFEField<dim, spacedim>::InternalData &data) const;
561 
562  /*
563  * Which components to use for the mapping.
564  */
565  const ComponentMask fe_mask;
566 
576  std::vector<unsigned int> fe_to_real;
577 
578  void
579  compute_data (const UpdateFlags update_flags,
580  const Quadrature<dim> &q,
581  const unsigned int n_original_q_points,
582  InternalData &data) const;
583 
584  void
585  compute_face_data (const UpdateFlags update_flags,
586  const Quadrature<dim> &q,
587  const unsigned int n_original_q_points,
588  InternalData &data) const;
589 
590 
594  template <int,int,class,class> friend class MappingFEField;
595 };
596 
599 /* -------------- declaration of explicit specializations ------------- */
600 
601 #ifndef DOXYGEN
602 
603 
604 
605 
606 #endif // DOXYGEN
607 
608 DEAL_II_NAMESPACE_CLOSE
609 
610 #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
std::vector< std::vector< Tensor< 1, dim > > > unit_tangentials
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:52
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const
std::vector< double > volume_elements
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 void compute_shapes_virtual(const std::vector< Point< dim > > &unit_points, typename MappingFEField< dim, spacedim >::InternalData &data) const
void update_internal_dofs(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim >::InternalData &data) 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
DeclException0(ExcInactiveCell)
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
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
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