Reference documentation for deal.II version Git 0b4b57890a 2021-10-16 12:51:50 +0200
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mapping_manifold.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2021 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_mapping_manifold_h
17 #define dealii_mapping_manifold_h
18 
19 
20 #include <deal.II/base/config.h>
21 
25 
26 #include <deal.II/fe/mapping.h>
27 
28 #include <cmath>
29 
31 
32 template <int, int>
33 class MappingQ;
34 
35 
38 
39 
57 template <int dim, int spacedim = dim>
58 class MappingManifold : public Mapping<dim, spacedim>
59 {
60 public:
64  MappingManifold() = default;
65 
70 
71  // for documentation, see the Mapping base class
72  virtual std::unique_ptr<Mapping<dim, spacedim>>
73  clone() const override;
74 
79  virtual bool
80  preserves_vertex_locations() const override;
81 
82  virtual bool
83  is_compatible_with(const ReferenceCell &cell_type) const override;
84 
90  // for documentation, see the Mapping base class
91  virtual Point<spacedim>
94  const Point<dim> &p) const override;
95 
96  // for documentation, see the Mapping base class
97  virtual Point<dim>
100  const Point<spacedim> &p) const override;
101 
111  // for documentation, see the Mapping base class
112  virtual void
113  transform(const ArrayView<const Tensor<1, dim>> & input,
114  const MappingKind kind,
115  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
116  const ArrayView<Tensor<1, spacedim>> &output) const override;
117 
118  // for documentation, see the Mapping base class
119  virtual void
121  const MappingKind kind,
122  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
123  const ArrayView<Tensor<2, spacedim>> &output) const override;
124 
125  // for documentation, see the Mapping base class
126  virtual void
127  transform(const ArrayView<const Tensor<2, dim>> & input,
128  const MappingKind kind,
129  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
130  const ArrayView<Tensor<2, spacedim>> &output) const override;
131 
132  // for documentation, see the Mapping base class
133  virtual void
135  const MappingKind kind,
136  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
137  const ArrayView<Tensor<3, spacedim>> &output) const override;
138 
139  // for documentation, see the Mapping base class
140  virtual void
141  transform(const ArrayView<const Tensor<3, dim>> & input,
142  const MappingKind kind,
143  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
144  const ArrayView<Tensor<3, spacedim>> &output) const override;
145 
166  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
167  {
168  public:
172  InternalData() = default;
173 
182  void
183  initialize(const UpdateFlags update_flags,
184  const Quadrature<dim> &quadrature,
185  const unsigned int n_original_q_points);
186 
192  void
193  initialize_face(const UpdateFlags update_flags,
194  const Quadrature<dim> &quadrature,
195  const unsigned int n_original_q_points);
196 
202  void
204 
208  void
210  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
211 
215  virtual std::size_t
216  memory_consumption() const override;
217 
223  mutable std::vector<Point<spacedim>> vertices;
224 
231 
238 
257  std::vector<std::vector<double>> cell_manifold_quadrature_weights;
258 
273  mutable std::vector<double> vertex_weights;
274 
288  std::array<std::vector<Tensor<1, dim>>,
291 
301  mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
302 
310  mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
311 
315  mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
316 
321  mutable std::vector<double> volume_elements;
322 
329  };
330 
331 private:
332  // documentation can be found in Mapping::requires_update_flags()
333  virtual UpdateFlags
334  requires_update_flags(const UpdateFlags update_flags) const override;
335 
336  // documentation can be found in Mapping::get_data()
337  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
338  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
339 
341 
342  // documentation can be found in Mapping::get_face_data()
343  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
344  get_face_data(const UpdateFlags flags,
345  const hp::QCollection<dim - 1> &quadrature) const override;
346 
347  // documentation can be found in Mapping::get_subface_data()
348  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
349  get_subface_data(const UpdateFlags flags,
350  const Quadrature<dim - 1> &quadrature) const override;
351 
352  // documentation can be found in Mapping::fill_fe_values()
355  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
356  const CellSimilarity::Similarity cell_similarity,
357  const Quadrature<dim> & quadrature,
358  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
360  &output_data) const override;
361 
363 
364  // documentation can be found in Mapping::fill_fe_face_values()
365  virtual void
367  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
368  const unsigned int face_no,
369  const hp::QCollection<dim - 1> & quadrature,
370  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
372  &output_data) const override;
373 
374  // documentation can be found in Mapping::fill_fe_subface_values()
375  virtual void
377  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
378  const unsigned int face_no,
379  const unsigned int subface_no,
380  const Quadrature<dim - 1> & quadrature,
381  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
383  &output_data) const override;
384 
388 };
389 
390 
391 
394 /*----------------------------------------------------------------------*/
395 
396 #ifndef DOXYGEN
397 
398 template <int dim, int spacedim>
399 inline void
402 {
404  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
405  vertices[i] = cell->vertex(i);
406  this->cell = cell;
407 }
408 
409 
410 template <int dim, int spacedim>
411 inline void
414 {
416  quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
417  for (unsigned int q = 0; q < quad.size(); ++q)
418  {
419  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
420  {
423  }
424  }
425 }
426 
427 
428 
429 template <int dim, int spacedim>
430 inline bool
432 {
433  return true;
434 }
435 
436 
437 template <int dim, int spacedim>
438 bool
440  const ReferenceCell &cell_type) const
441 {
442  if (cell_type.get_dimension() != dim)
443  return false; // TODO: or is this an error?
444 
445  if (cell_type.is_hyper_cube())
446  return true;
447 
448  return false;
449 }
450 
451 
452 
453 #endif // DOXYGEN
454 
455 /* -------------- declaration of explicit specializations ------------- */
456 
457 
459 
460 #endif
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
virtual bool preserves_vertex_locations() const override
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
SmartPointer< const Manifold< dim, spacedim > > manifold
bool is_hyper_cube() const
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
const Point< dim > & point(const unsigned int i) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
void store_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual std::size_t memory_consumption() const override
virtual bool is_compatible_with(const ReferenceCell &cell_type) const override
MappingKind
Definition: mapping.h:64
Triangulation< dim, spacedim >::cell_iterator cell
std::vector< Point< spacedim > > vertices
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
std::vector< std::vector< Tensor< 1, spacedim > > > aux
UpdateFlags
Abstract base class for mapping classes.
Definition: mapping.h:303
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:401
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< double > volume_elements
unsigned int size() const
std::vector< double > vertex_weights
unsigned int get_dimension() const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
void compute_manifold_quadrature_weights(const Quadrature< dim > &quadrature)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:400
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< std::vector< double > > cell_manifold_quadrature_weights
virtual void transform(const ArrayView< const Tensor< 1, dim >> &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim >> &output) const override
MappingManifold()=default