Reference documentation for deal.II version GIT relicensing-233-g802318d791 2024-03-28 20:20:02+00:00
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_mapping_manifold_h
16#define dealii_mapping_manifold_h
17
18
19#include <deal.II/base/config.h>
20
24
25#include <deal.II/fe/mapping.h>
26
27#include <cmath>
28
30
53template <int dim, int spacedim = dim>
54class MappingManifold : public Mapping<dim, spacedim>
55{
56public:
60 MappingManifold() = default;
61
66
67 // for documentation, see the Mapping base class
68 virtual std::unique_ptr<Mapping<dim, spacedim>>
69 clone() const override;
70
75 virtual bool
77
78 virtual bool
79 is_compatible_with(const ReferenceCell &cell_type) const override;
80
86 // for documentation, see the Mapping base class
87 virtual Point<spacedim>
90 const Point<dim> &p) const override;
91
92 // for documentation, see the Mapping base class
93 virtual Point<dim>
96 const Point<spacedim> &p) const override;
97
107 // for documentation, see the Mapping base class
108 virtual void
109 transform(const ArrayView<const Tensor<1, dim>> &input,
110 const MappingKind kind,
112 const ArrayView<Tensor<1, spacedim>> &output) const override;
113
114 // for documentation, see the Mapping base class
115 virtual void
117 const MappingKind kind,
119 const ArrayView<Tensor<2, spacedim>> &output) const override;
120
121 // for documentation, see the Mapping base class
122 virtual void
123 transform(const ArrayView<const Tensor<2, dim>> &input,
124 const MappingKind kind,
126 const ArrayView<Tensor<2, spacedim>> &output) const override;
127
128 // for documentation, see the Mapping base class
129 virtual void
131 const MappingKind kind,
133 const ArrayView<Tensor<3, spacedim>> &output) const override;
134
135 // for documentation, see the Mapping base class
136 virtual void
137 transform(const ArrayView<const Tensor<3, dim>> &input,
138 const MappingKind kind,
140 const ArrayView<Tensor<3, spacedim>> &output) const override;
141
162 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
163 {
164 public:
168 InternalData() = default;
169
178 void
179 initialize(const UpdateFlags update_flags,
180 const Quadrature<dim> &quadrature,
181 const unsigned int n_original_q_points);
182
188 void
189 initialize_face(const UpdateFlags update_flags,
190 const Quadrature<dim> &quadrature,
191 const unsigned int n_original_q_points);
192
198 void
200
204 void
207
211 virtual std::size_t
212 memory_consumption() const override;
213
219 mutable std::vector<Point<spacedim>> vertices;
220
227
234
253 std::vector<std::vector<double>> cell_manifold_quadrature_weights;
254
269 mutable std::vector<double> vertex_weights;
270
284 std::array<std::vector<Tensor<1, dim>>,
287
297 mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
298
306 mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
307
311 mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
312
317 mutable std::vector<double> volume_elements;
318
325 };
326
327private:
328 // documentation can be found in Mapping::requires_update_flags()
329 virtual UpdateFlags
330 requires_update_flags(const UpdateFlags update_flags) const override;
331
332 // documentation can be found in Mapping::get_data()
333 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
334 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
335
336 using Mapping<dim, spacedim>::get_face_data;
337
338 // documentation can be found in Mapping::get_face_data()
339 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
340 get_face_data(const UpdateFlags flags,
341 const hp::QCollection<dim - 1> &quadrature) const override;
342
343 // documentation can be found in Mapping::get_subface_data()
344 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
345 get_subface_data(const UpdateFlags flags,
346 const Quadrature<dim - 1> &quadrature) const override;
347
348 // documentation can be found in Mapping::fill_fe_values()
352 const CellSimilarity::Similarity cell_similarity,
353 const Quadrature<dim> &quadrature,
354 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
356 &output_data) const override;
357
358 using Mapping<dim, spacedim>::fill_fe_face_values;
359
360 // documentation can be found in Mapping::fill_fe_face_values()
361 virtual void
364 const unsigned int face_no,
365 const hp::QCollection<dim - 1> &quadrature,
366 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
368 &output_data) const override;
369
370 // documentation can be found in Mapping::fill_fe_subface_values()
371 virtual void
374 const unsigned int face_no,
375 const unsigned int subface_no,
376 const Quadrature<dim - 1> &quadrature,
377 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
379 &output_data) const override;
380
384};
385
386
387
390/*----------------------------------------------------------------------*/
391
392#ifndef DOXYGEN
393
394template <int dim, int spacedim>
395inline void
397 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
398{
400 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
401 vertices[i] = cell->vertex(i);
402 this->cell = cell;
403}
404
405
406template <int dim, int spacedim>
407inline void
410{
411 cell_manifold_quadrature_weights.resize(
412 quad.size(), std::vector<double>(GeometryInfo<dim>::vertices_per_cell));
413 for (unsigned int q = 0; q < quad.size(); ++q)
414 {
415 for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
416 {
417 cell_manifold_quadrature_weights[q][i] =
419 }
420 }
421}
422
423
424
425template <int dim, int spacedim>
426inline bool
428{
429 return true;
430}
431
432
433template <int dim, int spacedim>
434bool
436 const ReferenceCell &cell_type) const
437{
438 if (cell_type.get_dimension() != dim)
439 return false; // TODO: or is this an error?
440
441 if (cell_type.is_hyper_cube())
442 return true;
443
444 return false;
445}
446
447
448
449#endif // DOXYGEN
450
451/* -------------- declaration of explicit specializations ------------- */
452
453
455
456#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:316
Definition point.h:111
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:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
unsigned int vertex_indices[2]
UpdateFlags
MappingKind
Definition mapping.h:77
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)