Reference documentation for deal.II version 9.1.0-pre
fe_nedelec_sz.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2015 - 2017 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/derivative_form.h>
20 #include <deal.II/base/polynomials_integrated_legendre_sz.h>
21 #include <deal.II/base/qprojector.h>
22 
23 #include <deal.II/fe/fe.h>
24 #include <deal.II/fe/fe_values.h>
25 #include <deal.II/fe/mapping.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
31 
70 template <int dim, int spacedim = dim>
71 class FE_NedelecSZ : public FiniteElement<dim, dim>
72 {
73 public:
74  static_assert(dim == spacedim,
75  "FE_NedelecSZ is only implemented for dim==spacedim!");
76 
80  FE_NedelecSZ(const unsigned int degree);
81 
82  virtual UpdateFlags
83  requires_update_flags(const UpdateFlags update_flags) const override;
84 
85  virtual std::string
86  get_name() const override;
87 
88  virtual std::unique_ptr<FiniteElement<dim, dim>>
89  clone() const override;
90 
95  virtual double
96  shape_value(const unsigned int i, const Point<dim> &p) const override;
97 
101  virtual double
102  shape_value_component(const unsigned int i,
103  const Point<dim> & p,
104  const unsigned int component) const override;
105 
110  virtual Tensor<1, dim>
111  shape_grad(const unsigned int i, const Point<dim> &p) const override;
112 
116  virtual Tensor<1, dim>
117  shape_grad_component(const unsigned int i,
118  const Point<dim> & p,
119  const unsigned int component) const override;
120 
125  virtual Tensor<2, dim>
126  shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
127 
131  virtual Tensor<2, dim>
132  shape_grad_grad_component(const unsigned int i,
133  const Point<dim> & p,
134  const unsigned int component) const override;
135 
142  update_once(const UpdateFlags flags) const;
143 
150  update_each(const UpdateFlags flags) const;
151 
152 protected:
158 
159  virtual std::unique_ptr<
160  typename ::FiniteElement<dim, spacedim>::InternalDataBase>
161  get_data(
162  const UpdateFlags update_flags,
163  const Mapping<dim, spacedim> &mapping,
164  const Quadrature<dim> & quadrature,
166  spacedim>
167  &output_data) const override;
168 
174  virtual void
176  const typename Triangulation<dim, dim>::cell_iterator &cell,
177  const CellSimilarity::Similarity cell_similarity,
178  const Quadrature<dim> & quadrature,
179  const Mapping<dim, dim> & mapping,
180  const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
181  const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
182  & mapping_data,
183  const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
185  &data) const override;
186 
192  virtual void
194  const typename Triangulation<dim, dim>::cell_iterator &cell,
195  const unsigned int face_no,
196  const Quadrature<dim - 1> & quadrature,
197  const Mapping<dim, dim> & mapping,
198  const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
199  const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
200  & mapping_data,
201  const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
203  &data) const override;
204 
208  virtual void
210  const typename Triangulation<dim, dim>::cell_iterator &cell,
211  const unsigned int face_no,
212  const unsigned int sub_no,
213  const Quadrature<dim - 1> & quadrature,
214  const Mapping<dim, dim> & mapping,
215  const typename Mapping<dim, dim>::InternalDataBase & mapping_internal,
216  const ::internal::FEValuesImplementation::MappingRelatedData<dim, dim>
217  & mapping_data,
218  const typename FiniteElement<dim, dim>::InternalDataBase &fedata,
220  &data) const override;
221 
241  class InternalData : public FiniteElement<dim, dim>::InternalDataBase
242  {
243  public:
250  mutable std::vector<std::vector<Tensor<1, dim>>> shape_values;
251 
259  mutable std::vector<std::vector<DerivativeForm<1, dim, dim>>> shape_grads;
260 
276  std::vector<std::vector<std::vector<double>>> sigma_imj_values;
277 
293  std::vector<std::vector<std::vector<double>>> sigma_imj_grads;
294 
306  std::vector<std::vector<double>> edge_sigma_values;
307 
320  std::vector<std::vector<double>> edge_sigma_grads;
321 
337  std::vector<std::vector<double>> edge_lambda_values;
338 
348  std::vector<std::vector<double>> edge_lambda_grads_2d;
349 
359  std::vector<std::vector<std::vector<double>>> edge_lambda_grads_3d;
360 
371  std::vector<std::vector<std::vector<double>>> edge_lambda_gradgrads_3d;
372 
389  std::vector<std::vector<double>> face_lambda_values;
390 
399  std::vector<std::vector<double>> face_lambda_grads;
400  };
401 
402 private:
411  static std::vector<unsigned int>
412  get_dpo_vector(const unsigned int degree);
413 
417  std::vector<Polynomials::Polynomial<double>> IntegratedLegendrePolynomials;
418 
423  void
424  create_polynomials(const unsigned int degree);
425 
429  unsigned int
430  compute_num_dofs(const unsigned int degree) const;
431 
436  void
438  const Quadrature<dim> &quadrature,
439  const InternalData & fedata) const;
440 
445  void
447  const Quadrature<dim> &quadrature,
448  const InternalData & fedata) const;
449 };
450 
451 
452 
455 DEAL_II_NAMESPACE_CLOSE
456 
457 #endif
unsigned int compute_num_dofs(const unsigned int degree) const
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
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
MappingType
Definition: mapping.h:61
std::vector< std::vector< double > > edge_sigma_values
const unsigned int degree
Definition: fe_base.h:313
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
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
void fill_edge_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
MappingType mapping_type
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 std::string get_name() const override
std::vector< std::vector< DerivativeForm< 1, dim, dim > > > shape_grads
std::vector< std::vector< std::vector< double > > > edge_lambda_gradgrads_3d
UpdateFlags
UpdateFlags update_once(const UpdateFlags flags) const
void fill_face_values(const typename Triangulation< dim, dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, const InternalData &fedata) const
Abstract base class for mapping classes.
Definition: dof_tools.h:57
UpdateFlags update_each(const UpdateFlags flags) const
std::vector< Polynomials::Polynomial< double > > IntegratedLegendrePolynomials
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 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 Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
std::vector< std::vector< std::vector< double > > > edge_lambda_grads_3d
FE_NedelecSZ(const unsigned int degree)
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
std::vector< std::vector< std::vector< double > > > sigma_imj_grads