deal.II version GIT relicensing-2291-g668cd86249 2024-12-24 11:30: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_field.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_field_h
16#define dealii_mapping_fe_field_h
17
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/mutex.h>
23
25
26#include <deal.II/fe/fe.h>
28#include <deal.II/fe/mapping.h>
29
30#include <deal.II/lac/vector.h>
31
32#include <array>
33
34
36
37
78template <int dim, int spacedim = dim, typename VectorType = Vector<double>>
79class MappingFEField : public Mapping<dim, spacedim>
80{
81public:
115 const VectorType &euler_vector,
116 const ComponentMask &mask = {});
117
128 const std::vector<VectorType> &euler_vector,
129 const ComponentMask &mask = {});
130
140 const ComponentMask &mask = {});
141
146
151 virtual std::unique_ptr<Mapping<dim, spacedim>>
152 clone() const override;
153
159 virtual bool
160 preserves_vertex_locations() const override;
161
162 virtual bool
163 is_compatible_with(const ReferenceCell &reference_cell) const override;
164
172 virtual boost::container::small_vector<Point<spacedim>,
175 const override;
176
182 // for documentation, see the Mapping base class
183 virtual Point<spacedim>
186 const Point<dim> &p) const override;
187
188 // for documentation, see the Mapping base class
189 virtual Point<dim>
192 const Point<spacedim> &p) const override;
193
203 // for documentation, see the Mapping base class
204 virtual void
205 transform(const ArrayView<const Tensor<1, dim>> &input,
206 const MappingKind kind,
208 const ArrayView<Tensor<1, spacedim>> &output) const override;
209
210 // for documentation, see the Mapping base class
211 virtual void
213 const MappingKind kind,
215 const ArrayView<Tensor<2, spacedim>> &output) const override;
216
217 // for documentation, see the Mapping base class
218 virtual void
219 transform(const ArrayView<const Tensor<2, dim>> &input,
220 const MappingKind kind,
222 const ArrayView<Tensor<2, spacedim>> &output) const override;
223
224 // for documentation, see the Mapping base class
225 virtual void
227 const MappingKind kind,
229 const ArrayView<Tensor<3, spacedim>> &output) const override;
230
231 // for documentation, see the Mapping base class
232 virtual void
233 transform(const ArrayView<const Tensor<3, dim>> &input,
234 const MappingKind kind,
236 const ArrayView<Tensor<3, spacedim>> &output) const override;
237
246 unsigned int
247 get_degree() const;
248
254 get_component_mask() const;
255
260
261private:
267 // documentation can be found in Mapping::requires_update_flags()
268 virtual UpdateFlags
269 requires_update_flags(const UpdateFlags update_flags) const override;
270
271public:
283 class InternalData : public Mapping<dim, spacedim>::InternalDataBase
284 {
285 public:
290 const ComponentMask &mask);
291
292 // Documentation see Mapping::InternalDataBase.
293 virtual void
294 reinit(const UpdateFlags update_flags,
295 const Quadrature<dim> &quadrature) override;
296
301 const double &
302 shape(const unsigned int qpoint, const unsigned int shape_nr) const;
303
307 double &
308 shape(const unsigned int qpoint, const unsigned int shape_nr);
309
313 const Tensor<1, dim> &
314 derivative(const unsigned int qpoint, const unsigned int shape_nr) const;
315
320 derivative(const unsigned int qpoint, const unsigned int shape_nr);
321
325 const Tensor<2, dim> &
326 second_derivative(const unsigned int qpoint,
327 const unsigned int shape_nr) const;
328
333 second_derivative(const unsigned int qpoint, const unsigned int shape_nr);
334
338 const Tensor<3, dim> &
339 third_derivative(const unsigned int qpoint,
340 const unsigned int shape_nr) const;
341
346 third_derivative(const unsigned int qpoint, const unsigned int shape_nr);
347
351 const Tensor<4, dim> &
352 fourth_derivative(const unsigned int qpoint,
353 const unsigned int shape_nr) const;
354
359 fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr);
360
364 virtual std::size_t
365 memory_consumption() const override;
366
371
377 std::vector<double> shape_values;
378
384 std::vector<Tensor<1, dim>> shape_derivatives;
385
392 std::vector<Tensor<2, dim>> shape_second_derivatives;
393
400 std::vector<Tensor<3, dim>> shape_third_derivatives;
401
408 std::vector<Tensor<4, dim>> shape_fourth_derivatives;
409
423 std::array<std::vector<Tensor<1, dim>>,
426
433 unsigned int n_shape_functions;
434
448
458 mutable std::vector<DerivativeForm<1, dim, spacedim>> covariant;
459
467 mutable std::vector<DerivativeForm<1, dim, spacedim>> contravariant;
468
473 mutable std::vector<double> volume_elements;
474
478 mutable std::vector<double> quadrature_weights;
479
483 mutable std::vector<std::vector<Tensor<1, spacedim>>> aux;
484
488 mutable std::vector<types::global_dof_index> local_dof_indices;
489
493 mutable std::vector<double> local_dof_values;
494 };
495
496protected:
497 // documentation can be found in Mapping::get_data()
498 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
499 get_data(const UpdateFlags, const Quadrature<dim> &quadrature) const override;
500
501 using Mapping<dim, spacedim>::get_face_data;
502
503 // documentation can be found in Mapping::get_face_data()
504 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
505 get_face_data(const UpdateFlags flags,
506 const hp::QCollection<dim - 1> &quadrature) const override;
507
508 // documentation can be found in Mapping::get_subface_data()
509 virtual std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
510 get_subface_data(const UpdateFlags flags,
511 const Quadrature<dim - 1> &quadrature) const override;
512
513 // documentation can be found in Mapping::fill_fe_values()
517 const CellSimilarity::Similarity cell_similarity,
518 const Quadrature<dim> &quadrature,
519 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
521 &output_data) const override;
522
523 using Mapping<dim, spacedim>::fill_fe_face_values;
524
525 // documentation can be found in Mapping::fill_fe_face_values()
526 virtual void
529 const unsigned int face_no,
530 const hp::QCollection<dim - 1> &quadrature,
531 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
533 &output_data) const override;
534
535 // documentation can be found in Mapping::fill_fe_subface_values()
536 virtual void
539 const unsigned int face_no,
540 const unsigned int subface_no,
541 const Quadrature<dim - 1> &quadrature,
542 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
544 &output_data) const override;
545
546 virtual void
550 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
552 &output_data) const override;
553
563
568 const bool uses_level_dofs;
569
573 std::vector<ObserverPointer<const VectorType,
576
583
584private:
602
620 const Point<spacedim> &p,
621 const Point<dim> &starting_guess,
622 InternalData &mdata) const;
623
627 void
631 &data) const;
632
633 /*
634 * Which components to use for the mapping.
635 */
637
647 std::vector<unsigned int> fe_to_real;
648
654
659
660 void
661 compute_face_data(const unsigned int n_original_q_points,
662 InternalData &data) const;
663
664 // Declare other MappingFEField classes friends.
665 template <int, int, class>
666 friend class MappingFEField;
667};
668
671/* -------------- declaration of explicit specializations ------------- */
672
673#ifndef DOXYGEN
674
675#endif // DOXYGEN
676
678
679#endif
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
std::vector< double > quadrature_weights
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
std::vector< std::vector< Tensor< 1, spacedim > > > aux
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
ObserverPointer< const FiniteElement< dim, spacedim > > fe
std::vector< Tensor< 3, dim > > shape_third_derivatives
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< types::global_dof_index > local_dof_indices
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
virtual std::size_t memory_consumption() const override
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< Tensor< 1, dim > > shape_derivatives
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< double > local_dof_values
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< double > shape_values
std::vector< Tensor< 2, dim > > shape_second_derivatives
std::vector< double > volume_elements
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
std::vector< unsigned int > fe_to_real
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
const ComponentMask fe_mask
void update_internal_dofs(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim, VectorType >::InternalData &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 Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
virtual bool preserves_vertex_locations() const override
Point< dim > do_transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &starting_guess, InternalData &mdata) const
void compute_face_data(const unsigned int n_original_q_points, InternalData &data) const
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
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
friend class MappingFEField
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
std::vector< ObserverPointer< const VectorType, MappingFEField< dim, spacedim, VectorType > > > euler_vector
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
ReferenceCell reference_cell
unsigned int get_degree() const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Point< spacedim > do_transform_unit_to_real_cell(const InternalData &mdata) const
Threads::Mutex fe_values_mutex
const bool uses_level_dofs
ComponentMask get_component_mask() 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
ObserverPointer< const DoFHandler< dim, spacedim >, MappingFEField< dim, spacedim, VectorType > > euler_dof_handler
FEValues< dim, spacedim > fe_values
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() 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:318
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define DeclException0(Exception0)
Definition exceptions.h:466
static ::ExceptionBase & ExcInactiveCell()
UpdateFlags
MappingKind
Definition mapping.h:79
std::vector< index_type > data
Definition mpi.cc:735