Reference documentation for deal.II version 9.2.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\}}\)
fe_poly_face.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2019 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_poly_face_h
17 #define dealii_fe_poly_face_h
18 
19 
20 #include <deal.II/base/config.h>
21 
24 
25 #include <deal.II/fe/fe.h>
26 
27 
29 
32 
58 template <class PolynomialType,
59  int dim = PolynomialType::dimension + 1,
60  int spacedim = dim>
61 class FE_PolyFace : public FiniteElement<dim, spacedim>
62 {
63 public:
68  const FiniteElementData<dim> &fe_data,
69  const std::vector<bool> & restriction_is_additive_flags);
70 
75  unsigned int
76  get_degree() const;
77 
78  // for documentation, see the FiniteElement base class
79  virtual UpdateFlags
80  requires_update_flags(const UpdateFlags update_flags) const override;
81 
82 protected:
83  /*
84  * NOTE: The following functions have their definitions inlined into the class
85  * declaration because we otherwise run into a compiler error with MS Visual
86  * Studio.
87  */
88 
89 
90  virtual std::unique_ptr<
93  const UpdateFlags /*update_flags*/,
94  const Mapping<dim, spacedim> & /*mapping*/,
95  const Quadrature<dim> & /*quadrature*/,
97  spacedim>
98  & /*output_data*/) const override
99  {
100  return std_cxx14::make_unique<InternalData>();
101  }
102 
103  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
105  const UpdateFlags update_flags,
106  const Mapping<dim, spacedim> & /*mapping*/,
107  const Quadrature<dim - 1> &quadrature,
109  spacedim>
110  & /*output_data*/) const override
111  {
112  // generate a new data object and
113  // initialize some fields
114  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
115  data_ptr = std_cxx14::make_unique<InternalData>();
116  auto &data = dynamic_cast<InternalData &>(*data_ptr);
117  data.update_each = requires_update_flags(update_flags);
118 
119  const unsigned int n_q_points = quadrature.size();
120 
121  // some scratch arrays
122  std::vector<double> values(0);
123  std::vector<Tensor<1, dim - 1>> grads(0);
124  std::vector<Tensor<2, dim - 1>> grad_grads(0);
125  std::vector<Tensor<3, dim - 1>> empty_vector_of_3rd_order_tensors;
126  std::vector<Tensor<4, dim - 1>> empty_vector_of_4th_order_tensors;
127 
128  // initialize fields only if really
129  // necessary. otherwise, don't
130  // allocate memory
131  if (data.update_each & update_values)
132  {
133  values.resize(poly_space.n());
134  data.shape_values.resize(poly_space.n(),
135  std::vector<double>(n_q_points));
136  for (unsigned int i = 0; i < n_q_points; ++i)
137  {
138  poly_space.evaluate(quadrature.point(i),
139  values,
140  grads,
141  grad_grads,
142  empty_vector_of_3rd_order_tensors,
143  empty_vector_of_4th_order_tensors);
144 
145  for (unsigned int k = 0; k < poly_space.n(); ++k)
146  data.shape_values[k][i] = values[k];
147  }
148  }
149  // No derivatives of this element
150  // are implemented.
151  if (data.update_each & update_gradients ||
152  data.update_each & update_hessians)
153  {
154  Assert(false, ExcNotImplemented());
155  }
156 
157  return data_ptr;
158  }
159 
160  std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
162  const UpdateFlags update_flags,
163  const Mapping<dim, spacedim> &mapping,
164  const Quadrature<dim - 1> & quadrature,
166  spacedim>
167  &output_data) const override
168  {
169  return get_face_data(update_flags,
170  mapping,
172  quadrature),
173  output_data);
174  }
175 
176  virtual void
178  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
179  const CellSimilarity::Similarity cell_similarity,
180  const Quadrature<dim> & quadrature,
181  const Mapping<dim, spacedim> & mapping,
182  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
183  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
184  spacedim>
185  & mapping_data,
186  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
188  spacedim>
189  &output_data) const override;
190 
191  virtual void
193  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
194  const unsigned int face_no,
195  const Quadrature<dim - 1> & quadrature,
196  const Mapping<dim, spacedim> & mapping,
197  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
198  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
199  spacedim>
200  & mapping_data,
201  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
203  spacedim>
204  &output_data) const override;
205 
206  virtual void
208  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
209  const unsigned int face_no,
210  const unsigned int sub_no,
211  const Quadrature<dim - 1> & quadrature,
212  const Mapping<dim, spacedim> & mapping,
213  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
214  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
215  spacedim>
216  & mapping_data,
217  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
219  spacedim>
220  &output_data) const override;
221 
228  class InternalData : public FiniteElement<dim, spacedim>::InternalDataBase
229  {
230  public:
244  std::vector<std::vector<double>> shape_values;
245  };
246 
252 };
253 
257 
258 #endif
FE_PolyFace::get_subface_data
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
Definition: fe_poly_face.h:161
FE_PolyFace::poly_space
PolynomialType poly_space
Definition: fe_poly_face.h:251
FE_PolyFace
Definition: fe_poly_face.h:61
internal::FEValuesImplementation::FiniteElementRelatedData
Definition: fe_update_flags.h:524
FE_PolyFace::FE_PolyFace
FE_PolyFace(const PolynomialType &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
FiniteElementData
Definition: fe_base.h:147
FiniteElement< PolynomialType::dimension+1, PolynomialType::dimension+1 >::InternalDataBase
friend class InternalDataBase
Definition: fe.h:3048
CellSimilarity::Similarity
Similarity
Definition: fe_update_flags.h:378
FE_PolyFace::get_data
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Mapping< dim, spacedim > &, const Quadrature< dim > &, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &) const override
Definition: fe_poly_face.h:92
FE_PolyFace::InternalData
Definition: fe_poly_face.h:228
update_values
@ update_values
Shape function values.
Definition: fe_update_flags.h:78
QProjector
Definition: qprojector.h:78
Quadrature::point
const Point< dim > & point(const unsigned int i) const
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
FE_PolyFace::get_face_data
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &) const override
Definition: fe_poly_face.h:104
Tensor
Definition: tensor.h:450
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
update_gradients
@ update_gradients
Shape function gradients.
Definition: fe_update_flags.h:84
fe.h
FE_PolyFace::InternalData::shape_values
std::vector< std::vector< double > > shape_values
Definition: fe_poly_face.h:244
FiniteElement
Definition: fe.h:648
UpdateFlags
UpdateFlags
Definition: fe_update_flags.h:66
FE_PolyFace::requires_update_flags
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
FE_PolyFace::get_degree
unsigned int get_degree() const
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
update_hessians
@ update_hessians
Second derivatives of shape functions.
Definition: fe_update_flags.h:90
Mapping::InternalDataBase
Definition: mapping.h:597
FiniteElement::InternalDataBase
Definition: fe.h:682
FiniteElement< PolynomialType::dimension+1, PolynomialType::dimension+1 >::restriction_is_additive_flags
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2581
config.h
Quadrature::size
unsigned int size() const
memory.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
FE_PolyFace::fill_fe_face_values
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
Quadrature
Definition: quadrature.h:85
TriaIterator
Definition: tria_iterator.h:578
qprojector.h
FE_PolyFace::fill_fe_values
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
PolynomialType
FE_PolyFace::fill_fe_subface_values
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override