Reference documentation for deal.II version 8.5.1
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_nodal_values (const std::vector<Vector<double> > &support_point_values,
130  std::vector<double> &nodal_values) const;
131 
132  virtual void interpolate(std::vector<double> &local_dofs,
133  const std::vector<double> &values) const DEAL_II_DEPRECATED;
134 
135  virtual void interpolate(std::vector<double> &local_dofs,
136  const std::vector<Vector<double> > &values,
137  const unsigned int offset = 0) const DEAL_II_DEPRECATED;
138 
139  virtual void interpolate(std::vector<double> &local_dofs,
140  const VectorSlice<const std::vector<std::vector<double> > > &values) const DEAL_II_DEPRECATED;
141 
142  virtual std::size_t memory_consumption () const;
143  virtual FiniteElement<dim> *clone() const;
144 
145 private:
151  const unsigned int rt_order;
152 
159  static std::vector<unsigned int>
160  get_dpo_vector (const unsigned int degree);
161 
171  void initialize_support_points (const unsigned int rt_degree);
172 
179  void initialize_restriction ();
180 
187  class InternalData : public FiniteElement<dim>::InternalDataBase
188  {
189  public:
201  std::vector<std::vector<Tensor<1,dim> > > shape_values;
202 
212  std::vector<std::vector<Tensor<2,dim> > > shape_gradients;
213  };
214 
229 
230 
231 
246 
247 
251  template <int dim1> friend class FE_ABF;
252 };
253 
254 
255 
259 DEAL_II_NAMESPACE_CLOSE
260 
261 #endif
std::vector< std::vector< Tensor< 1, dim > > > shape_values
Definition: fe_abf.h:201
void initialize_restriction()
Definition: fe_abf.cc:306
Table< 3, double > interior_weights
Definition: fe_abf.h:228
Table< 2, double > boundary_weights
Definition: fe_abf.h:221
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
Definition: fe_abf.cc:575
const unsigned int degree
Definition: fe_base.h:311
virtual std::string get_name() const
Definition: fe_abf.cc:110
virtual std::size_t memory_consumption() const
Definition: fe_abf.cc:711
friend class FE_ABF
Definition: fe_abf.h:251
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_abf.cc:470
Definition: fe_abf.h:101
Table< 2, double > boundary_weights_abf
Definition: fe_abf.h:238
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
Definition: fe_abf.cc:514
Table< 3, double > interior_weights_abf
Definition: fe_abf.h:245
virtual FiniteElement< dim > * clone() const
Definition: fe_abf.cc:130
const unsigned int rt_order
Definition: fe_abf.h:151
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_abf.cc:436
void initialize_support_points(const unsigned int rt_degree)
Definition: fe_abf.cc:145
std::vector< std::vector< Tensor< 2, dim > > > shape_gradients
Definition: fe_abf.h:212