Reference documentation for deal.II version Git a7d922cd5a 2021-01-15 21:19:32 -0700
\(\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 - 2020 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 
87  // for documentation, see the Mapping base class
88  virtual Point<spacedim>
91  const Point<dim> &p) const override;
92 
93  // for documentation, see the Mapping base class
94  virtual Point<dim>
97  const Point<spacedim> &p) const override;
98 
108  // for documentation, see the Mapping base class
109  virtual void
110  transform(const ArrayView<const Tensor<1, dim>> & input,
111  const MappingKind kind,
112  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
113  const ArrayView<Tensor<1, spacedim>> &output) const override;
114 
115  // for documentation, see the Mapping base class
116  virtual void
118  const MappingKind kind,
119  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
120  const ArrayView<Tensor<2, spacedim>> &output) const override;
121 
122  // for documentation, see the Mapping base class
123  virtual void
124  transform(const ArrayView<const Tensor<2, dim>> & input,
125  const MappingKind kind,
126  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
127  const ArrayView<Tensor<2, spacedim>> &output) const override;
128 
129  // for documentation, see the Mapping base class
130  virtual void
132  const MappingKind kind,
133  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
134  const ArrayView<Tensor<3, spacedim>> &output) const override;
135 
136  // for documentation, see the Mapping base class
137  virtual void
138  transform(const ArrayView<const Tensor<3, dim>> & input,
139  const MappingKind kind,
140  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
141  const ArrayView<Tensor<3, spacedim>> &output) const override;
142 
152 public:
164  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
165  {
166  public:
170  InternalData() = default;
171 
180  void
181  initialize(const UpdateFlags update_flags,
182  const Quadrature<dim> &quadrature,
183  const unsigned int n_original_q_points);
184 
190  void
191  initialize_face(const UpdateFlags update_flags,
192  const Quadrature<dim> &quadrature,
193  const unsigned int n_original_q_points);
194 
195 
201  void
203 
207  void
209  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
210 
214  virtual std::size_t
215  memory_consumption() const override;
216 
222  mutable std::vector<Point<spacedim>> vertices;
223 
230 
237 
238 
257  std::vector<std::vector<double>> cell_manifold_quadrature_weights;
258 
259 
274  mutable std::vector<double> vertex_weights;
275 
289  std::array<std::vector<Tensor<1, dim>>,
292 
302  mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
303 
311  mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
312 
316  mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
317 
322  mutable std::vector<double> volume_elements;
323 
330  };
331 
332 
333  // documentation can be found in Mapping::requires_update_flags()
334  virtual UpdateFlags
335  requires_update_flags(const UpdateFlags update_flags) const override;
336 
337  // documentation can be found in Mapping::get_data()
338  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
339  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
340 
342 
343  // documentation can be found in Mapping::get_face_data()
344  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
345  get_face_data(const UpdateFlags flags,
346  const hp::QCollection<dim - 1> &quadrature) const override;
347 
348  // documentation can be found in Mapping::get_subface_data()
349  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
350  get_subface_data(const UpdateFlags flags,
351  const Quadrature<dim - 1> &quadrature) const override;
352 
353  // documentation can be found in Mapping::fill_fe_values()
356  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
357  const CellSimilarity::Similarity cell_similarity,
358  const Quadrature<dim> & quadrature,
359  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
361  &output_data) const override;
362 
364 
365  // documentation can be found in Mapping::fill_fe_face_values()
366  virtual void
368  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
369  const unsigned int face_no,
370  const hp::QCollection<dim - 1> & quadrature,
371  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
373  &output_data) const override;
374 
375  // documentation can be found in Mapping::fill_fe_subface_values()
376  virtual void
378  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
379  const unsigned int face_no,
380  const unsigned int subface_no,
381  const Quadrature<dim - 1> & quadrature,
382  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
384  &output_data) const override;
385 
389 };
390 
391 
392 
395 /*----------------------------------------------------------------------*/
396 
397 #ifndef DOXYGEN
398 
399 template <int dim, int spacedim>
400 inline void
403 {
405  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
406  vertices[i] = cell->vertex(i);
407  this->cell = cell;
408 }
409 
410 
411 template <int dim, int spacedim>
412 inline void
415 {
417  quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
418  for (unsigned int q = 0; q < quad.size(); ++q)
419  {
420  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
421  {
424  }
425  }
426 }
427 
428 
429 
430 template <int dim, int spacedim>
431 inline bool
433 {
434  return true;
435 }
436 
437 #endif // DOXYGEN
438 
439 /* -------------- declaration of explicit specializations ------------- */
440 
441 
443 
444 #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
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
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:380
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
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:379
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