Reference documentation for deal.II version 9.0.0
mapping_fe_field.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_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 #include <deal.II/fe/fe_values.h>
28 
29 #include <array>
30 
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 
37 
75 template <int dim, int spacedim=dim,
76  typename VectorType=Vector<double>,
77  typename DoFHandlerType=DoFHandler<dim,spacedim> >
78 class MappingFEField : public Mapping<dim,spacedim>
79 {
80 public:
113  MappingFEField (const DoFHandlerType &euler_dof_handler,
114  const VectorType &euler_vector,
115  const ComponentMask &mask = ComponentMask());
116 
121 
126  virtual
127  std::unique_ptr<Mapping<dim,spacedim>> clone () const override;
128 
134  virtual
135  bool preserves_vertex_locations () const override;
136 
144  virtual
145  std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
146  get_vertices (const typename Triangulation<dim,spacedim>::cell_iterator &cell) const override;
147 
153  // for documentation, see the Mapping base class
154  virtual
157  const Point<dim> &p) const override;
158 
159  // for documentation, see the Mapping base class
160  virtual
161  Point<dim>
163  const Point<spacedim> &p) const override;
164 
174  // for documentation, see the Mapping base class
175  virtual
176  void
177  transform (const ArrayView<const Tensor<1,dim> > &input,
178  const MappingType type,
179  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
180  const ArrayView<Tensor<1,spacedim> > &output) const override;
181 
182  // for documentation, see the Mapping base class
183  virtual
184  void
186  const MappingType type,
187  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
188  const ArrayView<Tensor<2,spacedim> > &output) const override;
189 
190  // for documentation, see the Mapping base class
191  virtual
192  void
193  transform (const ArrayView<const Tensor<2, dim> > &input,
194  const MappingType type,
195  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
196  const ArrayView<Tensor<2,spacedim> > &output) const override;
197 
198  // for documentation, see the Mapping base class
199  virtual
200  void
202  const MappingType type,
203  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
204  const ArrayView<Tensor<3,spacedim> > &output) const override;
205 
206  // for documentation, see the Mapping base class
207  virtual
208  void
209  transform (const ArrayView<const Tensor<3, dim> > &input,
210  const MappingType type,
211  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
212  const ArrayView<Tensor<3,spacedim> > &output) const override;
213 
222  unsigned int get_degree () const;
223 
229 
234 
235 private:
236 
242  // documentation can be found in Mapping::requires_update_flags()
243  virtual
245  requires_update_flags (const UpdateFlags update_flags) const override;
246 
247 public:
259  class InternalData : public Mapping<dim,spacedim>::InternalDataBase
260  {
261  public:
266  const ComponentMask &mask);
267 
272  const double &shape (const unsigned int qpoint,
273  const unsigned int shape_nr) const;
274 
278  double &shape (const unsigned int qpoint,
279  const unsigned int shape_nr);
280 
284  const Tensor<1,dim> &derivative (const unsigned int qpoint,
285  const unsigned int shape_nr) const;
286 
290  Tensor<1,dim> &derivative (const unsigned int qpoint,
291  const unsigned int shape_nr);
292 
296  const Tensor<2,dim> &second_derivative (const unsigned int qpoint,
297  const unsigned int shape_nr) const;
298 
302  Tensor<2,dim> &second_derivative (const unsigned int qpoint,
303  const unsigned int shape_nr);
304 
308  const Tensor<3,dim> &third_derivative (const unsigned int qpoint,
309  const unsigned int shape_nr) const;
310 
314  Tensor<3,dim> &third_derivative (const unsigned int qpoint,
315  const unsigned int shape_nr);
316 
320  const Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
321  const unsigned int shape_nr) const;
322 
326  Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
327  const unsigned int shape_nr);
328 
332  virtual
333  std::size_t memory_consumption () const override;
334 
340  std::vector<double> shape_values;
341 
347  std::vector<Tensor<1,dim> > shape_derivatives;
348 
355  std::vector<Tensor<2,dim> > shape_second_derivatives;
356 
363  std::vector<Tensor<3,dim> > shape_third_derivatives;
364 
371  std::vector<Tensor<4,dim> > shape_fourth_derivatives;
372 
386  std::array<std::vector<Tensor<1,dim> >, GeometryInfo<dim>::faces_per_cell *(dim-1)> unit_tangentials;
387 
394  unsigned int n_shape_functions;
395 
409 
419  mutable std::vector<DerivativeForm<1,dim, spacedim > > covariant;
420 
428  mutable std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
429 
434  mutable std::vector<double> volume_elements;
435 
439  mutable std::vector<std::vector<Tensor<1,spacedim> > > aux;
440 
444  mutable std::vector<types::global_dof_index> local_dof_indices;
445 
449  mutable std::vector<double> local_dof_values;
450  };
451 
452 private:
453 
454  // documentation can be found in Mapping::get_data()
455  virtual
456  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
457  get_data (const UpdateFlags,
458  const Quadrature<dim> &quadrature) const override;
459 
460  // documentation can be found in Mapping::get_face_data()
461  virtual
462  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
463  get_face_data (const UpdateFlags flags,
464  const Quadrature<dim-1>& quadrature) const override;
465 
466  // documentation can be found in Mapping::get_subface_data()
467  virtual
468  std::unique_ptr<typename Mapping<dim,spacedim>::InternalDataBase>
469  get_subface_data (const UpdateFlags flags,
470  const Quadrature<dim-1>& quadrature) const override;
471 
472  // documentation can be found in Mapping::fill_fe_values()
473  virtual
475  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
476  const CellSimilarity::Similarity cell_similarity,
477  const Quadrature<dim> &quadrature,
478  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
480 
481  // documentation can be found in Mapping::fill_fe_face_values()
482  virtual
483  void
484  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
485  const unsigned int face_no,
486  const Quadrature<dim-1> &quadrature,
487  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
489 
490  // documentation can be found in Mapping::fill_fe_subface_values()
491  virtual
492  void
493  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
494  const unsigned int face_no,
495  const unsigned int subface_no,
496  const Quadrature<dim-1> &quadrature,
497  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
499 
509 
514 
515 private:
532  do_transform_unit_to_real_cell (const InternalData &mdata) const;
533 
534 
549  Point<dim>
551  const Point<spacedim> &p,
552  const Point<dim> &initial_p_unit,
553  InternalData &mdata) const;
554 
560 
564  virtual
565  void
566  compute_shapes_virtual (const std::vector<Point<dim> > &unit_points,
568 
569  /*
570  * Which components to use for the mapping.
571  */
572  const ComponentMask fe_mask;
573 
583  std::vector<unsigned int> fe_to_real;
584 
590 
595 
596  void
597  compute_data (const UpdateFlags update_flags,
598  const Quadrature<dim> &q,
599  const unsigned int n_original_q_points,
600  InternalData &data) const;
601 
602  void
603  compute_face_data (const UpdateFlags update_flags,
604  const Quadrature<dim> &q,
605  const unsigned int n_original_q_points,
606  InternalData &data) const;
607 
608 
612  template <int,int,class,class> friend class MappingFEField;
613 };
614 
617 /* -------------- declaration of explicit specializations ------------- */
618 
619 #ifndef DOXYGEN
620 
621 #endif // DOXYGEN
622 
623 DEAL_II_NAMESPACE_CLOSE
624 
625 #endif
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
MappingType
Definition: mapping.h:51
std::vector< unsigned int > fe_to_real
Point< spacedim > do_transform_unit_to_real_cell(const InternalData &mdata) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
unsigned int get_degree() const
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
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
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< Tensor< 3, dim > > shape_third_derivatives
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim-1)> unit_tangentials
UpdateFlags
Threads::Mutex fe_values_mutex
Abstract base class for mapping classes.
Definition: dof_tools.h:46
#define DeclException0(Exception0)
Definition: exceptions.h:323
std::vector< double > volume_elements
SmartPointer< const DoFHandlerType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_dof_handler
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
std::vector< double > local_dof_values
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
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
FEValues< dim, spacedim > fe_values
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
Definition: fe.h:33
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
ComponentMask get_component_mask() const
std::vector< Tensor< 2, dim > > shape_second_derivatives
void update_internal_dofs(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
virtual bool preserves_vertex_locations() const override
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
virtual std::size_t memory_consumption() const override
InternalData(const FiniteElement< dim, spacedim > &fe, const ComponentMask &mask)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
std::vector< std::vector< Tensor< 1, spacedim > > > aux