version 3.8
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-FileCopyrightInfo: 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>
29{
30 static_assert(dim == 2 || dim == 3, "P1/Q1 bubble FE spaces only implemented for 2D and 3D grids");
31
32 // These are so-called non-conforming finite element spaces
33 // the local basis is only continuous at given points on the faces
34 using P1Bubble = Dumux::PQ1BubbleLocalFiniteElement<CoordScalar, Scalar, dim, Dune::GeometryTypes::simplex(dim).toId()>;
35 using Q1Bubble = Dumux::PQ1BubbleLocalFiniteElement<CoordScalar, Scalar, dim, Dune::GeometryTypes::cube(dim).toId()>;
36
37public:
38 using FiniteElementType = Dune::LocalFiniteElementVirtualInterface<typename P1Bubble::Traits::LocalBasisType::Traits>;
39
41 : p1BubbleBasis_(std::make_unique<Dune::LocalFiniteElementVirtualImp<P1Bubble>>(P1Bubble{}))
42 , q1BubbleBasis_(std::make_unique<Dune::LocalFiniteElementVirtualImp<Q1Bubble>>(Q1Bubble{}))
43 {}
44
46 const FiniteElementType& get(const Dune::GeometryType& gt) const
47 {
48 if (gt.isSimplex())
49 return *p1BubbleBasis_;
50 if (gt.isCube())
51 return *q1BubbleBasis_;
52 else
53 DUNE_THROW(Dune::NotImplemented,
54 "Lagrange bubble local finite element for geometry type " << gt
55 );
56 }
57
58private:
59 std::unique_ptr<FiniteElementType> p1BubbleBasis_;
60 std::unique_ptr<FiniteElementType> q1BubbleBasis_;
61};
62
63} // end namespace Dumux
64
65#endif
Definition: pq1bubblefecache.hh:29
const FiniteElementType & get(const Dune::GeometryType &gt) const
Get local finite element for given GeometryType.
Definition: pq1bubblefecache.hh:46
Dune::LocalFiniteElementVirtualInterface< typename P1Bubble::Traits::LocalBasisType::Traits > FiniteElementType
Definition: pq1bubblefecache.hh:38
PQ1BubbleFECache()
Definition: pq1bubblefecache.hh:40
P1/Q1 + Bubble finite element.
Definition: pq1bubblelocalfiniteelement.hh:285
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
Evaluate P1/Q1 basis with bubble function.