deal.II version GIT relicensing-2169-gec1b43f35b 2024-11-22 07:10: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_fe.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2015 - 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_fe_h
16#define dealii_mapping_fe_h
17
18
19#include <deal.II/base/config.h>
20
23#include <deal.II/base/table.h>
25
26#include <deal.II/fe/mapping.h>
27
29
30#include <array>
31#include <cmath>
32
34
35
58template <int dim, int spacedim = dim>
59class MappingFE : public Mapping<dim, spacedim>
60{
61public:
66
70 MappingFE(const MappingFE<dim, spacedim> &mapping);
71
72 // for documentation, see the Mapping base class
73 virtual std::unique_ptr<Mapping<dim, spacedim>>
74 clone() const override;
75
80 unsigned int
81 get_degree() const;
82
83 // for documentation, see the Mapping base class
86 &cell) const override;
87
88
89 virtual bool
90 is_compatible_with(const ReferenceCell &reference_cell) const override;
91
96 virtual bool
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
125 // for documentation, see the Mapping base class
126 virtual void
127 transform(const ArrayView<const Tensor<1, dim>> &input,
128 const MappingKind kind,
130 const ArrayView<Tensor<1, spacedim>> &output) const override;
131
132 // for documentation, see the Mapping base class
133 virtual void
135 const MappingKind kind,
137 const ArrayView<Tensor<2, spacedim>> &output) const override;
138
139 // for documentation, see the Mapping base class
140 virtual void
141 transform(const ArrayView<const Tensor<2, dim>> &input,
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
149 const MappingKind kind,
151 const ArrayView<Tensor<3, spacedim>> &output) const override;
152
153 // for documentation, see the Mapping base class
154 virtual void
155 transform(const ArrayView<const Tensor<3, dim>> &input,
156 const MappingKind kind,
158 const ArrayView<Tensor<3, spacedim>> &output) const override;
159
180 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
181 {
182 public:
187
188 // Documentation see Mapping::InternalDataBase.
189 virtual void
190 reinit(const UpdateFlags update_flags,
191 const Quadrature<dim> &quadrature) override;
192
198 void
199 initialize_face(const UpdateFlags update_flags,
200 const Quadrature<dim> &quadrature,
201 const unsigned int n_original_q_points);
202
207 void
208 compute_shape_function_values(const std::vector<Point<dim>> &unit_points);
209
210
215 const double &
216 shape(const unsigned int qpoint, const unsigned int shape_nr) const;
217
221 double &
222 shape(const unsigned int qpoint, const unsigned int shape_nr);
223
227 const Tensor<1, dim> &
228 derivative(const unsigned int qpoint, const unsigned int shape_nr) const;
229
234 derivative(const unsigned int qpoint, const unsigned int shape_nr);
235
239 const Tensor<2, dim> &
240 second_derivative(const unsigned int qpoint,
241 const unsigned int shape_nr) const;
242
247 second_derivative(const unsigned int qpoint, const unsigned int shape_nr);
248
252 const Tensor<3, dim> &
253 third_derivative(const unsigned int qpoint,
254 const unsigned int shape_nr) const;
255
260 third_derivative(const unsigned int qpoint, const unsigned int shape_nr);
261
265 const Tensor<4, dim> &
266 fourth_derivative(const unsigned int qpoint,
267 const unsigned int shape_nr) const;
268
273 fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr);
274
278 virtual std::size_t
279 memory_consumption() const override;
280
286 std::vector<double> shape_values;
287
293 std::vector<Tensor<1, dim>> shape_derivatives;
294
301 std::vector<Tensor<2, dim>> shape_second_derivatives;
302
309 std::vector<Tensor<3, dim>> shape_third_derivatives;
310
317 std::vector<Tensor<4, dim>> shape_fourth_derivatives;
318
325 std::array<std::vector<Tensor<1, dim>>,
328
333
337 const unsigned int polynomial_degree;
338
342 const unsigned int n_shape_functions;
343
353 mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
354
362 mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
363
367 mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
368
373 mutable std::vector<Point<spacedim>> mapping_support_points;
374
380
385 mutable std::vector<double> volume_elements;
386
390 mutable std::vector<double> quadrature_weights;
391 };
392
393
394 // documentation can be found in Mapping::requires_update_flags()
395 virtual UpdateFlags
396 requires_update_flags(const UpdateFlags update_flags) const override;
397
398 // documentation can be found in Mapping::get_data()
399 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
400 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
401
402 using Mapping<dim, spacedim>::get_face_data;
403
404 // documentation can be found in Mapping::get_face_data()
405 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
406 get_face_data(const UpdateFlags flags,
407 const hp::QCollection<dim - 1> &quadrature) const override;
408
409 // documentation can be found in Mapping::get_subface_data()
410 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
411 get_subface_data(const UpdateFlags flags,
412 const Quadrature<dim - 1> &quadrature) const override;
413
414 // documentation can be found in Mapping::fill_fe_values()
418 const CellSimilarity::Similarity cell_similarity,
419 const Quadrature<dim> &quadrature,
420 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
422 &output_data) const override;
423
424 using Mapping<dim, spacedim>::fill_fe_face_values;
425
426 // documentation can be found in Mapping::fill_fe_face_values()
427 virtual void
430 const unsigned int face_no,
431 const hp::QCollection<dim - 1> &quadrature,
432 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
434 &output_data) const override;
435
436 // documentation can be found in Mapping::fill_fe_subface_values()
437 virtual void
440 const unsigned int face_no,
441 const unsigned int subface_no,
442 const Quadrature<dim - 1> &quadrature,
443 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
445 &output_data) const override;
446
451protected:
452 const std::unique_ptr<FiniteElement<dim, spacedim>> fe;
453
458 const unsigned int polynomial_degree;
459
463 virtual std::vector<Point<spacedim>>
465 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;
466
467private:
469};
470
471
472
475/*----------------------------------------------------------------------*/
476
477#ifndef DOXYGEN
478
479template <int dim, int spacedim>
480inline const double &
481MappingFE<dim, spacedim>::InternalData::shape(const unsigned int qpoint,
482 const unsigned int shape_nr) const
483{
484 AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_values.size());
485 return shape_values[qpoint * n_shape_functions + shape_nr];
486}
487
488
489
490template <int dim, int spacedim>
491inline double &
492MappingFE<dim, spacedim>::InternalData::shape(const unsigned int qpoint,
493 const unsigned int shape_nr)
494{
495 AssertIndexRange(qpoint * n_shape_functions + shape_nr, shape_values.size());
496 return shape_values[qpoint * n_shape_functions + shape_nr];
497}
498
499
500template <int dim, int spacedim>
501inline const Tensor<1, dim> &
503 const unsigned int qpoint,
504 const unsigned int shape_nr) const
505{
506 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
507 shape_derivatives.size());
508 return shape_derivatives[qpoint * n_shape_functions + shape_nr];
509}
510
511
512
513template <int dim, int spacedim>
514inline Tensor<1, dim> &
516 const unsigned int shape_nr)
517{
518 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
519 shape_derivatives.size());
520 return shape_derivatives[qpoint * n_shape_functions + shape_nr];
521}
522
523
524template <int dim, int spacedim>
525inline const Tensor<2, dim> &
527 const unsigned int qpoint,
528 const unsigned int shape_nr) const
529{
530 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
531 shape_second_derivatives.size());
532 return shape_second_derivatives[qpoint * n_shape_functions + shape_nr];
533}
534
535
536template <int dim, int spacedim>
537inline Tensor<2, dim> &
539 const unsigned int qpoint,
540 const unsigned int shape_nr)
541{
542 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
543 shape_second_derivatives.size());
544 return shape_second_derivatives[qpoint * n_shape_functions + shape_nr];
545}
546
547template <int dim, int spacedim>
548inline const Tensor<3, dim> &
550 const unsigned int qpoint,
551 const unsigned int shape_nr) const
552{
553 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
554 shape_third_derivatives.size());
555 return shape_third_derivatives[qpoint * n_shape_functions + shape_nr];
556}
557
558
559template <int dim, int spacedim>
560inline Tensor<3, dim> &
562 const unsigned int qpoint,
563 const unsigned int shape_nr)
564{
565 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
566 shape_third_derivatives.size());
567 return shape_third_derivatives[qpoint * n_shape_functions + shape_nr];
568}
569
570
571template <int dim, int spacedim>
572inline const Tensor<4, dim> &
574 const unsigned int qpoint,
575 const unsigned int shape_nr) const
576{
577 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
578 shape_fourth_derivatives.size());
579 return shape_fourth_derivatives[qpoint * n_shape_functions + shape_nr];
580}
581
582
583template <int dim, int spacedim>
584inline Tensor<4, dim> &
586 const unsigned int qpoint,
587 const unsigned int shape_nr)
588{
589 AssertIndexRange(qpoint * n_shape_functions + shape_nr,
590 shape_fourth_derivatives.size());
591 return shape_fourth_derivatives[qpoint * n_shape_functions + shape_nr];
592}
593
594
595
596template <int dim, int spacedim>
597inline bool
599{
600 return true;
601}
602
603
604
605#endif // DOXYGEN
606
607/* -------------- declaration of explicit specializations ------------- */
608
609
611
612#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)
std::vector< Point< spacedim > > mapping_support_points
Definition mapping_fe.h:373
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:367
std::vector< Tensor< 3, dim > > shape_third_derivatives
Definition mapping_fe.h:309
std::vector< Tensor< 2, dim > > shape_second_derivatives
Definition mapping_fe.h:301
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:317
Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr)
Triangulation< dim, spacedim >::cell_iterator cell_of_current_support_points
Definition mapping_fe.h:379
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)
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
Definition mapping_fe.h:327
std::vector< double > quadrature_weights
Definition mapping_fe.h:390
const FiniteElement< dim, spacedim > & fe
Definition mapping_fe.h:332
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
Definition mapping_fe.h:362
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
Definition mapping_fe.h:353
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:286
std::vector< Tensor< 1, dim > > shape_derivatives
Definition mapping_fe.h:293
virtual std::size_t memory_consumption() const override
Definition mapping_fe.cc:60
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
const unsigned int n_shape_functions
Definition mapping_fe.h:342
std::vector< double > volume_elements
Definition mapping_fe.h:385
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
Definition mapping_fe.cc:81
const unsigned int polynomial_degree
Definition mapping_fe.h:337
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
const unsigned int polynomial_degree
Definition mapping_fe.h:458
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) 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 UpdateFlags requires_update_flags(const UpdateFlags update_flags) 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
Table< 2, double > mapping_support_point_weights
Definition mapping_fe.h:468
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) 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 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_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
unsigned int get_degree() const
virtual bool is_compatible_with(const ReferenceCell &reference_cell) 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 bool preserves_vertex_locations() const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
const std::unique_ptr< FiniteElement< dim, spacedim > > fe
Definition mapping_fe.h:452
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
#define AssertIndexRange(index, range)
UpdateFlags
MappingKind
Definition mapping.h:79