Reference documentation for deal.II version 9.5.0
\(\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\}}\)
Loading...
Searching...
No Matches
mapping_manifold.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2016 - 2022 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
54template <int dim, int spacedim = dim>
55class MappingManifold : public Mapping<dim, spacedim>
56{
57public:
61 MappingManifold() = default;
62
67
68 // for documentation, see the Mapping base class
69 virtual std::unique_ptr<Mapping<dim, spacedim>>
70 clone() const override;
71
76 virtual bool
78
79 virtual bool
80 is_compatible_with(const ReferenceCell &cell_type) 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,
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,
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,
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,
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,
141 const ArrayView<Tensor<3, spacedim>> &output) const override;
142
163 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
164 {
165 public:
169 InternalData() = default;
170
179 void
180 initialize(const UpdateFlags update_flags,
181 const Quadrature<dim> &quadrature,
182 const unsigned int n_original_q_points);
183
189 void
190 initialize_face(const UpdateFlags update_flags,
191 const Quadrature<dim> &quadrature,
192 const unsigned int n_original_q_points);
193
199 void
201
205 void
208
212 virtual std::size_t
213 memory_consumption() const override;
214
220 mutable std::vector<Point<spacedim>> vertices;
221
228
235
254 std::vector<std::vector<double>> cell_manifold_quadrature_weights;
255
270 mutable std::vector<double> vertex_weights;
271
285 std::array<std::vector<Tensor<1, dim>>,
288
298 mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
299
307 mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
308
312 mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
313
318 mutable std::vector<double> volume_elements;
319
326 };
327
328private:
329 // documentation can be found in Mapping::requires_update_flags()
330 virtual UpdateFlags
331 requires_update_flags(const UpdateFlags update_flags) const override;
332
333 // documentation can be found in Mapping::get_data()
334 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
335 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
336
337 using Mapping<dim, spacedim>::get_face_data;
338
339 // documentation can be found in Mapping::get_face_data()
340 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
341 get_face_data(const UpdateFlags flags,
342 const hp::QCollection<dim - 1> &quadrature) const override;
343
344 // documentation can be found in Mapping::get_subface_data()
345 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
346 get_subface_data(const UpdateFlags flags,
347 const Quadrature<dim - 1> &quadrature) const override;
348
349 // documentation can be found in Mapping::fill_fe_values()
353 const CellSimilarity::Similarity cell_similarity,
354 const Quadrature<dim> & quadrature,
355 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
357 &output_data) const override;
358
359 using Mapping<dim, spacedim>::fill_fe_face_values;
360
361 // documentation can be found in Mapping::fill_fe_face_values()
362 virtual void
365 const unsigned int face_no,
366 const hp::QCollection<dim - 1> & quadrature,
367 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
369 &output_data) const override;
370
371 // documentation can be found in Mapping::fill_fe_subface_values()
372 virtual void
375 const unsigned int face_no,
376 const unsigned int subface_no,
377 const Quadrature<dim - 1> & quadrature,
378 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
380 &output_data) const override;
381
385};
386
387
388
391/*----------------------------------------------------------------------*/
392
393#ifndef DOXYGEN
394
395template <int dim, int spacedim>
396inline void
398 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
399{
401 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
402 vertices[i] = cell->vertex(i);
403 this->cell = cell;
404}
405
406
407template <int dim, int spacedim>
408inline void
411{
412 cell_manifold_quadrature_weights.resize(
413 quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
414 for (unsigned int q = 0; q < quad.size(); ++q)
415 {
416 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
417 {
418 cell_manifold_quadrature_weights[q][i] =
420 }
421 }
422}
423
424
425
426template <int dim, int spacedim>
427inline bool
429{
430 return true;
431}
432
433
434template <int dim, int spacedim>
435bool
437 const ReferenceCell &cell_type) const
438{
439 if (cell_type.get_dimension() != dim)
440 return false; // TODO: or is this an error?
441
442 if (cell_type.is_hyper_cube())
443 return true;
444
445 return false;
446}
447
448
449
450#endif // DOXYGEN
451
452/* -------------- declaration of explicit specializations ------------- */
453
454
456
457#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 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 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 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
Abstract base class for mapping classes.
Definition mapping.h:317
Definition point.h:112
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:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
unsigned int vertex_indices[2]
UpdateFlags
MappingKind
Definition mapping.h:78
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)