Reference documentation for deal.II version 9.6.0
\(\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
fe_face.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2009 - 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_fe_face_h
16#define dealii_fe_face_h
17
18#include <deal.II/base/config.h>
19
22
24
26
27
53template <int dim, int spacedim = dim>
55 : public FE_PolyFace<TensorProductPolynomials<dim - 1>, dim, spacedim>
56{
57public:
63 FE_FaceQ(const unsigned int p);
64
65 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
66 clone() const override;
67
73 virtual std::string
74 get_name() const override;
75
83 virtual void
85 const std::vector<Vector<double>> &support_point_values,
86 std::vector<double> &nodal_values) const override;
87
96 virtual void
98 FullMatrix<double> &matrix,
99 const unsigned int face_no = 0) const override;
100
109 virtual void
111 const FiniteElement<dim, spacedim> &source,
112 const unsigned int subface,
113 FullMatrix<double> &matrix,
114 const unsigned int face_no = 0) const override;
115
120 virtual bool
121 has_support_on_face(const unsigned int shape_index,
122 const unsigned int face_index) const override;
123
146 virtual std::vector<std::pair<unsigned int, unsigned int>>
148 const FiniteElement<dim, spacedim> &fe_other) const override;
149
156 virtual std::vector<std::pair<unsigned int, unsigned int>>
158 const FiniteElement<dim, spacedim> &fe_other) const override;
159
166 virtual std::vector<std::pair<unsigned int, unsigned int>>
168 const unsigned int face_no = 0) const override;
169
174 virtual bool
175 hp_constraints_are_implemented() const override;
176
182 const unsigned int codim = 0) const override final;
183
192 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
193 get_constant_modes() const override;
194
195private:
199 static std::vector<unsigned int>
200 get_dpo_vector(const unsigned int deg);
201};
202
203
204
215template <int spacedim>
216class FE_FaceQ<1, spacedim> : public FiniteElement<1, spacedim>
217{
218public:
222 FE_FaceQ(const unsigned int p);
223
224 virtual std::unique_ptr<FiniteElement<1, spacedim>>
225 clone() const override;
226
232 virtual std::string
233 get_name() const override;
234
235 // for documentation, see the FiniteElement base class
236 virtual UpdateFlags
237 requires_update_flags(const UpdateFlags update_flags) const override;
238
247 virtual void
249 FullMatrix<double> &matrix,
250 const unsigned int face_no = 0) const override;
251
260 virtual void
262 const FiniteElement<1, spacedim> &source,
263 const unsigned int subface,
264 FullMatrix<double> &matrix,
265 const unsigned int face_no = 0) const override;
266
271 virtual bool
272 has_support_on_face(const unsigned int shape_index,
273 const unsigned int face_index) const override;
274
279 virtual bool
280 hp_constraints_are_implemented() const override;
281
299 virtual std::vector<std::pair<unsigned int, unsigned int>>
301 const FiniteElement<1, spacedim> &fe_other) const override;
302
309 virtual std::vector<std::pair<unsigned int, unsigned int>>
311 const FiniteElement<1, spacedim> &fe_other) const override;
312
319 virtual std::vector<std::pair<unsigned int, unsigned int>>
321 const unsigned int face_no = 0) const override;
322
327 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
328 get_constant_modes() const override;
329
330protected:
331 /*
332 * NOTE: The following functions have their definitions inlined into the class
333 * declaration because we otherwise run into a compiler error with MS Visual
334 * Studio.
335 */
336
337
338 virtual std::unique_ptr<typename FiniteElement<1, spacedim>::InternalDataBase>
340 const UpdateFlags /*update_flags*/,
341 const Mapping<1, spacedim> & /*mapping*/,
342 const Quadrature<1> & /*quadrature*/,
344 spacedim>
345 & /*output_data*/) const override
346 {
347 return std::make_unique<
349 }
350
351 using FiniteElement<1, spacedim>::get_face_data;
352
353 std::unique_ptr<typename FiniteElement<1, spacedim>::InternalDataBase>
355 const UpdateFlags update_flags,
356 const Mapping<1, spacedim> & /*mapping*/,
357 const hp::QCollection<0> &quadrature,
359 spacedim>
360 & /*output_data*/) const override
361 {
362 AssertDimension(quadrature.size(), 1);
363
364 // generate a new data object and initialize some fields
365 auto data_ptr =
366 std::make_unique<typename FiniteElement<1, spacedim>::InternalDataBase>();
367 data_ptr->update_each = requires_update_flags(update_flags);
368
369 const unsigned int n_q_points = quadrature[0].size();
370 AssertDimension(n_q_points, 1);
371 (void)n_q_points;
372
373 // No derivatives of this element are implemented.
374 if (data_ptr->update_each & update_gradients ||
375 data_ptr->update_each & update_hessians)
376 {
378 }
379
380 return data_ptr;
381 }
382
383 std::unique_ptr<typename FiniteElement<1, spacedim>::InternalDataBase>
385 const UpdateFlags update_flags,
386 const Mapping<1, spacedim> &mapping,
387 const Quadrature<0> &quadrature,
389 spacedim>
390 &output_data) const override
391 {
392 return get_face_data(update_flags,
393 mapping,
394 hp::QCollection<0>(quadrature),
395 output_data);
396 }
397
398 virtual void
401 const CellSimilarity::Similarity cell_similarity,
402 const Quadrature<1> &quadrature,
403 const Mapping<1, spacedim> &mapping,
404 const typename Mapping<1, spacedim>::InternalDataBase &mapping_internal,
406 &mapping_data,
407 const typename FiniteElement<1, spacedim>::InternalDataBase &fe_internal,
409 spacedim>
410 &output_data) const override;
411
412 using FiniteElement<1, spacedim>::fill_fe_face_values;
413
414 virtual void
417 const unsigned int face_no,
418 const hp::QCollection<0> &quadrature,
419 const Mapping<1, spacedim> &mapping,
420 const typename Mapping<1, spacedim>::InternalDataBase &mapping_internal,
422 &mapping_data,
423 const typename FiniteElement<1, spacedim>::InternalDataBase &fe_internal,
425 spacedim>
426 &output_data) const override;
427
428 virtual void
431 const unsigned int face_no,
432 const unsigned int sub_no,
433 const Quadrature<0> &quadrature,
434 const Mapping<1, spacedim> &mapping,
435 const typename Mapping<1, spacedim>::InternalDataBase &mapping_internal,
437 &mapping_data,
438 const typename FiniteElement<1, spacedim>::InternalDataBase &fe_internal,
440 spacedim>
441 &output_data) const override;
442
443private:
447 static std::vector<unsigned int>
448 get_dpo_vector(const unsigned int deg);
449};
450
451
452
475template <int dim, int spacedim = dim>
476class FE_FaceP : public FE_PolyFace<PolynomialSpace<dim - 1>, dim, spacedim>
477{
478public:
484 FE_FaceP(unsigned int p);
485
486 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
487 clone() const override;
488
494 virtual std::string
495 get_name() const override;
496
505 virtual void
507 FullMatrix<double> &matrix,
508 const unsigned int face_no = 0) const override;
509
518 virtual void
520 const FiniteElement<dim, spacedim> &source,
521 const unsigned int subface,
522 FullMatrix<double> &matrix,
523 const unsigned int face_no = 0) const override;
524
529 virtual bool
530 has_support_on_face(const unsigned int shape_index,
531 const unsigned int face_index) const override;
532
537 virtual bool
538 hp_constraints_are_implemented() const override;
539
545 const unsigned int codim = 0) const override final;
546
553 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
554 get_constant_modes() const override;
555
556private:
560 static std::vector<unsigned int>
561 get_dpo_vector(const unsigned int deg);
562};
563
564
565
570template <int spacedim>
571class FE_FaceP<1, spacedim> : public FE_FaceQ<1, spacedim>
572{
573public:
577 FE_FaceP(const unsigned int p);
578
582 std::string
583 get_name() const override;
584};
585
586
588
589#endif
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition fe_face.cc:993
virtual bool hp_constraints_are_implemented() const override
Definition fe_face.cc:823
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition fe_face.cc:796
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
Definition fe_face.cc:832
virtual std::string get_name() const override
Definition fe_face.cc:780
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
Definition fe_face.cc:807
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
Definition fe_face.cc:886
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_face.cc:771
FE_FaceP(unsigned int p)
Definition fe_face.cc:756
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
Definition fe_face.cc:871
std::unique_ptr< typename FiniteElement< 1, spacedim >::InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< 1, spacedim > &, const hp::QCollection< 0 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 1, spacedim > &) const override
Definition fe_face.h:354
virtual std::unique_ptr< typename FiniteElement< 1, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Mapping< 1, spacedim > &, const Quadrature< 1 > &, ::internal::FEValuesImplementation::FiniteElementRelatedData< 1, spacedim > &) const override
Definition fe_face.h:339
std::unique_ptr< typename FiniteElement< 1, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< 1, spacedim > &mapping, const Quadrature< 0 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 1, spacedim > &output_data) const override
Definition fe_face.h:384
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
Definition fe_face.cc:447
virtual std::string get_name() const override
Definition fe_face.cc:130
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
Definition fe_face.cc:265
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const override
Definition fe_face.cc:495
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition fe_face.cc:484
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
Definition fe_face.cc:146
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition fe_face.cc:298
FE_FaceQ(const unsigned int p)
Definition fe_face.cc:50
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition fe_face.cc:121
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
Definition fe_face.cc:161
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition fe_face.cc:254
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
Definition fe_face.cc:287
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const override
Definition fe_face.cc:368
virtual bool hp_constraints_are_implemented() const override
Definition fe_face.cc:278
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 Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const hp::QCollection< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
Abstract base class for mapping classes.
Definition mapping.h:318
unsigned int size() const
Definition collection.h:308
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
#define AssertDimension(dim1, dim2)
UpdateFlags
@ update_hessians
Second derivatives of shape functions.
@ update_gradients
Shape function gradients.
#define DEAL_II_NOT_IMPLEMENTED()