Reference documentation for deal.II version 9.0.0
fe_dgq.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_fe_dgq_h
17 #define dealii_fe_dgq_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/tensor_product_polynomials.h>
21 #include <deal.II/base/thread_management.h>
22 #include <deal.II/fe/fe_poly.h>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 template <int dim, int spacedim> class MappingQ;
27 template <int dim> class Quadrature;
28 
31 
104 template <int dim, int spacedim=dim>
105 class FE_DGQ : public FE_Poly<TensorProductPolynomials<dim>, dim, spacedim>
106 {
107 public:
114  FE_DGQ (const unsigned int p);
115 
121  virtual std::string get_name () const;
122 
132  virtual void
134  FullMatrix<double> &matrix) const;
135 
147  virtual void
149  FullMatrix<double> &matrix) const;
150 
162  virtual void
164  const unsigned int subface,
165  FullMatrix<double> &matrix) const;
166 
184  virtual const FullMatrix<double> &
185  get_restriction_matrix (const unsigned int child,
187 
209  virtual const FullMatrix<double> &
210  get_prolongation_matrix (const unsigned int child,
212 
236  virtual
237  std::vector<std::pair<unsigned int, unsigned int> >
239 
247  virtual
248  std::vector<std::pair<unsigned int, unsigned int> >
250 
258  virtual
259  std::vector<std::pair<unsigned int, unsigned int> >
261 
270  virtual bool hp_constraints_are_implemented () const;
271 
281  virtual
284 
293  virtual bool has_support_on_face (const unsigned int shape_index,
294  const unsigned int face_index) const;
295 
300  virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
301  get_constant_modes () const;
302 
310  virtual
311  void
312  convert_generalized_support_point_values_to_dof_values (const std::vector<Vector<double> > &support_point_values,
313  std::vector<double> &nodal_values) const;
314 
323  virtual std::size_t memory_consumption () const;
324 
325  virtual
326  std::unique_ptr<FiniteElement<dim,spacedim> >
327  clone() const;
328 
329 protected:
338  FE_DGQ (const std::vector<Polynomials::Polynomial<double> > &polynomials);
339 
340 private:
347  static std::vector<unsigned int> get_dpo_vector (const unsigned int degree);
348 
365  void rotate_indices (std::vector<unsigned int> &indices,
366  const char direction) const;
367 
368  /*
369  * Mutex for protecting initialization of restriction and embedding matrix.
370  */
371  mutable Threads::Mutex mutex;
372 
376  template <int dim1, int spacedim1> friend class FE_DGQ;
377 
381  template <int dim1, int spacedim1> friend class MappingQ;
382 };
383 
384 
385 
402 template <int dim,int spacedim=dim>
403 class FE_DGQArbitraryNodes : public FE_DGQ<dim,spacedim>
404 {
405 public:
412  FE_DGQArbitraryNodes (const Quadrature<1> &points);
413 
419  virtual std::string get_name () const;
420 
428  virtual
429  void
430  convert_generalized_support_point_values_to_dof_values (const std::vector<Vector<double> > &support_point_values,
431  std::vector<double> &nodal_values) const;
432  virtual
433  std::unique_ptr<FiniteElement<dim,spacedim> >
434  clone() const;
435 };
436 
437 
438 
449 template <int dim,int spacedim=dim>
450 class FE_DGQLegendre : public FE_DGQ<dim,spacedim>
451 {
452 public:
457  FE_DGQLegendre (const unsigned int degree);
458 
464  virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
465  get_constant_modes () const;
466 
473  virtual std::string get_name () const;
474 
475  virtual
476  std::unique_ptr<FiniteElement<dim,spacedim> >
477  clone() const;
478 };
479 
480 
481 
499 template <int dim,int spacedim=dim>
500 class FE_DGQHermite : public FE_DGQ<dim,spacedim>
501 {
502 public:
507  FE_DGQHermite (const unsigned int degree);
508 
515  virtual std::string get_name () const;
516 
517  virtual
518  std::unique_ptr<FiniteElement<dim,spacedim> >
519  clone() const;
520 };
521 
522 
525 DEAL_II_NAMESPACE_CLOSE
526 
527 #endif
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_dgq.cc:795
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_dgq.cc:352
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgq.cc:818
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgq.cc:565
FE_DGQLegendre(const unsigned int degree)
Definition: fe_dgq.cc:835
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgq.cc:867
Definition: fe_dgq.h:105
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_dgq.cc:116
const unsigned int degree
Definition: fe_base.h:311
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_dgq.cc:325
virtual std::size_t memory_consumption() const
Definition: fe_dgq.cc:683
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_dgq.cc:380
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgq.cc:138
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgq.cc:895
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_dgq.cc:670
virtual std::string get_name() const
Definition: fe_dgq.cc:885
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_dgq.cc:151
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_dgq.cc:577
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_dgq.cc:446
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_dgq.cc:242
void rotate_indices(std::vector< unsigned int > &indices, const char direction) const
Definition: fe_dgq.cc:164
FE_DGQArbitraryNodes(const Quadrature< 1 > &points)
Definition: fe_dgq.cc:694
FE_DGQHermite(const unsigned int degree)
Definition: fe_dgq.cc:877
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgq.cc:551
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgq.cc:536
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_dgq.cc:843
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgq.cc:521
virtual std::string get_name() const
Definition: fe_dgq.cc:857
friend class FE_DGQ
Definition: fe_dgq.h:376
virtual bool hp_constraints_are_implemented() const
Definition: fe_dgq.cc:511
virtual std::string get_name() const
Definition: fe_dgq.cc:706
virtual std::string get_name() const
Definition: fe_dgq.cc:98