Reference documentation for deal.II version 9.0.0
fe_abf.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_fe_abf_h
17 #define dealii_fe_abf_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/table.h>
21 #include <deal.II/base/polynomials_abf.h>
22 #include <deal.II/base/polynomial.h>
23 #include <deal.II/base/tensor_product_polynomials.h>
24 #include <deal.II/base/geometry_info.h>
25 #include <deal.II/fe/fe.h>
26 #include <deal.II/fe/fe_poly_tensor.h>
27 
28 #include <vector>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
35 
100 template <int dim>
101 class FE_ABF : public FE_PolyTensor<PolynomialsABF<dim>, dim>
102 {
103 public:
107  FE_ABF (const unsigned int p);
108 
114  virtual std::string get_name () const;
115 
123  virtual bool has_support_on_face (const unsigned int shape_index,
124  const unsigned int face_index) const;
125 
126  // documentation inherited from the base class
127  virtual
128  void
129  convert_generalized_support_point_values_to_dof_values (const std::vector<Vector<double> > &support_point_values,
130  std::vector<double> &nodal_values) const;
131 
132  virtual std::size_t memory_consumption () const;
133 
134  virtual
135  std::unique_ptr<FiniteElement<dim,dim> >
136  clone() const;
137 
138 private:
144  const unsigned int rt_order;
145 
152  static std::vector<unsigned int>
153  get_dpo_vector (const unsigned int degree);
154 
164  void initialize_support_points (const unsigned int rt_degree);
165 
172  void initialize_restriction ();
173 
180  class InternalData : public FiniteElement<dim>::InternalDataBase
181  {
182  public:
194  std::vector<std::vector<Tensor<1,dim> > > shape_values;
195 
205  std::vector<std::vector<Tensor<2,dim> > > shape_gradients;
206  };
207 
222 
223 
224 
239 
240 
244  template <int dim1> friend class FE_ABF;
245 };
246 
247 
248 
252 DEAL_II_NAMESPACE_CLOSE
253 
254 #endif
std::vector< std::vector< Tensor< 1, dim > > > shape_values
Definition: fe_abf.h:194
void initialize_restriction()
Definition: fe_abf.cc:308
Table< 3, double > interior_weights
Definition: fe_abf.h:221
Table< 2, double > boundary_weights
Definition: fe_abf.h:214
const unsigned int degree
Definition: fe_base.h:311
virtual std::string get_name() const
Definition: fe_abf.cc:112
virtual std::size_t memory_consumption() const
Definition: fe_abf.cc:577
friend class FE_ABF
Definition: fe_abf.h:244
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
Definition: fe_abf.cc:516
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_abf.cc:472
Definition: fe_abf.h:101
Table< 2, double > boundary_weights_abf
Definition: fe_abf.h:231
Table< 3, double > interior_weights_abf
Definition: fe_abf.h:238
const unsigned int rt_order
Definition: fe_abf.h:144
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_abf.cc:438
void initialize_support_points(const unsigned int rt_degree)
Definition: fe_abf.cc:147
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const
Definition: fe_abf.cc:132
std::vector< std::vector< Tensor< 2, dim > > > shape_gradients
Definition: fe_abf.h:205