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_q_hierarchical.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2002 - 2023 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_q_hierarchical_h
16#define dealii_fe_q_hierarchical_h
17
18#include <deal.II/base/config.h>
19
21
22#include <deal.II/fe/fe_poly.h>
23
25
539template <int dim>
540class FE_Q_Hierarchical : public FE_Poly<dim>
541{
542public:
546 FE_Q_Hierarchical(const unsigned int p);
547
553 virtual std::string
554 get_name() const override;
555
556 virtual std::unique_ptr<FiniteElement<dim, dim>>
557 clone() const override;
558
563 virtual bool
564 has_support_on_face(const unsigned int shape_index,
565 const unsigned int face_index) const override;
566
580 virtual bool
581 hp_constraints_are_implemented() const override;
582
587 virtual void
589 FullMatrix<double> &matrix) const override;
590
594 virtual const FullMatrix<double> &
596 const unsigned int child,
597 const RefinementCase<dim> &refinement_case =
599
615 virtual std::vector<std::pair<unsigned int, unsigned int>>
616 hp_vertex_dof_identities(const FiniteElement<dim> &fe_other) const override;
617
621 virtual std::vector<std::pair<unsigned int, unsigned int>>
622 hp_line_dof_identities(const FiniteElement<dim> &fe_other) const override;
623
627 virtual std::vector<std::pair<unsigned int, unsigned int>>
629 const unsigned int face_no = 0) const override;
630
636 const unsigned int codim = 0) const override final;
637
651 virtual void
653 FullMatrix<double> &matrix,
654 const unsigned int face_no = 0) const override;
655
667 virtual void
669 const FiniteElement<dim> &source,
670 const unsigned int subface,
671 FullMatrix<double> &matrix,
672 const unsigned int face_no = 0) const override;
673
682 virtual std::size_t
683 memory_consumption() const override;
684
690 std::vector<unsigned int>
691 get_embedding_dofs(const unsigned int sub_degree) const;
692
698 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
699 get_constant_modes() const override;
700
701private:
708 static std::vector<unsigned int>
709 get_dpo_vector(const unsigned int degree);
710
738 static std::vector<unsigned int>
740
744 static std::vector<unsigned int>
746
751 void
752 build_dofs_cell(std::vector<FullMatrix<double>> &dofs_cell,
753 std::vector<FullMatrix<double>> &dofs_subcell) const;
754
759 void
760 initialize_constraints(const std::vector<FullMatrix<double>> &dofs_subcell);
761
765 void
767 const std::vector<FullMatrix<double>> &dofs_cell,
768 const std::vector<FullMatrix<double>> &dofs_subcell);
769
774 void
776
781 void
783
787 const std::vector<unsigned int> face_renumber;
788
789 // Allow access from other dimensions. We need this since we want to call
790 // the functions @p get_dpo_vector and @p
791 // lexicographic_to_hierarchic_numbering for the faces of the finite element
792 // of dimension dim+1.
793 template <int dim1>
794 friend class FE_Q_Hierarchical;
795};
796
799/* -------------- declaration of explicit specializations ------------- */
800
801template <>
802void
804
805template <>
806bool
808 const unsigned int) const;
809
810template <>
811std::vector<unsigned int>
813 const unsigned int);
814
816
817#endif
const std::vector< unsigned int > face_renumber
void build_dofs_cell(std::vector< FullMatrix< double > > &dofs_cell, std::vector< FullMatrix< double > > &dofs_subcell) const
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual void get_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix) const override
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
void initialize_embedding_and_restriction(const std::vector< FullMatrix< double > > &dofs_cell, const std::vector< FullMatrix< double > > &dofs_subcell)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
friend class FE_Q_Hierarchical
std::vector< unsigned int > get_embedding_dofs(const unsigned int sub_degree) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim > &fe_other, const unsigned int face_no=0) const override
void initialize_generalized_support_points()
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim > &fe_other, const unsigned int codim=0) const override final
static std::vector< unsigned int > hierarchic_to_fe_q_hierarchical_numbering(const FiniteElementData< dim > &fe)
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
virtual void get_subface_interpolation_matrix(const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim > &fe_other) const override
static std::vector< unsigned int > face_fe_q_hierarchical_to_hierarchic_numbering(const unsigned int degree)
void initialize_generalized_face_support_points()
void initialize_constraints(const std::vector< FullMatrix< double > > &dofs_subcell)
virtual std::string get_name() const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::size_t memory_consumption() const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual bool hp_constraints_are_implemented() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim > &fe_other) const override
const unsigned int degree
Definition fe_data.h:452
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504