Dune Core Modules (unstable)

pq22d.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_PQ22DLOCALFINITEELEMENT_HH
6#define DUNE_PQ22DLOCALFINITEELEMENT_HH
7
9
10#include <dune/localfunctions/common/localfiniteelementvariant.hh>
11
12#include <dune/localfunctions/lagrange/lagrangesimplex.hh>
13#include <dune/localfunctions/lagrange/lagrangecube.hh>
14
15namespace Dune
16{
17 template<class D, class R>
18 class PQ22DLocalFiniteElement
19 {
20 using LFEVariant = LocalFiniteElementVariant<LagrangeSimplexLocalFiniteElement<D,R,2,2>,
21 LagrangeCubeLocalFiniteElement<D,R,2,2> >;
22 public:
23 using Traits = typename LFEVariant::Traits;
24
25 PQ22DLocalFiniteElement ( const GeometryType &gt )
26 {
27 if ( gt.isTriangle() )
28 lfeVariant_ = LagrangeSimplexLocalFiniteElement<D,R,2,2>();
29 else if ( gt.isQuadrilateral() )
30 lfeVariant_ = LagrangeCubeLocalFiniteElement<D,R,2,2>();
31 }
32
33 PQ22DLocalFiniteElement ( const GeometryType &gt, const std::vector<unsigned int> vertexmap )
34 {
35 if ( gt.isTriangle() )
36 lfeVariant_ = LagrangeSimplexLocalFiniteElement<D,R,2,2>(vertexmap);
37 else if ( gt.isQuadrilateral() )
38 lfeVariant_ = LagrangeCubeLocalFiniteElement<D,R,2,2>();
39 }
40
41 const typename Traits::LocalBasisType& localBasis () const
42 {
43 return lfeVariant_.localBasis();
44 }
45
46 const typename Traits::LocalCoefficientsType& localCoefficients () const
47 {
48 return lfeVariant_.localCoefficients();
49 }
50
51 const typename Traits::LocalInterpolationType& localInterpolation () const
52 {
53 return lfeVariant_.localInterpolation();
54 }
55
57 unsigned int size () const
58 {
59 return lfeVariant_.size();
60 }
61
62 GeometryType type () const
63 {
64 return lfeVariant_.type();
65 }
66
67 private:
68
69 LFEVariant lfeVariant_;
70 };
71
72}
73
74#endif
typename Dune::LocalFiniteElementTraits< LocalBasis, LocalCoefficients, LocalInterpolation > Traits
Export LocalFiniteElementTraits.
Definition: localfiniteelementvariant.hh:269
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
Implements a matrix constructed from a given type representing a field and compile-time given number ...
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)