3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
pq1bubble/elementsolution.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_PQ1BUBBLE_ELEMENT_SOLUTION_HH
25#define DUMUX_PQ1BUBBLE_ELEMENT_SOLUTION_HH
26
27#include <type_traits>
28#include <dune/common/reservedvector.hh>
30
31namespace Dumux {
32
37template<class FVElementGeometry, class PV>
39{
40 using GridGeometry = typename FVElementGeometry::GridGeometry;
41 using GridView = typename GridGeometry::GridView;
42 using Element = typename GridView::template Codim<0>::Entity;
43
44 static constexpr int dim = GridView::dimension;
45 // Dofs are located at corners and in the element center
46 static constexpr int numCubeDofs = (1<<dim) + 1;
47
48public:
50 using PrimaryVariables = PV;
51
54
56 template<class SolutionVector>
57 PQ1BubbleElementSolution(const Element& element, const SolutionVector& sol,
58 const GridGeometry& gridGeometry)
59 {
60 update(element, sol, gridGeometry);
61 }
62
64 template<class ElementVolumeVariables>
65 PQ1BubbleElementSolution(const Element& element, const ElementVolumeVariables& elemVolVars,
66 const FVElementGeometry& fvGeometry)
67 {
68 priVars_.resize(fvGeometry.numScv());
69 for (const auto& scv : scvs(fvGeometry))
70 priVars_[scv.localDofIndex()] = elemVolVars[scv].priVars();
71 }
72
74 template<class SolutionVector>
75 void update(const Element& element, const SolutionVector& sol,
76 const GridGeometry& gridGeometry)
77 {
78 const auto numVertices = element.subEntities(GridView::dimension);
79 priVars_.resize(numVertices+1);
80 for (int vIdx = 0; vIdx < numVertices; ++vIdx)
81 priVars_[vIdx] = sol[gridGeometry.dofMapper().subIndex(element, vIdx, GridView::dimension)];
82
83 priVars_[numVertices] = sol[gridGeometry.dofMapper().index(element)];
84 }
85
87 template<class SolutionVector>
88 void update(const Element& element, const SolutionVector& sol,
89 const FVElementGeometry& fvGeometry)
90 {
91 priVars_.resize(fvGeometry.numScv());
92 for (const auto& scv : scvs(fvGeometry))
93 priVars_[scv.localDofIndex()] = sol[scv.dofIndex()];
94 }
95
97 std::size_t size() const
98 { return priVars_.size(); }
99
101 template<typename IndexType>
102 const PrimaryVariables& operator [](IndexType i) const
103 { return priVars_[i]; }
104
106 template<typename IndexType>
108 { return priVars_[i]; }
109
110private:
111 Dune::ReservedVector<PrimaryVariables, numCubeDofs> priVars_;
112};
113
118template<class Element, class SolutionVector, class GridGeometry>
119auto elementSolution(const Element& element, const SolutionVector& sol, const GridGeometry& gg)
120-> std::enable_if_t<GridGeometry::discMethod == DiscretizationMethods::pq1bubble,
121 PQ1BubbleElementSolution<typename GridGeometry::LocalView,
122 std::decay_t<decltype(std::declval<SolutionVector>()[0])>>
123 >
124{
125 using PrimaryVariables = std::decay_t<decltype(std::declval<SolutionVector>()[0])>;
126 return PQ1BubbleElementSolution<typename GridGeometry::LocalView, PrimaryVariables>(element, sol, gg);
127}
128
133template<class Element, class ElementVolumeVariables, class FVElementGeometry>
134auto elementSolution(const Element& element, const ElementVolumeVariables& elemVolVars, const FVElementGeometry& gg)
135-> std::enable_if_t<FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::pq1bubble,
136 PQ1BubbleElementSolution<FVElementGeometry,
137 typename ElementVolumeVariables::VolumeVariables::PrimaryVariables>>
138{
139 using PrimaryVariables = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables;
140 return PQ1BubbleElementSolution<FVElementGeometry, PrimaryVariables>(element, elemVolVars, gg);
141}
142
143} // end namespace Dumux
144
145#endif
The available discretization methods in Dumux.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr PQ1Bubble pq1bubble
Definition: method.hh:137
The element solution vector.
Definition: pq1bubble/elementsolution.hh:39
PV PrimaryVariables
export the primary variables type
Definition: pq1bubble/elementsolution.hh:50
PQ1BubbleElementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
Constructor with element and solution and grid geometry.
Definition: pq1bubble/elementsolution.hh:57
void update(const Element &element, const SolutionVector &sol, const FVElementGeometry &fvGeometry)
extract the element solution from the solution vector using a local fv geometry
Definition: pq1bubble/elementsolution.hh:88
PQ1BubbleElementSolution(const Element &element, const ElementVolumeVariables &elemVolVars, const FVElementGeometry &fvGeometry)
Constructor with element and elemVolVars and fvGeometry.
Definition: pq1bubble/elementsolution.hh:65
PQ1BubbleElementSolution()=default
Default constructor.
const PrimaryVariables & operator[](IndexType i) const
bracket operator const access
Definition: pq1bubble/elementsolution.hh:102
void update(const Element &element, const SolutionVector &sol, const GridGeometry &gridGeometry)
extract the element solution from the solution vector using a mapper
Definition: pq1bubble/elementsolution.hh:75
std::size_t size() const
return the size of the element solution
Definition: pq1bubble/elementsolution.hh:97