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_nedelec.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_nedelec_h
16#define dealii_fe_nedelec_h
17
18#include <deal.II/base/config.h>
19
21#include <deal.II/base/mutex.h>
24#include <deal.II/base/table.h>
25#include <deal.II/base/tensor.h>
27
28#include <deal.II/fe/fe.h>
30
31#include <vector>
32
34
399template <int dim>
400class FE_Nedelec : public FE_PolyTensor<dim>
401{
402public:
417 FE_Nedelec(const unsigned int order);
418
424 virtual std::string
425 get_name() const override;
426
427
432 virtual bool
433 has_support_on_face(const unsigned int shape_index,
434 const unsigned int face_index) const override;
435
444 virtual bool
445 hp_constraints_are_implemented() const override;
446
452 const unsigned int codim = 0) const override final;
453
469 virtual std::vector<std::pair<unsigned int, unsigned int>>
470 hp_vertex_dof_identities(const FiniteElement<dim> &fe_other) const override;
471
476 virtual std::vector<std::pair<unsigned int, unsigned int>>
477 hp_line_dof_identities(const FiniteElement<dim> &fe_other) const override;
478
483 virtual std::vector<std::pair<unsigned int, unsigned int>>
485 const unsigned int face_no = 0) const override;
486
498 virtual void
500 FullMatrix<double> &matrix,
501 const unsigned int face_no = 0) const override;
502
514 virtual void
516 const FiniteElement<dim> &source,
517 const unsigned int subface,
518 FullMatrix<double> &matrix,
519 const unsigned int face_no = 0) const override;
520
535 virtual const FullMatrix<double> &
537 const unsigned int child,
538 const RefinementCase<dim> &refinement_case =
540
561 virtual const FullMatrix<double> &
563 const unsigned int child,
564 const RefinementCase<dim> &refinement_case =
566
567 // documentation inherited from the base class
568 virtual void
570 const std::vector<Vector<double>> &support_point_values,
571 std::vector<double> &nodal_values) const override;
572
576 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
577 get_constant_modes() const override;
578
579 virtual std::size_t
580 memory_consumption() const override;
581
582 virtual std::unique_ptr<FiniteElement<dim, dim>>
583 clone() const override;
584
594 std::vector<unsigned int>
595 get_embedding_dofs(const unsigned int sub_degree) const;
596
597private:
608 static std::vector<unsigned int>
609 get_dpo_vector(const unsigned int degree, bool dg = false);
610
616 void
617 initialize_support_points(const unsigned int order);
618
624 void
626
636
643
651 void
653
654 // Allow access from other dimensions.
655 template <int dim1>
656 friend class FE_Nedelec;
657};
658
659/* -------------- declaration of explicit specializations ------------- */
660
661#ifndef DOXYGEN
662
663template <>
664void
666
667#endif // DOXYGEN
668
672
673#endif
void initialize_quad_dof_index_permutation_and_sign_change()
std::vector< unsigned int > get_embedding_dofs(const unsigned int sub_degree) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim > &fe_other) const override
virtual bool hp_constraints_are_implemented() 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
void initialize_restriction()
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
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree, bool dg=false)
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) 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::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
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
Threads::Mutex restriction_matrix_mutex
Definition fe_nedelec.h:641
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) 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 std::string get_name() const override
friend class FE_Nedelec
Definition fe_nedelec.h:656
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim > &fe_other, const unsigned int codim=0) const override final
Threads::Mutex prolongation_matrix_mutex
Definition fe_nedelec.h:642
Table< 2, double > boundary_weights
Definition fe_nedelec.h:635
void initialize_support_points(const unsigned int order)
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_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