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
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 virtual bool
97 is_compatible_with(const ReferenceCell &reference_cell) const override;
98
104 // for documentation, see the Mapping base class
105 virtual Point<spacedim>
108 const Point<dim> &p) const override;
109
110 // for documentation, see the Mapping base class
111 virtual Point<dim>
114 const Point<spacedim> &p) const override;
115
116 // for documentation, see the Mapping base class
117 virtual void
120 const ArrayView<const Point<spacedim>> &real_points,
121 const ArrayView<Point<dim>> &unit_points) const override;
122
132 // for documentation, see the Mapping base class
133 virtual void
134 transform(const ArrayView<const Tensor<1, dim>> &input,
135 const MappingKind kind,
137 const ArrayView<Tensor<1, spacedim>> &output) const override;
138
139 // for documentation, see the Mapping base class
140 virtual void
142 const MappingKind kind,
144 const ArrayView<Tensor<2, spacedim>> &output) const override;
145
146 // for documentation, see the Mapping base class
147 virtual void
148 transform(const ArrayView<const Tensor<2, dim>> &input,
149 const MappingKind kind,
151 const ArrayView<Tensor<2, spacedim>> &output) const override;
152
153 // for documentation, see the Mapping base class
154 virtual void
156 const MappingKind kind,
158 const ArrayView<Tensor<3, spacedim>> &output) const override;
159
160 // for documentation, see the Mapping base class
161 virtual void
162 transform(const ArrayView<const Tensor<3, dim>> &input,
163 const MappingKind kind,
165 const ArrayView<Tensor<3, spacedim>> &output) const override;
166
187 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
188 {
189 public:
193 InternalData() = default;
194
198 InternalData(const Quadrature<dim> &quadrature);
199
200 // Documentation see Mapping::InternalDataBase.
201 virtual void
202 reinit(const UpdateFlags update_flags,
203 const Quadrature<dim> &quadrature) override;
204
208 virtual std::size_t
209 memory_consumption() const override;
210
216
221
228
232 mutable double volume_element;
233
239 std::vector<Point<dim>> quadrature_points;
240 };
241
242private:
243 // documentation can be found in Mapping::requires_update_flags()
244 virtual UpdateFlags
245 requires_update_flags(const UpdateFlags update_flags) const override;
246
247 // documentation can be found in Mapping::get_data()
248 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
249 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
250
251 using Mapping<dim, spacedim>::get_face_data;
252
253 // documentation can be found in Mapping::get_subface_data()
254 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
255 get_subface_data(const UpdateFlags flags,
256 const Quadrature<dim - 1> &quadrature) const override;
257
258 // documentation can be found in Mapping::fill_fe_values()
262 const CellSimilarity::Similarity cell_similarity,
263 const Quadrature<dim> &quadrature,
264 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
266 &output_data) const override;
267
268 using Mapping<dim, spacedim>::fill_fe_face_values;
269
270 // documentation can be found in Mapping::fill_fe_subface_values()
271 virtual void
274 const unsigned int face_no,
275 const unsigned int subface_no,
276 const Quadrature<dim - 1> &quadrature,
277 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
279 &output_data) const override;
280
281 // documentation can be found in Mapping::fill_fe_immersed_surface_values()
282 virtual void
286 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
288 &output_data) const override;
289
298 void
301 const CellSimilarity::Similarity cell_similarity,
302 const InternalData &data) const;
303
310 void
313 const InternalData &data,
314 const ArrayView<const Point<dim>> &unit_quadrature_points,
315 std::vector<Point<dim>> &quadrature_points) const;
316
321 void
323 const unsigned int face_no,
324 const InternalData &data,
325 std::vector<Tensor<1, dim>> &normal_vectors) const;
326
332 void
334 const InternalData &data,
335 const CellSimilarity::Similarity cell_similarity,
337 &output_data) const;
338
339
344 void
346
351 void
353 const InternalData &data,
354 const CellSimilarity::Similarity cell_similarity,
356 &output_data) const;
357
362 void
364 const InternalData &data,
365 const CellSimilarity::Similarity cell_similarity,
367 &output_data) const;
368
372 std::vector<BoundingBox<dim>> boxes;
373
377 std::map<types::global_cell_index, types::global_cell_index>
379};
380
384
385#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
UpdateFlags
MappingKind
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)