deal.II version GIT relicensing-2169-gec1b43f35b 2024-11-22 07:10:00+00:00
\(\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
75template <int dim, int spacedim = dim>
76class FE_NedelecSZ : public FiniteElement<dim, dim>
77{
78public:
79 static_assert(dim == spacedim,
80 "FE_NedelecSZ is only implemented for dim==spacedim!");
81
96 FE_NedelecSZ(const unsigned int order);
97
98 virtual UpdateFlags
99 requires_update_flags(const UpdateFlags update_flags) const override;
100
101 virtual std::string
102 get_name() const override;
103
104 virtual std::unique_ptr<FiniteElement<dim, dim>>
105 clone() const override;
106
115 virtual double
116 shape_value(const unsigned int i, const Point<dim> &p) const override;
117
124 virtual double
125 shape_value_component(const unsigned int i,
126 const Point<dim> &p,
127 const unsigned int component) const override;
128
137 virtual Tensor<1, dim>
138 shape_grad(const unsigned int i, const Point<dim> &p) const override;
139
146 virtual Tensor<1, dim>
147 shape_grad_component(const unsigned int i,
148 const Point<dim> &p,
149 const unsigned int component) const override;
150
159 virtual Tensor<2, dim>
160 shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
161
168 virtual Tensor<2, dim>
169 shape_grad_grad_component(const unsigned int i,
170 const Point<dim> &p,
171 const unsigned int component) const override;
172
193 virtual const FullMatrix<double> &
195 const unsigned int child,
196 const RefinementCase<dim> &refinement_case =
198
199protected:
205
210 void
211 evaluate(const std::vector<Point<dim>> &p_list,
212 const UpdateFlags update_flags,
213 std::unique_ptr<
214 typename ::FiniteElement<dim, spacedim>::InternalDataBase>
215 &data_ptr) const;
216
217 virtual std::unique_ptr<
218 typename ::FiniteElement<dim, spacedim>::InternalDataBase>
219 get_data(
220 const UpdateFlags update_flags,
221 const Mapping<dim, spacedim> &mapping,
222 const Quadrature<dim> &quadrature,
224 spacedim>
225 &output_data) const override;
226
232 virtual void
234 const typename Triangulation<dim, dim>::cell_iterator &cell,
235 const CellSimilarity::Similarity cell_similarity,
236 const Quadrature<dim> &quadrature,
237 const Mapping<dim, dim> &mapping,
238 const typename Mapping<dim, dim>::InternalDataBase &mapping_internal,
240 &mapping_data,
241 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
243 &data) const override;
244
245 using FiniteElement<dim, spacedim>::fill_fe_face_values;
246
252 virtual void
254 const typename Triangulation<dim, dim>::cell_iterator &cell,
255 const unsigned int face_no,
256 const hp::QCollection<dim - 1> &quadrature,
257 const Mapping<dim, dim> &mapping,
258 const typename Mapping<dim, dim>::InternalDataBase &mapping_internal,
260 &mapping_data,
261 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
263 &data) const override;
264
268 virtual void
270 const typename Triangulation<dim, dim>::cell_iterator &cell,
271 const unsigned int face_no,
272 const unsigned int sub_no,
273 const Quadrature<dim - 1> &quadrature,
274 const Mapping<dim, dim> &mapping,
275 const typename Mapping<dim, dim>::InternalDataBase &mapping_internal,
277 &mapping_data,
278 const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
280 &data) const override;
281
302 {
303 public:
310 mutable std::vector<std::vector<Tensor<1, dim>>> shape_values;
311
319 mutable std::vector<std::vector<DerivativeForm<1, dim, dim>>> shape_grads;
320
328 mutable std::vector<std::vector<DerivativeForm<2, dim, dim>>>
330
346 std::vector<std::vector<std::vector<double>>> sigma_imj_values;
347
363 std::vector<std::vector<std::vector<double>>> sigma_imj_grads;
364
376 std::vector<std::vector<double>> edge_sigma_values;
377
389 std::vector<std::vector<double>> edge_sigma_grads;
390
406 std::vector<std::vector<double>> edge_lambda_values;
407
417 std::vector<std::vector<double>> edge_lambda_grads_2d;
418
428 std::vector<std::vector<std::vector<double>>> edge_lambda_grads_3d;
429
440 std::vector<std::vector<std::vector<double>>> edge_lambda_gradgrads_3d;
441
458 std::vector<std::vector<double>> face_lambda_values;
459
468 std::vector<std::vector<double>> face_lambda_grads;
469 };
470
471private:
480 static std::vector<unsigned int>
481 get_dpo_vector(const unsigned int degree);
482
486 std::vector<Polynomials::Polynomial<double>> IntegratedLegendrePolynomials;
487
492 void
493 create_polynomials(const unsigned int degree);
494
498 unsigned int
499 compute_num_dofs(const unsigned int degree) const;
500
505 void
507 const Quadrature<dim> &quadrature,
508 const InternalData &fedata) const;
509
514 void
516 const Quadrature<dim> &quadrature,
517 const InternalData &fedata) const;
518
524};
525
526
527
531
532#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
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:3093
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
UpdateFlags
MappingKind
Definition mapping.h:79
std::vector< index_type > data
Definition mpi.cc:735