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_sz.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2018 - 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_nedelec_sz_h
16#define dealii_fe_nedelec_sz_h
17
18#include <deal.II/base/config.h>
19
24
25#include <deal.II/fe/fe.h>
27#include <deal.II/fe/mapping.h>
28
30
83template <int dim, int spacedim = dim>
84class FE_NedelecSZ : public FiniteElement<dim, dim>
85{
86public:
87 static_assert(dim == spacedim,
88 "FE_NedelecSZ is only implemented for dim==spacedim!");
89
104 FE_NedelecSZ(const unsigned int order);
105
106 virtual UpdateFlags
107 requires_update_flags(const UpdateFlags update_flags) const override;
108
109 virtual std::string
110 get_name() const override;
111
112 virtual std::unique_ptr<FiniteElement<dim, dim>>
113 clone() const override;
114
123 virtual double
124 shape_value(const unsigned int i, const Point<dim> &p) const override;
125
132 virtual double
133 shape_value_component(const unsigned int i,
134 const Point<dim> &p,
135 const unsigned int component) const override;
136
145 virtual Tensor<1, dim>
146 shape_grad(const unsigned int i, const Point<dim> &p) const override;
147
154 virtual Tensor<1, dim>
155 shape_grad_component(const unsigned int i,
156 const Point<dim> &p,
157 const unsigned int component) const override;
158
167 virtual Tensor<2, dim>
168 shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
169
176 virtual Tensor<2, dim>
177 shape_grad_grad_component(const unsigned int i,
178 const Point<dim> &p,
179 const unsigned int component) const override;
180
201 virtual const FullMatrix<double> &
203 const unsigned int child,
204 const RefinementCase<dim> &refinement_case =
206
207protected:
213
218 void
219 evaluate(const std::vector<Point<dim>> &p_list,
220 const UpdateFlags update_flags,
221 std::unique_ptr<
222 typename ::FiniteElement<dim, spacedim>::InternalDataBase>
223 &data_ptr) const;
224
225 virtual std::unique_ptr<
226 typename ::FiniteElement<dim, spacedim>::InternalDataBase>
227 get_data(
228 const UpdateFlags update_flags,
229 const Mapping<dim, spacedim> &mapping,
230 const Quadrature<dim> &quadrature,
232 spacedim>
233 &output_data) const override;
234
240 virtual void
242 const typename Triangulation<dim, dim>::cell_iterator &cell,
243 const CellSimilarity::Similarity cell_similarity,
244 const Quadrature<dim> &quadrature,
245 const Mapping<dim, dim> &mapping,
246 const typename Mapping<dim, dim>::InternalDataBase &mapping_internal,
248 &mapping_data,
249 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
251 &data) const override;
252
253 using FiniteElement<dim, spacedim>::fill_fe_face_values;
254
260 virtual void
262 const typename Triangulation<dim, dim>::cell_iterator &cell,
263 const unsigned int face_no,
264 const hp::QCollection<dim - 1> &quadrature,
265 const Mapping<dim, dim> &mapping,
266 const typename Mapping<dim, dim>::InternalDataBase &mapping_internal,
268 &mapping_data,
269 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
271 &data) const override;
272
276 virtual void
278 const typename Triangulation<dim, dim>::cell_iterator &cell,
279 const unsigned int face_no,
280 const unsigned int sub_no,
281 const Quadrature<dim - 1> &quadrature,
282 const Mapping<dim, dim> &mapping,
283 const typename Mapping<dim, dim>::InternalDataBase &mapping_internal,
285 &mapping_data,
286 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
288 &data) const override;
289
310 {
311 public:
318 mutable std::vector<std::vector<Tensor<1, dim>>> shape_values;
319
327 mutable std::vector<std::vector<DerivativeForm<1, dim, dim>>> shape_grads;
328
336 mutable std::vector<std::vector<DerivativeForm<2, dim, dim>>>
338
354 std::vector<std::vector<std::vector<double>>> sigma_imj_values;
355
371 std::vector<std::vector<std::vector<double>>> sigma_imj_grads;
372
384 std::vector<std::vector<double>> edge_sigma_values;
385
397 std::vector<std::vector<double>> edge_sigma_grads;
398
414 std::vector<std::vector<double>> edge_lambda_values;
415
425 std::vector<std::vector<double>> edge_lambda_grads_2d;
426
436 std::vector<std::vector<std::vector<double>>> edge_lambda_grads_3d;
437
448 std::vector<std::vector<std::vector<double>>> edge_lambda_gradgrads_3d;
449
466 std::vector<std::vector<double>> face_lambda_values;
467
476 std::vector<std::vector<double>> face_lambda_grads;
477 };
478
479private:
488 static std::vector<unsigned int>
489 get_dpo_vector(const unsigned int degree);
490
494 std::vector<Polynomials::Polynomial<double>> IntegratedLegendrePolynomials;
495
500 void
501 create_polynomials(const unsigned int degree);
502
506 unsigned int
507 compute_num_dofs(const unsigned int degree) const;
508
513 void
515 const Quadrature<dim> &quadrature,
516 const InternalData &fedata) const;
517
522 void
524 const Quadrature<dim> &quadrature,
525 const InternalData &fedata) const;
526
532};
533
534
535
539
540#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
void evaluate(const std::vector< Point< dim > > &p_list, const UpdateFlags update_flags, std::unique_ptr< typename ::FiniteElement< dim, spacedim >::InternalDataBase > &data_ptr) const
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
FE_NedelecSZ(const unsigned int order)
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 const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) 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
Threads::Mutex prolongation_matrix_mutex
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:452
friend class InternalDataBase
Definition fe.h:3075
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
UpdateFlags
MappingKind
Definition mapping.h:79