3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_FECACHE_HH
25#define DUMUX_DISCRETIZATION_PQ1BUBBLE_FECACHE_HH
26
27#include <memory>
28
29#include <dune/common/exceptions.hh>
30#include <dune/geometry/type.hh>
31
32#include <dune/localfunctions/common/virtualinterface.hh>
33#include <dune/localfunctions/common/virtualwrappers.hh>
34
36
37namespace Dumux {
38
39template< class CoordScalar, class Scalar, unsigned int dim>
41{
42 static_assert(dim == 2 || dim == 3, "P1/Q1 bubble FE spaces only implemented for 2D and 3D grids");
43
44 // These are so-called non-conforming finite element spaces
45 // the local basis is only continuous at given points on the faces
46 using P1Bubble = Dumux::PQ1BubbleLocalFiniteElement<CoordScalar, Scalar, dim, Dune::GeometryTypes::simplex(dim).toId()>;
47 using Q1Bubble = Dumux::PQ1BubbleLocalFiniteElement<CoordScalar, Scalar, dim, Dune::GeometryTypes::cube(dim).toId()>;
48
49public:
50 using FiniteElementType = Dune::LocalFiniteElementVirtualInterface<typename P1Bubble::Traits::LocalBasisType::Traits>;
51
53 : p1BubbleBasis_(std::make_unique<Dune::LocalFiniteElementVirtualImp<P1Bubble>>(P1Bubble{}))
54 , q1BubbleBasis_(std::make_unique<Dune::LocalFiniteElementVirtualImp<Q1Bubble>>(Q1Bubble{}))
55 {}
56
58 const FiniteElementType& get(const Dune::GeometryType& gt) const
59 {
60 if (gt.isSimplex())
61 return *p1BubbleBasis_;
62 if (gt.isCube())
63 return *q1BubbleBasis_;
64 else
65 DUNE_THROW(Dune::NotImplemented,
66 "Lagrange bubble local finite element for geometry type " << gt
67 );
68 }
69
70private:
71 std::unique_ptr<FiniteElementType> p1BubbleBasis_;
72 std::unique_ptr<FiniteElementType> q1BubbleBasis_;
73};
74
75} // end namespace Dumux
76
77#endif
Evaluate P1/Q1 basis with bubble function.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
Definition: pq1bubblefecache.hh:41
const FiniteElementType & get(const Dune::GeometryType &gt) const
Get local finite element for given GeometryType.
Definition: pq1bubblefecache.hh:58
Dune::LocalFiniteElementVirtualInterface< typename P1Bubble::Traits::LocalBasisType::Traits > FiniteElementType
Definition: pq1bubblefecache.hh:50
PQ1BubbleFECache()
Definition: pq1bubblefecache.hh:52
P1/Q1 + Bubble finite element.
Definition: pq1bubblelocalfiniteelement.hh:297