deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+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_q.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_q_h
16#define dealii_mapping_q_h
17
18
19#include <deal.II/base/config.h>
20
24#include <deal.II/base/table.h>
26
27#include <deal.II/fe/mapping.h>
28
30
32
33#include <array>
34#include <cmath>
35
37
38#ifndef DOXYGEN
39template <int, int>
40class MappingQCache;
41#endif
42
108template <int dim, int spacedim = dim>
109class MappingQ : public Mapping<dim, spacedim>
110{
111public:
117 MappingQ(const unsigned int polynomial_degree);
118
122 MappingQ(const MappingQ<dim, spacedim> &mapping);
123
124 // for documentation, see the Mapping base class
125 virtual std::unique_ptr<Mapping<dim, spacedim>>
126 clone() const override;
127
132 unsigned int
133 get_degree() const;
134
139 virtual bool
141
142 // for documentation, see the Mapping base class
145 &cell) const override;
146
147 virtual bool
148 is_compatible_with(const ReferenceCell &reference_cell) const override;
149
155 // for documentation, see the Mapping base class
156 virtual Point<spacedim>
159 const Point<dim> &p) const override;
160
161 // for documentation, see the Mapping base class
162 virtual Point<dim>
165 const Point<spacedim> &p) const override;
166
167 // for documentation, see the Mapping base class
168 virtual void
171 const ArrayView<const Point<spacedim>> &real_points,
172 const ArrayView<Point<dim>> &unit_points) const override;
173
183 // for documentation, see the Mapping base class
184 virtual void
185 transform(const ArrayView<const Tensor<1, dim>> &input,
186 const MappingKind kind,
188 const ArrayView<Tensor<1, spacedim>> &output) const override;
189
190 // for documentation, see the Mapping base class
191 virtual void
193 const MappingKind kind,
195 const ArrayView<Tensor<2, spacedim>> &output) const override;
196
197 // for documentation, see the Mapping base class
198 virtual void
199 transform(const ArrayView<const Tensor<2, dim>> &input,
200 const MappingKind kind,
202 const ArrayView<Tensor<2, spacedim>> &output) const override;
203
204 // for documentation, see the Mapping base class
205 virtual void
207 const MappingKind kind,
209 const ArrayView<Tensor<3, spacedim>> &output) const override;
210
211 // for documentation, see the Mapping base class
212 virtual void
213 transform(const ArrayView<const Tensor<3, dim>> &input,
214 const MappingKind kind,
216 const ArrayView<Tensor<3, spacedim>> &output) const override;
217
239 void
242 const ArrayView<const Point<dim>> &unit_points,
243 const UpdateFlags update_flags,
245 &output_data) const;
246
269 void
272 const unsigned int face_number,
273 const Quadrature<dim - 1> &face_quadrature,
274 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
276 &output_data) const;
277
294 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
295 {
296 public:
301 InternalData(const unsigned int polynomial_degree);
302
303 // Documentation see Mapping::InternalDataBase.
304 virtual void
305 reinit(const UpdateFlags update_flags,
306 const Quadrature<dim> &quadrature) override;
307
313 void
314 initialize_face(const UpdateFlags update_flags,
315 const Quadrature<dim> &quadrature,
316 const unsigned int n_original_q_points);
317
321 virtual std::size_t
322 memory_consumption() const override;
323
329 std::vector<Point<dim>> quadrature_points;
330
344 std::array<std::vector<Tensor<1, dim>>,
347
352 const unsigned int polynomial_degree;
353
363 const unsigned int n_shape_functions;
364
365 /*
366 * The default line support points. Is used in when the shape function
367 * values are computed.
368 *
369 * The number of quadrature points depends on the degree of this
370 * class, and it matches the number of degrees of freedom of an
371 * FE_Q<1>(this->degree).
372 */
374
390 VectorizedArray<double,
391 std::min<std::size_t>(VectorizedArray<double>::size(),
392 (dim <= 2 ? 2 : 4))>;
393
400
406
411
415 mutable std::vector<AlignedVector<Tensor<1, spacedim>>> aux;
416
421 mutable std::vector<Point<spacedim>> mapping_support_points;
422
428
434
442 };
443
444protected:
445 // documentation can be found in Mapping::requires_update_flags()
446 virtual UpdateFlags
447 requires_update_flags(const UpdateFlags update_flags) const override;
448
449 // documentation can be found in Mapping::get_data()
450 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
451 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
452
453 using Mapping<dim, spacedim>::get_face_data;
454
455 // documentation can be found in Mapping::get_face_data()
456 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
457 get_face_data(const UpdateFlags flags,
458 const hp::QCollection<dim - 1> &quadrature) const override;
459
460 // documentation can be found in Mapping::get_subface_data()
461 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
462 get_subface_data(const UpdateFlags flags,
463 const Quadrature<dim - 1> &quadrature) const override;
464
465 // documentation can be found in Mapping::fill_fe_values()
469 const CellSimilarity::Similarity cell_similarity,
470 const Quadrature<dim> &quadrature,
471 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
473 &output_data) const override;
474
475 using Mapping<dim, spacedim>::fill_fe_face_values;
476
477 // documentation can be found in Mapping::fill_fe_face_values()
478 virtual void
481 const unsigned int face_no,
482 const hp::QCollection<dim - 1> &quadrature,
483 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
485 &output_data) const override;
486
487 // documentation can be found in Mapping::fill_fe_subface_values()
488 virtual void
491 const unsigned int face_no,
492 const unsigned int subface_no,
493 const Quadrature<dim - 1> &quadrature,
494 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
496 &output_data) const override;
497
498 // documentation can be found in Mapping::fill_fe_immersed_surface_values()
499 virtual void
503 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
505 &output_data) const override;
506
515 const unsigned int polynomial_degree;
516
517 /*
518 * The default line support points. These are used when computing the
519 * location in real space of the support points on lines and quads, which
520 * are needed by the Manifold<dim,spacedim> class.
521 *
522 * The number of points depends on the degree of this class, and it matches
523 * the number of degrees of freedom of an FE_Q<1>(this->degree).
524 */
525 const std::vector<Point<1>> line_support_points;
526
527 /*
528 * The one-dimensional polynomials defined as Lagrange polynomials from the
529 * line support points. These are used for point evaluations and match the
530 * polynomial space of an FE_Q<1>(this->degree).
531 */
532 const std::vector<Polynomials::Polynomial<double>> polynomials_1d;
533
534 /*
535 * The numbering from the lexicographic to the hierarchical ordering used
536 * when expanding the tensor product with the mapping support points (which
537 * come in hierarchical numbers).
538 */
539 const std::vector<unsigned int> renumber_lexicographic_to_hierarchic;
540
541 /*
542 * The support points in reference coordinates. These are used for
543 * constructing approximations of the output of
544 * compute_mapping_support_points() when evaluating the mapping on the fly,
545 * rather than going through the FEValues interface provided by
546 * InternalData.
547 *
548 * The number of points depends on the degree of this class, and it matches
549 * the number of degrees of freedom of an FE_Q<dim>(this->degree).
550 */
551 const std::vector<Point<dim>> unit_cell_support_points;
552
572 const std::vector<Table<2, double>>
574
588
618 virtual std::vector<Point<spacedim>>
620 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
621
629 const Point<spacedim> &p,
630 const Point<dim> &initial_p_unit) const;
631
645 virtual void
648 std::vector<Point<spacedim>> &a) const;
649
664 virtual void
667 std::vector<Point<spacedim>> &a) const;
668
669 // Make MappingQCache a friend since it needs to call the
670 // compute_mapping_support_points() function.
671 template <int, int>
672 friend class MappingQCache;
673};
674
675
676
682template <int dim, int spacedim = dim>
684
688/*----------------------------------------------------------------------*/
689
690#ifndef DOXYGEN
691
692template <int dim, int spacedim>
693inline bool
695{
696 return true;
697}
698
699#endif // DOXYGEN
700
701/* -------------- declaration of explicit specializations ------------- */
702
703
705
706#endif
QGaussLobatto< 1 > line_support_points
Definition mapping_q.h:373
Triangulation< dim, spacedim >::cell_iterator cell_of_current_support_points
Definition mapping_q.h:427
virtual std::size_t memory_consumption() const override
Definition mapping_q.cc:63
std::vector< Point< spacedim > > mapping_support_points
Definition mapping_q.h:421
const unsigned int polynomial_degree
Definition mapping_q.h:352
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
Definition mapping_q.h:346
internal::MatrixFreeFunctions::ShapeInfo< double > shape_info
Definition mapping_q.h:399
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
Definition mapping_q.cc:81
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > * output_data
Definition mapping_q.h:441
std::vector< AlignedVector< Tensor< 1, spacedim > > > aux
Definition mapping_q.h:415
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
Definition mapping_q.cc:159
AlignedVector< double > volume_elements
Definition mapping_q.h:433
std::vector< Point< dim > > quadrature_points
Definition mapping_q.h:329
const unsigned int n_shape_functions
Definition mapping_q.h:363
AlignedVector< VectorizedArrayType > scratch
Definition mapping_q.h:405
const std::vector< unsigned int > renumber_lexicographic_to_hierarchic
Definition mapping_q.h:539
const Table< 2, double > support_point_weights_cell
Definition mapping_q.h:587
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Definition mapping_q.cc:774
void fill_mapping_data_for_face_quadrature(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_number, const Quadrature< dim - 1 > &face_quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
Definition mapping_q.cc:808
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 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 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
Definition mapping_q.cc:827
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
const unsigned int polynomial_degree
Definition mapping_q.h:515
virtual bool preserves_vertex_locations() 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
Definition mapping_q.cc:626
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
Definition mapping_q.cc:787
Point< dim > transform_real_to_unit_cell_internal(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit) const
Definition mapping_q.cc:321
const std::vector< Table< 2, double > > support_point_weights_perimeter_to_interior
Definition mapping_q.h:573
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
const std::vector< Point< 1 > > line_support_points
Definition mapping_q.h:525
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
Definition mapping_q.cc:280
const std::vector< Polynomials::Polynomial< double > > polynomials_1d
Definition mapping_q.h:532
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition mapping_q.cc:262
const std::vector< Point< dim > > unit_cell_support_points
Definition mapping_q.h:551
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
Definition mapping_q.cc:506
virtual void add_line_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
virtual void add_quad_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Definition mapping_q.cc:718
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 is_compatible_with(const ReferenceCell &reference_cell) const override
unsigned int get_degree() const
Definition mapping_q.cc:271
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
UpdateFlags
MappingKind
Definition mapping.h:79