Reference documentation for deal.II version 8.5.1
fe_q_bubbles.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 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 
17 #ifndef dealii__fe_q_bubbles_h
18 #define dealii__fe_q_bubbles_h
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/tensor_product_polynomials_bubbles.h>
22 #include <deal.II/fe/fe_q_base.h>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 
29 
81 template <int dim, int spacedim=dim>
82 class FE_Q_Bubbles : public FE_Q_Base<TensorProductPolynomialsBubbles<dim>,dim,spacedim>
83 {
84 public:
90  FE_Q_Bubbles (const unsigned int p);
91 
98  FE_Q_Bubbles (const Quadrature<1> &points);
99 
105  virtual std::string get_name () const;
106 
107  // documentation inherited from the base class
108  virtual
109  void
110  convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
111  std::vector<double> &nodal_values) const;
112 
113  virtual void interpolate(std::vector<double> &local_dofs,
114  const std::vector<double> &values) const DEAL_II_DEPRECATED;
115 
116  virtual void interpolate(std::vector<double> &local_dofs,
117  const std::vector<Vector<double> > &values,
118  const unsigned int offset = 0) const DEAL_II_DEPRECATED;
119 
120  virtual void interpolate(std::vector<double> &local_dofs,
121  const VectorSlice<const std::vector<std::vector<double> > > &values) const DEAL_II_DEPRECATED;
122 
132  virtual void
134  FullMatrix<double> &matrix) const;
135 
136  virtual const FullMatrix<double> &
137  get_prolongation_matrix (const unsigned int child,
138  const RefinementCase<dim> &refinement_case) const;
139 
140  virtual const FullMatrix<double> &
141  get_restriction_matrix (const unsigned int child,
142  const RefinementCase<dim> &refinement_case) const;
143 
152  virtual bool has_support_on_face (const unsigned int shape_index,
153  const unsigned int face_index) const;
154 
155 protected:
161  virtual FiniteElement<dim,spacedim> *clone() const;
162 
163 private:
164 
169  static std::vector<bool> get_riaf_vector(const unsigned int degree);
170 
177  static std::vector<unsigned int> get_dpo_vector(const unsigned int degree);
178 
182  const unsigned int n_bubbles;
183 };
184 
185 
186 
190 DEAL_II_NAMESPACE_CLOSE
191 
192 #endif
const unsigned int n_bubbles
Definition: fe_q_bubbles.h:182
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
const unsigned int degree
Definition: fe_base.h:311
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case) const
virtual FiniteElement< dim, spacedim > * clone() const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case) const
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual std::string get_name() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
static std::vector< bool > get_riaf_vector(const unsigned int degree)
FE_Q_Bubbles(const unsigned int p)
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const