PolyDEAL
Loading...
Searching...
No Matches
mapping_box.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 2024 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_box_h
16#define dealii_mapping_box_h
17
18
19#include <deal.II/base/config.h>
20
23
24#include <deal.II/fe/mapping.h>
25
26#include <cmath>
27
28
30
35
78template <int dim, int spacedim = dim>
79class MappingBox : public Mapping<dim, spacedim>
80{
81public:
82 MappingBox(const std::vector<BoundingBox<dim>> &local_boxes,
83 const std::map<types::global_cell_index, types::global_cell_index>
85 // for documentation, see the Mapping base class
86 virtual std::unique_ptr<Mapping<dim, spacedim>>
87 clone() const override;
88
93 virtual bool
94 preserves_vertex_locations() const override;
95
96#if DEAL_II_VERSION_GTE(9, 8, 0)
97 virtual bool
99#else
100 virtual bool
101 is_compatible_with(const ReferenceCell &reference_cell) const override;
102#endif
103
108
109 // for documentation, see the Mapping base class
110 virtual Point<spacedim>
113 const Point<dim> &p) const override;
114
115 // for documentation, see the Mapping base class
116 virtual Point<dim>
119 const Point<spacedim> &p) const override;
120
121 // for documentation, see the Mapping base class
122 virtual void
125 const ArrayView<const Point<spacedim>> &real_points,
126 const ArrayView<Point<dim>> &unit_points) const override;
127
131
136
137 // for documentation, see the Mapping base class
138 virtual void
139 transform(const ArrayView<const Tensor<1, dim>> &input,
140 const MappingKind kind,
142 const ArrayView<Tensor<1, spacedim>> &output) const override;
143
144 // for documentation, see the Mapping base class
145 virtual void
147 const MappingKind kind,
149 const ArrayView<Tensor<2, spacedim>> &output) const override;
150
151 // for documentation, see the Mapping base class
152 virtual void
153 transform(const ArrayView<const Tensor<2, dim>> &input,
154 const MappingKind kind,
156 const ArrayView<Tensor<2, spacedim>> &output) const override;
157
158 // for documentation, see the Mapping base class
159 virtual void
161 const MappingKind kind,
163 const ArrayView<Tensor<3, spacedim>> &output) const override;
164
165 // for documentation, see the Mapping base class
166 virtual void
167 transform(const ArrayView<const Tensor<3, dim>> &input,
168 const MappingKind kind,
170 const ArrayView<Tensor<3, spacedim>> &output) const override;
171
175
180
192 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
193 {
194 public:
198 InternalData() = default;
199
203 InternalData(const Quadrature<dim> &quadrature);
204
205 // Documentation see Mapping::InternalDataBase.
206 virtual void
207 reinit(const UpdateFlags update_flags,
208 const Quadrature<dim> &quadrature) override;
209
213 virtual std::size_t
214 memory_consumption() const override;
215
221
226
233
237 mutable double volume_element;
238
244 std::vector<Point<dim>> quadrature_points;
245 };
246
247private:
248 // documentation can be found in Mapping::requires_update_flags()
249 virtual UpdateFlags
250 requires_update_flags(const UpdateFlags update_flags) const override;
251
252 // documentation can be found in Mapping::get_data()
253 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
254 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
255
256 using Mapping<dim, spacedim>::get_face_data;
257
258 // documentation can be found in Mapping::get_subface_data()
259 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
260 get_subface_data(const UpdateFlags flags,
261 const Quadrature<dim - 1> &quadrature) const override;
262
263 // documentation can be found in Mapping::fill_fe_values()
267 const CellSimilarity::Similarity cell_similarity,
268 const Quadrature<dim> &quadrature,
269 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
271 &output_data) const override;
272
273 using Mapping<dim, spacedim>::fill_fe_face_values;
274
275 // documentation can be found in Mapping::fill_fe_subface_values()
276 virtual void
279 const unsigned int face_no,
280 const unsigned int subface_no,
281 const Quadrature<dim - 1> &quadrature,
282 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
284 &output_data) const override;
285
286 // documentation can be found in Mapping::fill_fe_immersed_surface_values()
287 virtual void
291 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
293 &output_data) const override;
294
298
303 void
306 const CellSimilarity::Similarity cell_similarity,
307 const InternalData &data) const;
308
315 void
318 const InternalData &data,
319 const ArrayView<const Point<dim>> &unit_quadrature_points,
320 std::vector<Point<dim>> &quadrature_points) const;
321
326 void
328 const unsigned int face_no,
329 const InternalData &data,
330 std::vector<Tensor<1, dim>> &normal_vectors) const;
331
337 void
339 const InternalData &data,
340 const CellSimilarity::Similarity cell_similarity,
342 &output_data) const;
343
344
349 void
351
356 void
358 const InternalData &data,
359 const CellSimilarity::Similarity cell_similarity,
361 &output_data) const;
362
367 void
369 const InternalData &data,
370 const CellSimilarity::Similarity cell_similarity,
372 &output_data) const;
373
377 std::vector<BoundingBox<dim>> boxes;
378
382 std::map<types::global_cell_index, types::global_cell_index>
384};
385
387
389
390#endif
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
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const hp::QCollection< dim - 1 > &quadrature) const
virtual std::size_t memory_consumption() const override
Tensor< 1, dim > traslation
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
Tensor< 1, dim > cell_extents
std::vector< Point< dim > > quadrature_points
Tensor< 1, dim > inverse_cell_extents
void maybe_update_volume_elements(const InternalData &data) const
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
MappingBox(const std::vector< BoundingBox< dim > > &local_boxes, const std::map< types::global_cell_index, types::global_cell_index > &polytope_translator)
void maybe_update_jacobian_derivatives(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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 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
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const override
void update_cell_extents(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const InternalData &data) const
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_immersed_surface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const NonMatching::ImmersedSurfaceQuadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) 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 std::unique_ptr< Mapping< dim, spacedim > > clone() const override
void maybe_update_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
std::vector< BoundingBox< dim > > boxes
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
void maybe_update_cell_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const ArrayView< const Point< dim > > &unit_quadrature_points, std::vector< Point< dim > > &quadrature_points) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
void maybe_update_inverse_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
void maybe_update_normal_vectors(const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, dim > > &normal_vectors) const
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
std::map< types::global_cell_index, types::global_cell_index > polytope_translator
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
UpdateFlags
MappingKind
std::vector< index_type > data
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell< dim > &reference_cell)