Reference documentation for deal.II version 9.3.3
\(\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_fe.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2020 - 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_fe_h
17#define dealii_mapping_fe_h
18
19
20#include <deal.II/base/config.h>
21
24#include <deal.II/base/table.h>
26
27#include <deal.II/fe/mapping.h>
28
30
31#include <array>
32#include <cmath>
33
35
36
39
40
57template <int dim, int spacedim = dim>
58class MappingFE : public Mapping<dim, spacedim>
59{
60public:
65
69 MappingFE(const MappingFE<dim, spacedim> &mapping);
70
71 // for documentation, see the Mapping base class
72 virtual std::unique_ptr<Mapping<dim, spacedim>>
73 clone() const override;
74
79 unsigned int
80 get_degree() const;
81
82 // for documentation, see the Mapping base class
85 &cell) const override;
86
87
88 virtual bool
89 is_compatible_with(const ReferenceCell &reference_cell) const override;
90
95 virtual bool
97
103 // for documentation, see the Mapping base class
104 virtual Point<spacedim>
107 const Point<dim> &p) const override;
108
109 // for documentation, see the Mapping base class
110 virtual Point<dim>
113 const Point<spacedim> &p) const override;
114
124 // for documentation, see the Mapping base class
125 virtual void
126 transform(const ArrayView<const Tensor<1, dim>> & input,
127 const MappingKind kind,
129 const ArrayView<Tensor<1, spacedim>> &output) const override;
130
131 // for documentation, see the Mapping base class
132 virtual void
134 const MappingKind kind,
136 const ArrayView<Tensor<2, spacedim>> &output) const override;
137
138 // for documentation, see the Mapping base class
139 virtual void
140 transform(const ArrayView<const Tensor<2, dim>> & input,
141 const MappingKind kind,
143 const ArrayView<Tensor<2, spacedim>> &output) const override;
144
145 // for documentation, see the Mapping base class
146 virtual void
148 const MappingKind kind,
150 const ArrayView<Tensor<3, spacedim>> &output) const override;
151
152 // for documentation, see the Mapping base class
153 virtual void
154 transform(const ArrayView<const Tensor<3, dim>> & input,
155 const MappingKind kind,
157 const ArrayView<Tensor<3, spacedim>> &output) const override;
158
179 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
180 {
181 public:
186
195 void
196 initialize(const UpdateFlags update_flags,
197 const Quadrature<dim> &quadrature,
198 const unsigned int n_original_q_points);
199
205 void
206 initialize_face(const UpdateFlags update_flags,
207 const Quadrature<dim> &quadrature,
208 const unsigned int n_original_q_points);
209
214 void
215 compute_shape_function_values(const std::vector<Point<dim>> &unit_points);
216
217
222 const double &
223 shape(const unsigned int qpoint, const unsigned int shape_nr) const;
224
228 double &
229 shape(const unsigned int qpoint, const unsigned int shape_nr);
230
234 const Tensor<1, dim> &
235 derivative(const unsigned int qpoint, const unsigned int shape_nr) const;
236
241 derivative(const unsigned int qpoint, const unsigned int shape_nr);
242
246 const Tensor<2, dim> &
247 second_derivative(const unsigned int qpoint,
248 const unsigned int shape_nr) const;
249
254 second_derivative(const unsigned int qpoint, const unsigned int shape_nr);
255
259 const Tensor<3, dim> &
260 third_derivative(const unsigned int qpoint,
261 const unsigned int shape_nr) const;
262
267 third_derivative(const unsigned int qpoint, const unsigned int shape_nr);
268
272 const Tensor<4, dim> &
273 fourth_derivative(const unsigned int qpoint,
274 const unsigned int shape_nr) const;
275
280 fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr);
281
285 virtual std::size_t
286 memory_consumption() const override;
287
293 std::vector<double> shape_values;
294
300 std::vector<Tensor<1, dim>> shape_derivatives;
301
308 std::vector<Tensor<2, dim>> shape_second_derivatives;
309
316 std::vector<Tensor<3, dim>> shape_third_derivatives;
317
324 std::vector<Tensor<4, dim>> shape_fourth_derivatives;
325
332 std::array<std::vector<Tensor<1, dim>>,
335
340
344 const unsigned int polynomial_degree;
345
349 const unsigned int n_shape_functions;
350
360 mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
361
369 mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
370
374 mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
375
380 mutable std::vector<Point<spacedim>> mapping_support_points;
381
387
392 mutable std::vector<double> volume_elements;
393
397 mutable std::vector<double> quadrature_weights;
398 };
399
400
401 // documentation can be found in Mapping::requires_update_flags()
402 virtual UpdateFlags
403 requires_update_flags(const UpdateFlags update_flags) const override;
404
405 // documentation can be found in Mapping::get_data()
406 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
407 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
408
409 using Mapping<dim, spacedim>::get_face_data;
410
411 // documentation can be found in Mapping::get_face_data()
412 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
413 get_face_data(const UpdateFlags flags,
414 const hp::QCollection<dim - 1> &quadrature) const override;
415
416 // documentation can be found in Mapping::get_subface_data()
417 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
418 get_subface_data(const UpdateFlags flags,
419 const Quadrature<dim - 1> &quadrature) const override;
420
421 // documentation can be found in Mapping::fill_fe_values()
425 const CellSimilarity::Similarity cell_similarity,
426 const Quadrature<dim> & quadrature,
427 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
429 &output_data) const override;
430
431 using Mapping<dim, spacedim>::fill_fe_face_values;
432
433 // documentation can be found in Mapping::fill_fe_face_values()
434 virtual void
437 const unsigned int face_no,
438 const hp::QCollection<dim - 1> & quadrature,
439 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
441 &output_data) const override;
442
443 // documentation can be found in Mapping::fill_fe_subface_values()
444 virtual void
447 const unsigned int face_no,
448 const unsigned int subface_no,
449 const Quadrature<dim - 1> & quadrature,
450 const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
452 &output_data) const override;
453
458protected:
459 const std::unique_ptr<FiniteElement<dim, spacedim>> fe;
460
465 const unsigned int polynomial_degree;
466
470 virtual std::vector<Point<spacedim>>
472 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
473
474private:
476};
477
478
479
482/*----------------------------------------------------------------------*/
483
484#ifndef DOXYGEN
485
486template <int dim, int spacedim>
487inline const double &
488MappingFE<dim, spacedim>::InternalData::shape(const unsigned int qpoint,
489 const unsigned int shape_nr) const
490{
491 AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_values.size());
492 return shape_values[qpoint * n_shape_functions + shape_nr];
493}
494
495
496
497template <int dim, int spacedim>
498inline double &
499MappingFE<dim, spacedim>::InternalData::shape(const unsigned int qpoint,
500 const unsigned int shape_nr)
501{
502 AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_values.size());
503 return shape_values[qpoint * n_shape_functions + shape_nr];
504}
505
506
507template <int dim, int spacedim>
508inline const Tensor<1, dim> &
510 const unsigned int qpoint,
511 const unsigned int shape_nr) const
512{
513 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
514 shape_derivatives.size());
515 return shape_derivatives[qpoint * n_shape_functions + shape_nr];
516}
517
518
519
520template <int dim, int spacedim>
521inline Tensor<1, dim> &
523 const unsigned int shape_nr)
524{
525 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
526 shape_derivatives.size());
527 return shape_derivatives[qpoint * n_shape_functions + shape_nr];
528}
529
530
531template <int dim, int spacedim>
532inline const Tensor<2, dim> &
534 const unsigned int qpoint,
535 const unsigned int shape_nr) const
536{
537 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
538 shape_second_derivatives.size());
539 return shape_second_derivatives[qpoint * n_shape_functions + shape_nr];
540}
541
542
543template <int dim, int spacedim>
544inline Tensor<2, dim> &
546 const unsigned int qpoint,
547 const unsigned int shape_nr)
548{
549 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
550 shape_second_derivatives.size());
551 return shape_second_derivatives[qpoint * n_shape_functions + shape_nr];
552}
553
554template <int dim, int spacedim>
555inline const Tensor<3, dim> &
557 const unsigned int qpoint,
558 const unsigned int shape_nr) const
559{
560 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
561 shape_third_derivatives.size());
562 return shape_third_derivatives[qpoint * n_shape_functions + shape_nr];
563}
564
565
566template <int dim, int spacedim>
567inline Tensor<3, dim> &
569 const unsigned int qpoint,
570 const unsigned int shape_nr)
571{
572 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
573 shape_third_derivatives.size());
574 return shape_third_derivatives[qpoint * n_shape_functions + shape_nr];
575}
576
577
578template <int dim, int spacedim>
579inline const Tensor<4, dim> &
581 const unsigned int qpoint,
582 const unsigned int shape_nr) const
583{
584 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
585 shape_fourth_derivatives.size());
586 return shape_fourth_derivatives[qpoint * n_shape_functions + shape_nr];
587}
588
589
590template <int dim, int spacedim>
591inline Tensor<4, dim> &
593 const unsigned int qpoint,
594 const unsigned int shape_nr)
595{
596 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
597 shape_fourth_derivatives.size());
598 return shape_fourth_derivatives[qpoint * n_shape_functions + shape_nr];
599}
600
601
602
603template <int dim, int spacedim>
604inline bool
606{
607 return true;
608}
609
610
611
612#endif // DOXYGEN
613
614/* -------------- declaration of explicit specializations ------------- */
615
616
618
619#endif
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr)
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
Definition: mapping_fe.cc:141
std::vector< Point< spacedim > > mapping_support_points
Definition: mapping_fe.h:380
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< std::vector< Tensor< 1, spacedim > > > aux
Definition: mapping_fe.h:374
std::vector< Tensor< 3, dim > > shape_third_derivatives
Definition: mapping_fe.h:316
std::vector< Tensor< 2, dim > > shape_second_derivatives
Definition: mapping_fe.h:308
Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr)
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
Definition: mapping_fe.h:324
Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr)
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
Definition: mapping_fe.cc:83
Triangulation< dim, spacedim >::cell_iterator cell_of_current_support_points
Definition: mapping_fe.h:386
Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr)
void compute_shape_function_values(const std::vector< Point< dim > > &unit_points)
Definition: mapping_fe.cc:181
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
Definition: mapping_fe.h:334
std::vector< double > quadrature_weights
Definition: mapping_fe.h:397
const FiniteElement< dim, spacedim > & fe
Definition: mapping_fe.h:339
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
Definition: mapping_fe.h:369
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
Definition: mapping_fe.h:360
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
double & shape(const unsigned int qpoint, const unsigned int shape_nr)
std::vector< double > shape_values
Definition: mapping_fe.h:293
std::vector< Tensor< 1, dim > > shape_derivatives
Definition: mapping_fe.h:300
virtual std::size_t memory_consumption() const override
Definition: mapping_fe.cc:63
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
InternalData(const FiniteElement< dim, spacedim > &fe)
Definition: mapping_fe.cc:52
const unsigned int n_shape_functions
Definition: mapping_fe.h:349
std::vector< double > volume_elements
Definition: mapping_fe.h:392
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
const unsigned int polynomial_degree
Definition: mapping_fe.h:344
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
Definition: mapping_fe.cc:1636
const unsigned int polynomial_degree
Definition: mapping_fe.h:465
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
Definition: mapping_fe.cc:1590
MappingFE(const FiniteElement< dim, spacedim > &fe)
Definition: mapping_fe.cc:846
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_fe.cc:1076
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
Definition: mapping_fe.cc:924
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Definition: mapping_fe.cc:1005
Table< 2, double > mapping_support_point_weights
Definition: mapping_fe.h:475
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Definition: mapping_fe.cc:2279
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
Definition: mapping_fe.cc:2092
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_fe.cc:1114
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
Definition: mapping_fe.cc:2306
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Definition: mapping_fe.cc:1061
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
Definition: mapping_fe.cc:1095
unsigned int get_degree() const
Definition: mapping_fe.cc:896
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
Definition: mapping_fe.cc:2316
virtual bool preserves_vertex_locations() const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_fe.cc:887
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
Definition: mapping_fe.cc:905
const std::unique_ptr< FiniteElement< dim, spacedim > > fe
Definition: mapping_fe.h:459
Abstract base class for mapping classes.
Definition: mapping.h:304
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
UpdateFlags
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
MappingKind
Definition: mapping.h:65
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)