Reference documentation for deal.II version 9.3.3
\(\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
32template <int, int>
33class MappingQ;
34
35
38
39
57template <int dim, int spacedim = dim>
58class MappingManifold : public Mapping<dim, spacedim>
59{
60public:
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
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,
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,
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,
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,
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,
144 const ArrayView<Tensor<3, spacedim>> &output) const override;
145
155public:
167 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
168 {
169 public:
173 InternalData() = default;
174
183 void
184 initialize(const UpdateFlags update_flags,
185 const Quadrature<dim> &quadrature,
186 const unsigned int n_original_q_points);
187
193 void
194 initialize_face(const UpdateFlags update_flags,
195 const Quadrature<dim> &quadrature,
196 const unsigned int n_original_q_points);
197
198
204 void
206
210 void
213
217 virtual std::size_t
218 memory_consumption() const override;
219
225 mutable std::vector<Point<spacedim>> vertices;
226
233
240
241
260 std::vector<std::vector<double>> cell_manifold_quadrature_weights;
261
262
277 mutable std::vector<double> vertex_weights;
278
292 std::array<std::vector<Tensor<1, dim>>,
295
305 mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
306
314 mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
315
319 mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
320
325 mutable std::vector<double> volume_elements;
326
333 };
334
335
336 // documentation can be found in Mapping::requires_update_flags()
337 virtual UpdateFlags
338 requires_update_flags(const UpdateFlags update_flags) const override;
339
340 // documentation can be found in Mapping::get_data()
341 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
342 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
343
344 using Mapping<dim, spacedim>::get_face_data;
345
346 // documentation can be found in Mapping::get_face_data()
347 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
348 get_face_data(const UpdateFlags flags,
349 const hp::QCollection<dim - 1> &quadrature) const override;
350
351 // documentation can be found in Mapping::get_subface_data()
352 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
353 get_subface_data(const UpdateFlags flags,
354 const Quadrature<dim - 1> &quadrature) const override;
355
356 // documentation can be found in Mapping::fill_fe_values()
360 const CellSimilarity::Similarity cell_similarity,
361 const Quadrature<dim> & quadrature,
362 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
364 &output_data) const override;
365
366 using Mapping<dim, spacedim>::fill_fe_face_values;
367
368 // documentation can be found in Mapping::fill_fe_face_values()
369 virtual void
372 const unsigned int face_no,
373 const hp::QCollection<dim - 1> & quadrature,
374 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
376 &output_data) const override;
377
378 // documentation can be found in Mapping::fill_fe_subface_values()
379 virtual void
382 const unsigned int face_no,
383 const unsigned int subface_no,
384 const Quadrature<dim - 1> & quadrature,
385 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
387 &output_data) const override;
388
392};
393
394
395
398/*----------------------------------------------------------------------*/
399
400#ifndef DOXYGEN
401
402template <int dim, int spacedim>
403inline void
405 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
406{
408 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
409 vertices[i] = cell->vertex(i);
410 this->cell = cell;
411}
412
413
414template <int dim, int spacedim>
415inline void
418{
419 cell_manifold_quadrature_weights.resize(
420 quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
421 for (unsigned int q = 0; q < quad.size(); ++q)
422 {
423 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
424 {
425 cell_manifold_quadrature_weights[q][i] =
427 }
428 }
429}
430
431
432
433template <int dim, int spacedim>
434inline bool
436{
437 return true;
438}
439
440
441template <int dim, int spacedim>
442bool
444 const ReferenceCell &cell_type) const
445{
446 if (cell_type.get_dimension() != dim)
447 return false; // TODO: or is this an error?
448
449 if (cell_type.is_hyper_cube())
450 return true;
451
452 return false;
453}
454
455
456
457#endif // DOXYGEN
458
459/* -------------- declaration of explicit specializations ------------- */
460
461
463
464#endif
void store_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
std::vector< double > volume_elements
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
SmartPointer< const Manifold< dim, spacedim > > manifold
std::vector< std::vector< Tensor< 1, spacedim > > > aux
virtual std::size_t memory_consumption() const override
std::vector< std::vector< double > > cell_manifold_quadrature_weights
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< Point< spacedim > > vertices
std::vector< double > vertex_weights
void compute_manifold_quadrature_weights(const Quadrature< dim > &quadrature)
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
Triangulation< dim, spacedim >::cell_iterator cell
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
MappingManifold()=default
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
virtual bool is_compatible_with(const ReferenceCell &cell_type) const override
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
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 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 bool preserves_vertex_locations() const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
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
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) 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
Abstract base class for mapping classes.
Definition: mapping.h:304
const Point< dim > & point(const unsigned int i) const
unsigned int size() const
bool is_hyper_cube() const
unsigned int get_dimension() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
UpdateFlags
MappingKind
Definition: mapping.h:65
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)