Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
fe_abf.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2018 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_abf_h
17 #define dealii_fe_abf_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/geometry_info.h>
22 #include <deal.II/base/polynomial.h>
23 #include <deal.II/base/polynomials_abf.h>
24 #include <deal.II/base/table.h>
25 #include <deal.II/base/tensor_product_polynomials.h>
26 
27 #include <deal.II/fe/fe.h>
28 #include <deal.II/fe/fe_poly_tensor.h>
29 
30 #include <vector>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 
37 
102 template <int dim>
103 class FE_ABF : public FE_PolyTensor<PolynomialsABF<dim>, dim>
104 {
105 public:
109  FE_ABF(const unsigned int p);
110 
116  virtual std::string
117  get_name() const override;
118 
126  virtual bool
127  has_support_on_face(const unsigned int shape_index,
128  const unsigned int face_index) const override;
129 
130  // documentation inherited from the base class
131  virtual void
133  const std::vector<Vector<double>> &support_point_values,
134  std::vector<double> & nodal_values) const override;
135 
136  virtual std::size_t
137  memory_consumption() const override;
138 
139  virtual std::unique_ptr<FiniteElement<dim, dim>>
140  clone() const override;
141 
142 private:
148  const unsigned int rt_order;
149 
156  static std::vector<unsigned int>
157  get_dpo_vector(const unsigned int degree);
158 
168  void
169  initialize_support_points(const unsigned int rt_degree);
170 
177  void
179 
186  class InternalData : public FiniteElement<dim>::InternalDataBase
187  {
188  public:
200  std::vector<std::vector<Tensor<1, dim>>> shape_values;
201 
211  std::vector<std::vector<Tensor<2, dim>>> shape_gradients;
212  };
213 
228 
229 
230 
245 
246 
250  template <int dim1>
251  friend class FE_ABF;
252 };
253 
254 
255 
259 DEAL_II_NAMESPACE_CLOSE
260 
261 #endif
void initialize_restriction()
Definition: fe_abf.cc:321
Table< 3, double > interior_weights
Definition: fe_abf.h:227
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) const override
Definition: fe_abf.cc:535
Table< 2, double > boundary_weights
Definition: fe_abf.h:220
virtual std::string get_name() const override
Definition: fe_abf.cc:116
const unsigned int degree
Definition: fe_base.h:298
friend class FE_ABF
Definition: fe_abf.h:251
std::vector< std::vector< Tensor< 1, dim > > > shape_values
Definition: fe_abf.h:200
Definition: fe_abf.h:103
Table< 2, double > boundary_weights_abf
Definition: fe_abf.h:237
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
Definition: fe_abf.cc:136
Table< 3, double > interior_weights_abf
Definition: fe_abf.h:244
std::vector< std::vector< Tensor< 2, dim > > > shape_gradients
Definition: fe_abf.h:211
const unsigned int rt_order
Definition: fe_abf.h:148
virtual std::size_t memory_consumption() const override
Definition: fe_abf.cc:612
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_abf.cc:458
void initialize_support_points(const unsigned int rt_degree)
Definition: fe_abf.cc:151
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition: fe_abf.cc:491