dune-pdelab
2.7-git
dune
pdelab
finiteelementmap
ap/pk1d.hh
Go to the documentation of this file.
1
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2
// vi: set et ts=4 sw=2 sts=2:
3
4
// Pk in one dimension with k as runtime variable
5
6
#ifndef DUNE_PDELAB_FINITEELEMENTMAP_PK1D_HH
7
#define DUNE_PDELAB_FINITEELEMENTMAP_PK1D_HH
8
9
#include<dune/localfunctions/lagrange.hh>
10
#include<dune/localfunctions/lagrange/equidistantpoints.hh>
11
12
#include<
dune/pdelab/finiteelementmap/finiteelementmap.hh
>
13
14
namespace
Dune
{
15
16
namespace
PDELab {
17
25
template
<
class
D,
class
R>
26
class
Pk1dLocalFiniteElementMap
27
:
public
Dune::PDELab::SimpleLocalFiniteElementMap
<LagrangeLocalFiniteElement<EquidistantPointSet,1,D,R>,1>
28
{
29
public
:
30
31
Pk1dLocalFiniteElementMap
(std::size_t k)
32
:
Dune
::PDELab::
SimpleLocalFiniteElementMap
<LagrangeLocalFiniteElement<EquidistantPointSet,1,D,R>,1>(LagrangeLocalFiniteElement<EquidistantPointSet,1,D,R>(GeometryTypes::cube(1),k))
33
, _k(k)
34
{}
35
36
static
constexpr
bool
fixedSize
()
37
{
38
return
true
;
39
}
40
41
bool
hasDOFs
(
int
codim)
const
42
{
43
switch
(codim)
44
{
45
case
0:
46
return
_k != 1;
47
case
1:
48
return
_k > 0;
49
}
50
return
false
;
51
}
52
53
std::size_t
size
(GeometryType gt)
const
54
{
55
if
(gt.isVertex())
56
return
_k > 0 ? 1 : 0;
57
if
(gt.isLine())
58
return
_k > 0 ? _k - 1 : 1;
59
return
0;
60
}
61
62
std::size_t
maxLocalSize
()
const
63
{
64
return
_k + 1;
65
}
66
67
private
:
68
const
std::size_t _k;
69
};
70
}
71
}
72
#endif // DUNE_PDELAB_FINITEELEMENTMAP_PK1D_HH
Dune::PDELab::SimpleLocalFiniteElementMap
simple implementation where all entities have the same finite element
Definition:
finiteelementmap.hh:98
Dune
For backward compatibility – Do not use this!
Definition:
adaptivity.hh:28
finiteelementmap.hh
Dune::PDELab::Pk1dLocalFiniteElementMap::maxLocalSize
std::size_t maxLocalSize() const
Definition:
ap/pk1d.hh:62
Dune::PDELab::Pk1dLocalFiniteElementMap::fixedSize
static constexpr bool fixedSize()
Definition:
ap/pk1d.hh:36
Dune::PDELab::Pk1dLocalFiniteElementMap::Pk1dLocalFiniteElementMap
Pk1dLocalFiniteElementMap(std::size_t k)
Definition:
ap/pk1d.hh:31
Dune::PDELab::Pk1dLocalFiniteElementMap::size
std::size_t size(GeometryType gt) const
Definition:
ap/pk1d.hh:53
Dune::PDELab::Pk1dLocalFiniteElementMap
FiniteElementMap for the Pk basis in 1d.
Definition:
ap/pk1d.hh:26
Dune::PDELab::Pk1dLocalFiniteElementMap::hasDOFs
bool hasDOFs(int codim) const
Definition:
ap/pk1d.hh:41
Generated by
1.8.17