Doxygen API reference documentation for ideal.II
fe_dg.hh
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2022 - 2023 by the ideal.II authors
4 //
5 // This file is part of the ideal.II library.
6 //
7 // The ideal.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 3.0 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of ideal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef INCLUDE_IDEAL_II_FE_FE_DG_HH_
17 #define INCLUDE_IDEAL_II_FE_FE_DG_HH_
18 
19 #include <deal.II/fe/fe.h>
20 #include <deal.II/fe/fe_dgq.h>
21 
22 #include <memory>
23 
24 namespace idealii::spacetime
25 {
26  template <int dim>
36  {
37  public:
38  // Info/Todo: With proper Quadrature class could be extended to RaudauLeft
39  // and RadauRight Would need Radau based on Legendre polynomial roots
56  {
65  };
69  DG_FiniteElement(std::shared_ptr<dealii::FiniteElement<dim>> fe_space,
70  const unsigned int r,
71  support_type type = support_type::Lobatto);
76  std::shared_ptr<dealii::FiniteElement<dim>>
82  std::shared_ptr<dealii::FiniteElement<1>>
91  const unsigned int dofs_per_cell;
101  type();
102 
103  private:
104  std::shared_ptr<dealii::FiniteElement<dim>> _fe_space;
105  std::shared_ptr<dealii::FiniteElement<1>> _fe_time;
106  support_type _type;
107  };
108 } // namespace idealii::spacetime
109 
110 #endif /* INCLUDE_IDEAL_II_FE_FE_DG_HH_ */
A class for dG elements in time and arbitrary elements in space.
Definition: fe_dg.hh:36
std::shared_ptr< dealii::FiniteElement< dim > > spatial()
The underlying spatial finite element.
support_type type()
The support_type used for construction.
std::shared_ptr< dealii::FiniteElement< 1 > > temporal()
The underlying temporal finite element.
DG_FiniteElement(std::shared_ptr< dealii::FiniteElement< dim >> fe_space, const unsigned int r, support_type type=support_type::Lobatto)
Constructor for the finite element class.
support_type
Choice of underlying temporal support points.
Definition: fe_dg.hh:56
@ Lobatto
Definition: fe_dg.hh:60
@ Legendre
Definition: fe_dg.hh:58
@ RadauLeft
Definition: fe_dg.hh:62
@ RadauRight
Definition: fe_dg.hh:64
const unsigned int dofs_per_cell
The number of degrees of freedom per space-time element.
Definition: fe_dg.hh:91
Namespace for general spacetime object and collections of slab objects.
Definition: idealii.hh:89