version 3.11-dev
pq1bubblefecache.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_FECACHE_HH
13#define DUMUX_DISCRETIZATION_PQ1BUBBLE_FECACHE_HH
14
15#include <memory>
16
17#include <dune/common/exceptions.hh>
18#include <dune/geometry/type.hh>
19
20#include <dune/localfunctions/common/virtualinterface.hh>
21#include <dune/localfunctions/common/virtualwrappers.hh>
22
24
25namespace Dumux {
26
27template< class CoordScalar, class Scalar, unsigned int dim, std::size_t numCubeBubbleDofs = 1>
29{
30 static_assert(dim == 2 || dim == 3, "P1/Q1 bubble FE spaces only implemented for 2D and 3D grids");
31 static_assert(numCubeBubbleDofs == 1 || numCubeBubbleDofs == 2,
32 "P1/Q1 bubble FE spaces supports numCubeBubbleDofs = 1 or 2");
33
34 // These are so-called non-conforming finite element spaces
35 // the local basis is only continuous at given points on the faces
36 using P1Bubble = Dumux::PQ1BubbleLocalFiniteElement<CoordScalar, Scalar, dim, Dune::GeometryTypes::simplex(dim).toId(), 1>;
37 using Q1Bubble = Dumux::PQ1BubbleLocalFiniteElement<CoordScalar, Scalar, dim, Dune::GeometryTypes::cube(dim).toId(), numCubeBubbleDofs>;
38
39public:
40 using FiniteElementType = Dune::LocalFiniteElementVirtualInterface<typename P1Bubble::Traits::LocalBasisType::Traits>;
41
43 : p1BubbleBasis_(std::make_unique<Dune::LocalFiniteElementVirtualImp<P1Bubble>>(P1Bubble{}))
44 , q1BubbleBasis_(std::make_unique<Dune::LocalFiniteElementVirtualImp<Q1Bubble>>(Q1Bubble{}))
45 {}
46
48 const FiniteElementType& get(const Dune::GeometryType& gt) const
49 {
50 if (gt.isSimplex())
51 return *p1BubbleBasis_;
52 if (gt.isCube())
53 return *q1BubbleBasis_;
54 else
55 DUNE_THROW(Dune::NotImplemented,
56 "Lagrange bubble local finite element for geometry type " << gt
57 );
58 }
59
60private:
61 std::unique_ptr<FiniteElementType> p1BubbleBasis_;
62 std::unique_ptr<FiniteElementType> q1BubbleBasis_;
63};
64
65} // end namespace Dumux
66
67#endif
Definition: pq1bubblefecache.hh:29
PQ1BubbleFECache()
Definition: pq1bubblefecache.hh:42
const FiniteElementType & get(const Dune::GeometryType &gt) const
Get local finite element for given GeometryType.
Definition: pq1bubblefecache.hh:48
Dune::LocalFiniteElementVirtualInterface< typename P1Bubble::Traits::LocalBasisType::Traits > FiniteElementType
Definition: pq1bubblefecache.hh:40
P1/Q1 + Bubble finite element.
Definition: pq1bubblelocalfiniteelement.hh:322
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
Evaluate P1/Q1 basis with bubble function.