Reference documentation for deal.II version Git 193422c69f 2020-07-08 17:07:46 +0200
\(\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\}}\)
fe_nedelec_sz.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2015 - 2020 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>
26 #include <deal.II/fe/fe_values.h>
27 #include <deal.II/fe/mapping.h>
28 
30 
33 
72 template <int dim, int spacedim = dim>
73 class FE_NedelecSZ : public FiniteElement<dim, dim>
74 {
75 public:
76  static_assert(dim == spacedim,
77  "FE_NedelecSZ is only implemented for dim==spacedim!");
78 
91  FE_NedelecSZ(const unsigned int order);
92 
93  virtual UpdateFlags
94  requires_update_flags(const UpdateFlags update_flags) const override;
95 
96  virtual std::string
97  get_name() const override;
98 
99  virtual std::unique_ptr<FiniteElement<dim, dim>>
100  clone() const override;
101 
106  virtual double
107  shape_value(const unsigned int i, const Point<dim> &p) const override;
108 
112  virtual double
113  shape_value_component(const unsigned int i,
114  const Point<dim> & p,
115  const unsigned int component) const override;
116 
121  virtual Tensor<1, dim>
122  shape_grad(const unsigned int i, const Point<dim> &p) const override;
123 
127  virtual Tensor<1, dim>
128  shape_grad_component(const unsigned int i,
129  const Point<dim> & p,
130  const unsigned int component) const override;
131 
136  virtual Tensor<2, dim>
137  shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
138 
142  virtual Tensor<2, dim>
143  shape_grad_grad_component(const unsigned int i,
144  const Point<dim> & p,
145  const unsigned int component) const override;
146 
147 protected:
153 
154  virtual std::unique_ptr<
155  typename ::FiniteElement<dim, spacedim>::InternalDataBase>
156  get_data(
157  const UpdateFlags update_flags,
158  const Mapping<dim, spacedim> &mapping,
159  const Quadrature<dim> & quadrature,
161  spacedim>
162  &output_data) const override;
163 
169  virtual void
171  const typename Triangulation<dim, dim>::cell_iterator &cell,
172  const CellSimilarity::Similarity cell_similarity,
173  const Quadrature<dim> & quadrature,
174  const Mapping<dim, dim> & mapping,
175  const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
176  const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
177  & mapping_data,
178  const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
180  &data) const override;
181 
187  virtual void
189  const typename Triangulation<dim, dim>::cell_iterator &cell,
190  const unsigned int face_no,
191  const Quadrature<dim - 1> & quadrature,
192  const Mapping<dim, dim> & mapping,
193  const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
194  const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
195  & mapping_data,
196  const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
198  &data) const override;
199 
203  virtual void
205  const typename Triangulation<dim, dim>::cell_iterator &cell,
206  const unsigned int face_no,
207  const unsigned int sub_no,
208  const Quadrature<dim - 1> & quadrature,
209  const Mapping<dim, dim> & mapping,
210  const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
211  const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
212  & mapping_data,
213  const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
215  &data) const override;
216 
237  {
238  public:
245  mutable std::vector<std::vector<Tensor<1, dim>>> shape_values;
246 
254  mutable std::vector<std::vector<DerivativeForm<1, dim, dim>>> shape_grads;
255 
271  std::vector<std::vector<std::vector<double>>> sigma_imj_values;
272 
288  std::vector<std::vector<std::vector<double>>> sigma_imj_grads;
289 
301  std::vector<std::vector<double>> edge_sigma_values;
302 
314  std::vector<std::vector<double>> edge_sigma_grads;
315 
331  std::vector<std::vector<double>> edge_lambda_values;
332 
342  std::vector<std::vector<double>> edge_lambda_grads_2d;
343 
353  std::vector<std::vector<std::vector<double>>> edge_lambda_grads_3d;
354 
365  std::vector<std::vector<std::vector<double>>> edge_lambda_gradgrads_3d;
366 
383  std::vector<std::vector<double>> face_lambda_values;
384 
393  std::vector<std::vector<double>> face_lambda_grads;
394  };
395 
396 private:
405  static std::vector<unsigned int>
406  get_dpo_vector(const unsigned int degree);
407 
411  std::vector<Polynomials::Polynomial<double>> IntegratedLegendrePolynomials;
412 
417  void
418  create_polynomials(const unsigned int degree);
419 
423  unsigned int
424  compute_num_dofs(const unsigned int degree) const;
425 
430  void
432  const Quadrature<dim> &quadrature,
433  const InternalData & fedata) const;
434 
439  void
441  const Quadrature<dim> &quadrature,
442  const InternalData & fedata) const;
443 };
444 
445 
446 
450 
451 #endif
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
std::vector< std::vector< double > > edge_lambda_grads_2d
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
std::vector< std::vector< double > > face_lambda_values
void fill_edge_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
std::vector< std::vector< double > > edge_sigma_values
const unsigned int degree
Definition: fe_base.h:296
MappingKind
Definition: mapping.h:62
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
MappingKind mapping_kind
void fill_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::string get_name() const override
virtual void fill_fe_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const unsigned int face_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
std::vector< std::vector< DerivativeForm< 1, dim, dim > > > shape_grads
std::vector< std::vector< std::vector< double > > > edge_lambda_gradgrads_3d
UpdateFlags
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
Abstract base class for mapping classes.
Definition: mapping.h:301
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:363
std::vector< Polynomials::Polynomial< double > > IntegratedLegendrePolynomials
friend class InternalDataBase
Definition: fe.h:3028
std::vector< std::vector< Tensor< 1, dim > > > shape_values
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
void create_polynomials(const unsigned int degree)
std::vector< std::vector< double > > edge_lambda_values
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:362
std::vector< std::vector< std::vector< double > > > edge_lambda_grads_3d
std::vector< std::vector< double > > edge_sigma_grads
std::vector< std::vector< std::vector< double > > > sigma_imj_values
std::vector< std::vector< double > > face_lambda_grads
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
std::vector< std::vector< std::vector< double > > > sigma_imj_grads
FE_NedelecSZ(const unsigned int order)
unsigned int compute_num_dofs(const unsigned int degree) const