Reference documentation for deal.II version GIT 9c182271f7 2023-03-28 14:30:01+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\}}\)
mapping_cartesian.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2021 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_cartesian_h
17 #define dealii_mapping_cartesian_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 
24 #include <deal.II/fe/mapping.h>
25 
26 #include <cmath>
27 
28 
30 
78 template <int dim, int spacedim = dim>
79 class MappingCartesian : public Mapping<dim, spacedim>
80 {
81 public:
82  // for documentation, see the Mapping base class
83  virtual std::unique_ptr<Mapping<dim, spacedim>>
84  clone() const override;
85 
90  virtual bool
91  preserves_vertex_locations() const override;
92 
93  virtual bool
94  is_compatible_with(const ReferenceCell &reference_cell) const override;
95 
101  // for documentation, see the Mapping base class
102  virtual Point<spacedim>
104  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
105  const Point<dim> &p) const override;
106 
107  // for documentation, see the Mapping base class
108  virtual Point<dim>
110  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
111  const Point<spacedim> &p) const override;
112 
122  // for documentation, see the Mapping base class
123  virtual void
124  transform(const ArrayView<const Tensor<1, dim>> & input,
125  const MappingKind kind,
126  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
127  const ArrayView<Tensor<1, spacedim>> &output) const override;
128 
129  // for documentation, see the Mapping base class
130  virtual void
132  const MappingKind kind,
133  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
134  const ArrayView<Tensor<2, spacedim>> &output) const override;
135 
136  // for documentation, see the Mapping base class
137  virtual void
138  transform(const ArrayView<const Tensor<2, dim>> & input,
139  const MappingKind kind,
140  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
141  const ArrayView<Tensor<2, spacedim>> &output) const override;
142 
143  // for documentation, see the Mapping base class
144  virtual void
146  const MappingKind kind,
147  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
148  const ArrayView<Tensor<3, spacedim>> &output) const override;
149 
150  // for documentation, see the Mapping base class
151  virtual void
152  transform(const ArrayView<const Tensor<3, dim>> & input,
153  const MappingKind kind,
154  const typename Mapping<dim, spacedim>::InternalDataBase &internal,
155  const ArrayView<Tensor<3, spacedim>> &output) const override;
156 
178  void
180  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
181  const ArrayView<const Point<dim>> & unit_points,
182  const UpdateFlags update_flags,
184  &output_data) const;
185 
202  class InternalData : public Mapping<dim, spacedim>::InternalDataBase
203  {
204  public:
208  InternalData() = default;
209 
213  InternalData(const Quadrature<dim> &quadrature);
214 
218  virtual std::size_t
219  memory_consumption() const override;
220 
226 
230  mutable double volume_element;
231 
235  std::vector<Point<dim>> quadrature_points;
236  };
237 
238 private:
239  // documentation can be found in Mapping::requires_update_flags()
240  virtual UpdateFlags
241  requires_update_flags(const UpdateFlags update_flags) const override;
242 
243  // documentation can be found in Mapping::get_data()
244  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
245  get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
246 
248 
249  // documentation can be found in Mapping::get_face_data()
250  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
251  get_face_data(const UpdateFlags flags,
252  const hp::QCollection<dim - 1> &quadrature) const override;
253 
254  // documentation can be found in Mapping::get_subface_data()
255  virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
256  get_subface_data(const UpdateFlags flags,
257  const Quadrature<dim - 1> &quadrature) const override;
258 
259  // documentation can be found in Mapping::fill_fe_values()
262  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
263  const CellSimilarity::Similarity cell_similarity,
264  const Quadrature<dim> & quadrature,
265  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
267  &output_data) const override;
268 
270 
271  // documentation can be found in Mapping::fill_fe_face_values()
272  virtual void
274  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
275  const unsigned int face_no,
276  const hp::QCollection<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_subface_values()
282  virtual void
284  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
285  const unsigned int face_no,
286  const unsigned int subface_no,
287  const Quadrature<dim - 1> & quadrature,
288  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
290  &output_data) const override;
291 
292  // documentation can be found in Mapping::fill_fe_immersed_surface_values()
293  virtual void
295  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
297  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
299  &output_data) const override;
300 
309  void
311  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
312  const CellSimilarity::Similarity cell_similarity,
313  const InternalData & data) const;
314 
321  void
323  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
324  const InternalData & data,
325  std::vector<Point<dim>> &quadrature_points) const;
326 
333  void
335  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
336  const unsigned int face_no,
337  const InternalData & data,
338  std::vector<Point<dim>> &quadrature_points) const;
339 
346  void
348  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
349  const unsigned int face_no,
350  const unsigned int sub_no,
351  const InternalData & data,
352  std::vector<Point<dim>> &quadrature_points) const;
353 
360  void
362  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
363  const InternalData & data,
364  const typename QProjector<dim>::DataSetDescriptor & offset,
365  std::vector<Point<dim>> &quadrature_points) const;
366 
371  void
373  const unsigned int face_no,
374  const InternalData & data,
375  std::vector<Tensor<1, dim>> &normal_vectors) const;
376 
382  void
384  const InternalData & data,
385  const CellSimilarity::Similarity cell_similarity,
387  &output_data) const;
388 
389 
394  void
395  maybe_update_volume_elements(const InternalData &data) const;
396 
401  void
403  const InternalData & data,
404  const CellSimilarity::Similarity cell_similarity,
406  &output_data) const;
407 
412  void
414  const InternalData & data,
415  const CellSimilarity::Similarity cell_similarity,
417  &output_data) const;
418 };
419 
423 
424 #endif
std::vector< Point< dim > > quadrature_points
virtual std::size_t memory_consumption() 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 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
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 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
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 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 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
virtual bool is_compatible_with(const ReferenceCell &reference_cell) 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 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 maybe_update_jacobian_derivatives(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 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 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
Abstract base class for mapping classes.
Definition: mapping.h:314
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
UpdateFlags
MappingKind
Definition: mapping.h:75
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim >>> &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const std::vector< std::vector< double >> &properties={})
Definition: generators.cc:493