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_cartesian.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 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_cartesian_h
16#define dealii_mapping_cartesian_h
17
18
19#include <deal.II/base/config.h>
20
22
23#include <deal.II/fe/mapping.h>
24
25#include <cmath>
26
27
29
77template <int dim, int spacedim = dim>
78class MappingCartesian : public Mapping<dim, spacedim>
79{
80public:
81 // for documentation, see the Mapping base class
82 virtual std::unique_ptr<Mapping<dim, spacedim>>
83 clone() const override;
84
89 virtual bool
90 preserves_vertex_locations() const override;
91
92 virtual bool
93 is_compatible_with(const ReferenceCell &reference_cell) const override;
94
100 // for documentation, see the Mapping base class
101 virtual Point<spacedim>
104 const Point<dim> &p) const override;
105
106 // for documentation, see the Mapping base class
107 virtual Point<dim>
110 const Point<spacedim> &p) const override;
111
121 // for documentation, see the Mapping base class
122 virtual void
123 transform(const ArrayView<const Tensor<1, dim>> &input,
124 const MappingKind kind,
126 const ArrayView<Tensor<1, spacedim>> &output) const override;
127
128 // for documentation, see the Mapping base class
129 virtual void
131 const MappingKind kind,
133 const ArrayView<Tensor<2, spacedim>> &output) const override;
134
135 // for documentation, see the Mapping base class
136 virtual void
137 transform(const ArrayView<const Tensor<2, dim>> &input,
138 const MappingKind kind,
140 const ArrayView<Tensor<2, spacedim>> &output) const override;
141
142 // for documentation, see the Mapping base class
143 virtual void
145 const MappingKind kind,
147 const ArrayView<Tensor<3, spacedim>> &output) const override;
148
149 // for documentation, see the Mapping base class
150 virtual void
151 transform(const ArrayView<const Tensor<3, dim>> &input,
152 const MappingKind kind,
154 const ArrayView<Tensor<3, spacedim>> &output) const override;
155
177 void
180 const ArrayView<const Point<dim>> &unit_points,
181 const UpdateFlags update_flags,
183 &output_data) const;
184
201 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
202 {
203 public:
207 InternalData() = default;
208
212 InternalData(const Quadrature<dim> &quadrature);
213
217 virtual std::size_t
218 memory_consumption() const override;
219
225
229 mutable double volume_element;
230
234 std::vector<Point<dim>> quadrature_points;
235 };
236
237private:
238 // documentation can be found in Mapping::requires_update_flags()
239 virtual UpdateFlags
240 requires_update_flags(const UpdateFlags update_flags) const override;
241
242 // documentation can be found in Mapping::get_data()
243 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
244 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
245
246 using Mapping<dim, spacedim>::get_face_data;
247
248 // documentation can be found in Mapping::get_face_data()
249 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
250 get_face_data(const UpdateFlags flags,
251 const hp::QCollection<dim - 1> &quadrature) const override;
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_face_values()
271 virtual void
274 const unsigned int face_no,
275 const hp::QCollection<dim - 1> &quadrature,
276 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
278 &output_data) const override;
279
280 // documentation can be found in Mapping::fill_fe_subface_values()
281 virtual void
284 const unsigned int face_no,
285 const unsigned int subface_no,
286 const Quadrature<dim - 1> &quadrature,
287 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
289 &output_data) const override;
290
291 // documentation can be found in Mapping::fill_fe_immersed_surface_values()
292 virtual void
296 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
298 &output_data) const override;
299
308 void
311 const CellSimilarity::Similarity cell_similarity,
312 const InternalData &data) const;
313
320 void
323 const InternalData &data,
324 std::vector<Point<dim>> &quadrature_points) const;
325
332 void
335 const unsigned int face_no,
336 const InternalData &data,
337 std::vector<Point<dim>> &quadrature_points) const;
338
345 void
348 const unsigned int face_no,
349 const unsigned int sub_no,
350 const InternalData &data,
351 std::vector<Point<dim>> &quadrature_points) const;
352
359 void
362 const InternalData &data,
363 const typename QProjector<dim>::DataSetDescriptor &offset,
364 std::vector<Point<dim>> &quadrature_points) const;
365
370 void
372 const unsigned int face_no,
373 const InternalData &data,
374 std::vector<Tensor<1, dim>> &normal_vectors) const;
375
381 void
383 const InternalData &data,
384 const CellSimilarity::Similarity cell_similarity,
386 &output_data) const;
387
388
393 void
395
400 void
402 const InternalData &data,
403 const CellSimilarity::Similarity cell_similarity,
405 &output_data) const;
406
411 void
413 const InternalData &data,
414 const CellSimilarity::Similarity cell_similarity,
416 &output_data) const;
417};
418
422
423#endif
std::vector< Point< dim > > quadrature_points
virtual std::size_t memory_consumption() 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
void maybe_update_volume_elements(const InternalData &data) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
void update_cell_extents(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const InternalData &data) const
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
void maybe_update_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
void maybe_update_subface_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const InternalData &data, std::vector< Point< dim > > &quadrature_points) const
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
void maybe_update_face_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const InternalData &data, std::vector< Point< dim > > &quadrature_points) const
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
void maybe_update_cell_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, std::vector< Point< dim > > &quadrature_points) const
virtual bool preserves_vertex_locations() const override
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
void fill_mapping_data_for_generic_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points, const UpdateFlags update_flags, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
void transform_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const typename QProjector< dim >::DataSetDescriptor &offset, std::vector< Point< dim > > &quadrature_points) const
void maybe_update_jacobian_derivatives(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) 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
void maybe_update_normal_vectors(const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, dim > > &normal_vectors) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
void maybe_update_inverse_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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
Abstract base class for mapping classes.
Definition mapping.h:316
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
UpdateFlags
MappingKind
Definition mapping.h:77