Reference documentation for deal.II version 8.4.2
fe_q_bubbles.h
1 // ---------------------------------------------------------------------
2 // @f$Id@f$
3 //
4 // Copyright (C) 2012 - 2016 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 
18 #ifndef dealii__fe_q_bubbles_h
19 #define dealii__fe_q_bubbles_h
20 
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/tensor_product_polynomials_bubbles.h>
23 #include <deal.II/fe/fe_q_base.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 
30 
82 template <int dim, int spacedim=dim>
83 class FE_Q_Bubbles : public FE_Q_Base<TensorProductPolynomialsBubbles<dim>,dim,spacedim>
84 {
85 public:
91  FE_Q_Bubbles (const unsigned int p);
92 
99  FE_Q_Bubbles (const Quadrature<1> &points);
100 
106  virtual std::string get_name () const;
107 
112  virtual void interpolate(std::vector<double> &local_dofs,
113  const std::vector<double> &values) const;
114 
124  virtual void interpolate(std::vector<double> &local_dofs,
125  const std::vector<Vector<double> > &values,
126  unsigned int offset = 0) const;
127 
132  virtual void interpolate(
133  std::vector<double> &local_dofs,
134  const VectorSlice<const std::vector<std::vector<double> > > &values) const;
135 
145  virtual void
147  FullMatrix<double> &matrix) const;
148 
149  virtual const FullMatrix<double> &
150  get_prolongation_matrix (const unsigned int child,
151  const RefinementCase<dim> &refinement_case) const;
152 
153  virtual const FullMatrix<double> &
154  get_restriction_matrix (const unsigned int child,
155  const RefinementCase<dim> &refinement_case) const;
156 
165  virtual bool has_support_on_face (const unsigned int shape_index,
166  const unsigned int face_index) const;
167 
168 protected:
174  virtual FiniteElement<dim,spacedim> *clone() const;
175 
176 private:
177 
182  static std::vector<bool> get_riaf_vector(const unsigned int degree);
183 
190  static std::vector<unsigned int> get_dpo_vector(const unsigned int degree);
191 
195  const unsigned int n_bubbles;
196 };
197 
198 
199 
203 DEAL_II_NAMESPACE_CLOSE
204 
205 #endif
const unsigned int n_bubbles
Definition: fe_q_bubbles.h:195
const unsigned int degree
Definition: fe_base.h:299
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 void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) 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)