Reference documentation for deal.II version 9.5.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_sz.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2015 - 2023 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_fe_nedelec_sz_h
17#define dealii_fe_nedelec_sz_h
18
19#include <deal.II/base/config.h>
20
24
25#include <deal.II/fe/fe.h>
27#include <deal.II/fe/mapping.h>
28
30
74template <int dim, int spacedim = dim>
75class FE_NedelecSZ : public FiniteElement<dim, dim>
76{
77public:
78 static_assert(dim == spacedim,
79 "FE_NedelecSZ is only implemented for dim==spacedim!");
80
95 FE_NedelecSZ(const unsigned int order);
96
97 virtual UpdateFlags
98 requires_update_flags(const UpdateFlags update_flags) const override;
99
100 virtual std::string
101 get_name() const override;
102
103 virtual std::unique_ptr<FiniteElement<dim, dim>>
104 clone() const override;
105
110 virtual double
111 shape_value(const unsigned int i, const Point<dim> &p) const override;
112
116 virtual double
117 shape_value_component(const unsigned int i,
118 const Point<dim> & p,
119 const unsigned int component) const override;
120
125 virtual Tensor<1, dim>
126 shape_grad(const unsigned int i, const Point<dim> &p) const override;
127
131 virtual Tensor<1, dim>
132 shape_grad_component(const unsigned int i,
133 const Point<dim> & p,
134 const unsigned int component) const override;
135
140 virtual Tensor<2, dim>
141 shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
142
146 virtual Tensor<2, dim>
147 shape_grad_grad_component(const unsigned int i,
148 const Point<dim> & p,
149 const unsigned int component) const override;
150
151protected:
157
158 virtual std::unique_ptr<
159 typename ::FiniteElement<dim, spacedim>::InternalDataBase>
160 get_data(
161 const UpdateFlags update_flags,
162 const Mapping<dim, spacedim> &mapping,
163 const Quadrature<dim> & quadrature,
165 spacedim>
166 &output_data) const override;
167
173 virtual void
175 const typename Triangulation<dim, dim>::cell_iterator &cell,
176 const CellSimilarity::Similarity cell_similarity,
177 const Quadrature<dim> & quadrature,
178 const Mapping<dim, dim> & mapping,
179 const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
181 & mapping_data,
182 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
184 &data) const override;
185
186 using FiniteElement<dim, spacedim>::fill_fe_face_values;
187
193 virtual void
195 const typename Triangulation<dim, dim>::cell_iterator &cell,
196 const unsigned int face_no,
197 const hp::QCollection<dim - 1> & quadrature,
198 const Mapping<dim, dim> & mapping,
199 const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
201 & mapping_data,
202 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
204 &data) const override;
205
209 virtual void
211 const typename Triangulation<dim, dim>::cell_iterator &cell,
212 const unsigned int face_no,
213 const unsigned int sub_no,
214 const Quadrature<dim - 1> & quadrature,
215 const Mapping<dim, dim> & mapping,
216 const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
218 & mapping_data,
219 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
221 &data) const override;
222
243 {
244 public:
251 mutable std::vector<std::vector<Tensor<1, dim>>> shape_values;
252
260 mutable std::vector<std::vector<DerivativeForm<1, dim, dim>>> shape_grads;
261
269 mutable std::vector<std::vector<DerivativeForm<2, dim, dim>>>
271
287 std::vector<std::vector<std::vector<double>>> sigma_imj_values;
288
304 std::vector<std::vector<std::vector<double>>> sigma_imj_grads;
305
317 std::vector<std::vector<double>> edge_sigma_values;
318
330 std::vector<std::vector<double>> edge_sigma_grads;
331
347 std::vector<std::vector<double>> edge_lambda_values;
348
358 std::vector<std::vector<double>> edge_lambda_grads_2d;
359
369 std::vector<std::vector<std::vector<double>>> edge_lambda_grads_3d;
370
381 std::vector<std::vector<std::vector<double>>> edge_lambda_gradgrads_3d;
382
399 std::vector<std::vector<double>> face_lambda_values;
400
409 std::vector<std::vector<double>> face_lambda_grads;
410 };
411
412private:
421 static std::vector<unsigned int>
422 get_dpo_vector(const unsigned int degree);
423
427 std::vector<Polynomials::Polynomial<double>> IntegratedLegendrePolynomials;
428
433 void
434 create_polynomials(const unsigned int degree);
435
439 unsigned int
440 compute_num_dofs(const unsigned int degree) const;
441
446 void
448 const Quadrature<dim> &quadrature,
449 const InternalData & fedata) const;
450
455 void
457 const Quadrature<dim> &quadrature,
458 const InternalData & fedata) const;
459};
460
461
462
466
467#endif
std::vector< std::vector< double > > edge_lambda_grads_2d
std::vector< std::vector< Tensor< 1, dim > > > shape_values
std::vector< std::vector< std::vector< double > > > sigma_imj_values
std::vector< std::vector< double > > edge_sigma_grads
std::vector< std::vector< double > > edge_sigma_values
std::vector< std::vector< DerivativeForm< 2, dim, dim > > > shape_hessians
std::vector< std::vector< std::vector< double > > > sigma_imj_grads
std::vector< std::vector< double > > face_lambda_values
std::vector< std::vector< DerivativeForm< 1, dim, dim > > > shape_grads
std::vector< std::vector< std::vector< double > > > edge_lambda_grads_3d
std::vector< std::vector< double > > edge_lambda_values
std::vector< std::vector< std::vector< double > > > edge_lambda_gradgrads_3d
std::vector< std::vector< double > > face_lambda_grads
unsigned int compute_num_dofs(const unsigned int degree) const
virtual void fill_fe_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
void create_polynomials(const unsigned int degree)
virtual std::unique_ptr< typename ::FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual void fill_fe_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, dim > &mapping, const typename Mapping< dim, dim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, dim > &mapping_data, const typename FiniteElement< dim, dim >::InternalDataBase &fedata, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, dim > &data) const override
virtual std::string get_name() const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
std::vector< Polynomials::Polynomial< double > > IntegratedLegendrePolynomials
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
void fill_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
void fill_edge_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
MappingKind mapping_kind
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
const unsigned int degree
Definition fe_data.h:453
friend class InternalDataBase
Definition fe.h:3061
Abstract base class for mapping classes.
Definition mapping.h:317
Definition point.h:112
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
UpdateFlags
MappingKind
Definition mapping.h:78