Reference documentation for deal.II version 8.5.1
fe_poly_face.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2016 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__fe_poly_face_h
17 #define dealii__fe_poly_face_h
18 
19 
20 #include <deal.II/base/qprojector.h>
21 #include <deal.II/fe/fe.h>
22 
23 
24 DEAL_II_NAMESPACE_OPEN
25 
28 
54 template <class PolynomialType, int dim=PolynomialType::dimension+1, int spacedim=dim>
55 class FE_PolyFace : public FiniteElement<dim,spacedim>
56 {
57 public:
61  FE_PolyFace (const PolynomialType &poly_space,
62  const FiniteElementData<dim> &fe_data,
63  const std::vector<bool> &restriction_is_additive_flags);
64 
69  unsigned int get_degree () const;
70 
71  // for documentation, see the FiniteElement base class
72  virtual
74  requires_update_flags (const UpdateFlags update_flags) const;
75 
76 protected:
77  /*
78  * NOTE: The following functions have their definitions inlined into the class declaration
79  * because we otherwise run into a compiler error with MS Visual Studio.
80  */
81 
82 
83  virtual
85  get_data (const UpdateFlags /*update_flags*/,
86  const Mapping<dim,spacedim> &/*mapping*/,
87  const Quadrature<dim> &/*quadrature*/,
89  {
90  InternalData *data = new InternalData;
91  return data;
92  }
93 
95  get_face_data(const UpdateFlags update_flags,
96  const Mapping<dim,spacedim> &/*mapping*/,
97  const Quadrature<dim-1> &quadrature,
99  {
100  // generate a new data object and
101  // initialize some fields
102  InternalData *data = new InternalData;
103  data->update_each = requires_update_flags(update_flags);
104 
105  const unsigned int n_q_points = quadrature.size();
106 
107  // some scratch arrays
108  std::vector<double> values(0);
109  std::vector<Tensor<1,dim-1> > grads(0);
110  std::vector<Tensor<2,dim-1> > grad_grads(0);
111  std::vector<Tensor<3,dim-1> > empty_vector_of_3rd_order_tensors;
112  std::vector<Tensor<4,dim-1> > empty_vector_of_4th_order_tensors;
113 
114  // initialize fields only if really
115  // necessary. otherwise, don't
116  // allocate memory
117  if (data->update_each & update_values)
118  {
119  values.resize (poly_space.n());
120  data->shape_values.resize (poly_space.n(),
121  std::vector<double> (n_q_points));
122  for (unsigned int i=0; i<n_q_points; ++i)
123  {
124  poly_space.compute(quadrature.point(i),
125  values, grads, grad_grads,
126  empty_vector_of_3rd_order_tensors,
127  empty_vector_of_4th_order_tensors);
128 
129  for (unsigned int k=0; k<poly_space.n(); ++k)
130  data->shape_values[k][i] = values[k];
131  }
132  }
133  // No derivatives of this element
134  // are implemented.
135  if (data->update_each & update_gradients || data->update_each & update_hessians)
136  {
137  Assert(false, ExcNotImplemented());
138  }
139 
140  return data;
141  }
142 
144  get_subface_data(const UpdateFlags update_flags,
145  const Mapping<dim,spacedim> &mapping,
146  const Quadrature<dim-1> &quadrature,
148  {
149  return get_face_data(update_flags, mapping,
151  output_data);
152  }
153 
154  virtual
155  void
156  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
157  const CellSimilarity::Similarity cell_similarity,
158  const Quadrature<dim> &quadrature,
159  const Mapping<dim,spacedim> &mapping,
160  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
161  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
162  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
164 
165  virtual
166  void
167  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
168  const unsigned int face_no,
169  const Quadrature<dim-1> &quadrature,
170  const Mapping<dim,spacedim> &mapping,
171  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
172  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
173  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
175 
176  virtual
177  void
178  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
179  const unsigned int face_no,
180  const unsigned int sub_no,
181  const Quadrature<dim-1> &quadrature,
182  const Mapping<dim,spacedim> &mapping,
183  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
184  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
185  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
187 
194  class InternalData : public FiniteElement<dim,spacedim>::InternalDataBase
195  {
196  public:
210  std::vector<std::vector<double> > shape_values;
211  };
212 
217  PolynomialType poly_space;
218 };
219 
222 DEAL_II_NAMESPACE_CLOSE
223 
224 #endif
Shape function values.
std::vector< std::vector< double > > shape_values
Definition: fe_poly_face.h:210
const Point< dim > & point(const unsigned int i) const
virtual FiniteElement< dim, spacedim >::InternalDataBase * get_data(const UpdateFlags, const Mapping< dim, spacedim > &, const Quadrature< dim > &, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &) const
Definition: fe_poly_face.h:85
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
Second derivatives of shape functions.
PolynomialType poly_space
Definition: fe_poly_face.h:217
unsigned int size() const
Definition: mpi.h:41
unsigned int get_degree() const
FiniteElement< dim, spacedim >::InternalDataBase * get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe_poly_face.h:144
Shape function gradients.
static ::ExceptionBase & ExcNotImplemented()
FE_PolyFace(const PolynomialType &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2346
UpdateFlags update_each
Definition: fe.h:638
FiniteElement< dim, spacedim >::InternalDataBase * get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &) const
Definition: fe_poly_face.h:95