dune-localfunctions 2.10
Loading...
Searching...
No Matches
pq22d.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// 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
8#include <dune/common/fmatrix.hh>
9
11
14
15namespace Dune
16{
17 template<class D, class R>
19 {
22 public:
23 using Traits = typename LFEVariant::Traits;
24
25 PQ22DLocalFiniteElement ( const GeometryType &gt )
26 {
27 if ( gt.isTriangle() )
29 else if ( gt.isQuadrilateral() )
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() )
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
Definition bdfmcube.hh:18
typename Dune::LocalFiniteElementTraits< LocalBasis, LocalCoefficients, LocalInterpolation > Traits
Export LocalFiniteElementTraits.
Definition localfiniteelementvariant.hh:269
unsigned int size() const
Number of shape functions.
Definition localfiniteelementvariant.hh:374
constexpr GeometryType type() const
Number of shape functions.
Definition localfiniteelementvariant.hh:382
const Traits::LocalBasisType & localBasis() const
Provide access to LocalBasis implementation of this LocalFiniteElement.
Definition localfiniteelementvariant.hh:350
const Traits::LocalCoefficientsType & localCoefficients() const
Provide access to LocalCoefficients implementation of this LocalFiniteElement.
Definition localfiniteelementvariant.hh:358
const Traits::LocalInterpolationType & localInterpolation() const
Provide access to LocalInterpolation implementation of this LocalFiniteElement.
Definition localfiniteelementvariant.hh:366
Lagrange finite element for cubes with arbitrary compile-time dimension and polynomial order.
Definition lagrangecube.hh:709
Lagrange finite element for simplices with arbitrary compile-time dimension and polynomial order.
Definition lagrangesimplex.hh:836
Definition pq22d.hh:19
typename LFEVariant::Traits Traits
Definition pq22d.hh:23
const Traits::LocalCoefficientsType & localCoefficients() const
Definition pq22d.hh:46
PQ22DLocalFiniteElement(const GeometryType &gt, const std::vector< unsigned int > vertexmap)
Definition pq22d.hh:33
PQ22DLocalFiniteElement(const GeometryType &gt)
Definition pq22d.hh:25
unsigned int size() const
Number of shape functions in this finite element.
Definition pq22d.hh:57
const Traits::LocalInterpolationType & localInterpolation() const
Definition pq22d.hh:51
GeometryType type() const
Definition pq22d.hh:62
const Traits::LocalBasisType & localBasis() const
Definition pq22d.hh:41