3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/fvgridvariables.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_MULTIDOMAIN_FV_GRID_VARIABLES_HH
25#define DUMUX_MULTIDOMAIN_FV_GRID_VARIABLES_HH
26
27#include <tuple>
28#include <memory>
29#include <utility>
30
31#include <dune/common/hybridutilities.hh>
32#include <dune/common/indices.hh>
33
36
37namespace Dumux {
38
44template<class MDTraits>
46{
47 using SolutionVector = typename MDTraits::SolutionVector;
48 static constexpr std::size_t numSubDomains = MDTraits::numSubDomains;
49
50 template<std::size_t i>
51 using GridGeometry = typename MDTraits::template SubDomain<i>::GridGeometry;
52 using GridGeometries = typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
53
54 template<std::size_t i>
55 using Problem = typename MDTraits::template SubDomain<i>::Problem;
56 using Problems = typename MDTraits::template TupleOfSharedPtrConst<Problem>;
57
58public:
60 template<std::size_t i>
61 using Type = typename MDTraits::template SubDomain<i>::GridVariables;
62
64 template<std::size_t i>
65 using PtrType = std::shared_ptr<Type<i>>;
66
68 using TupleType = typename MDTraits::template Tuple<PtrType>;
69
76 {
77 using namespace Dune::Hybrid;
78 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
79 {
80 constexpr auto i = std::decay_t<decltype(id)>::value;
81 std::get<i>(gridVars_) = std::make_shared<Type<i>>(
82 problems.template get<i>(), gridGeometries.template get<i>()
83 );
84 });
85 }
86
92 : gridVars_(std::move(ggTuple))
93 {}
94
96 void init(const SolutionVector& sol)
97 {
98 using namespace Dune::Hybrid;
99 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
100 {
101 elementAt(gridVars_, id)->init(sol[id]);
102 });
103 }
104
106 void update(const SolutionVector& sol, bool forceFluxCacheUpdate = false)
107 {
108 using namespace Dune::Hybrid;
109 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
110 {
111 elementAt(gridVars_, id)->update(sol[id], forceFluxCacheUpdate);
112 });
113 }
114
116 void updateAfterGridAdaption(const SolutionVector& sol)
117 {
118 using namespace Dune::Hybrid;
119 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
120 {
121 elementAt(gridVars_, id)->updateAfterGridAdaption(sol[id]);
122 });
123 }
124
130 {
131 using namespace Dune::Hybrid;
132 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
133 {
134 elementAt(gridVars_, id)->advanceTimeStep();
135 });
136 }
137
139 void resetTimeStep(const SolutionVector& sol)
140 {
141 using namespace Dune::Hybrid;
142 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
143 {
144 elementAt(gridVars_, id)->resetTimeStep(sol[id]);
145 });
146 }
147
149 template<std::size_t i>
150 const Type<i>& operator[] (Dune::index_constant<i> id) const
151 { return *std::get<i>(gridVars_); }
152
154 template<std::size_t i>
155 Type<i>& operator[] (Dune::index_constant<i> id)
156 { return *std::get<i>(gridVars_); }
157
159 template<std::size_t i>
160 const PtrType<i>& get(Dune::index_constant<i> id = Dune::index_constant<i>{}) const
161 { return std::get<i>(gridVars_); }
162
164 template<std::size_t i>
165 PtrType<i>& get(Dune::index_constant<i> id = Dune::index_constant<i>{})
166 { return std::get<i>(gridVars_); }
167
172 { return gridVars_; }
173
177 const TupleType& asTuple() const
178 { return gridVars_; }
179
180private:
181
183 TupleType gridVars_;
184};
185
186} // end namespace Dumux
187
188#endif
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
A multidomain wrapper for multiple grid geometries.
Definition: multidomain/fvgridgeometry.hh:58
A multidomain wrapper for multiple grid variables.
Definition: multidomain/fvgridvariables.hh:46
typename MDTraits::template Tuple< PtrType > TupleType
export type of tuple of pointers
Definition: multidomain/fvgridvariables.hh:68
MultiDomainFVGridVariables(TupleType ggTuple)
Construct wrapper from a tuple of grid variables.
Definition: multidomain/fvgridvariables.hh:91
typename MDTraits::template SubDomain< i >::GridVariables Type
export base types of the stored type
Definition: multidomain/fvgridvariables.hh:61
const TupleType & asTuple() const
Access the underlying tuple representation.
Definition: multidomain/fvgridvariables.hh:177
void update(const SolutionVector &sol, bool forceFluxCacheUpdate=false)
update all variables
Definition: multidomain/fvgridvariables.hh:106
MultiDomainFVGridVariables(MultiDomainFVGridGeometry< MDTraits > gridGeometries, MultiDomainFVProblem< MDTraits > problems)
Construct the grid variables.
Definition: multidomain/fvgridvariables.hh:75
const Type< i > & operator[](Dune::index_constant< i > id) const
return the grid variables for domain with index i
Definition: multidomain/fvgridvariables.hh:150
void updateAfterGridAdaption(const SolutionVector &sol)
update all variables after grid adaption
Definition: multidomain/fvgridvariables.hh:116
const PtrType< i > & get(Dune::index_constant< i > id=Dune::index_constant< i >{}) const
access the ith grid variables pointer we are wrapping
Definition: multidomain/fvgridvariables.hh:160
std::shared_ptr< Type< i > > PtrType
export pointer types the stored type
Definition: multidomain/fvgridvariables.hh:65
void init(const SolutionVector &sol)
initialize all variables
Definition: multidomain/fvgridvariables.hh:96
void resetTimeStep(const SolutionVector &sol)
resets state to the one before time integration
Definition: multidomain/fvgridvariables.hh:139
PtrType< i > & get(Dune::index_constant< i > id=Dune::index_constant< i >{})
access the ith grid variables pointer we are wrapping
Definition: multidomain/fvgridvariables.hh:165
void advanceTimeStep()
Sets the current state as the previous for next time step.
Definition: multidomain/fvgridvariables.hh:129
TupleType & asTuple()
Access the underlying tuple representation.
Definition: multidomain/fvgridvariables.hh:171
A multidomain wrapper for multiple problems.
Definition: multidomain/fvproblem.hh:45
Multidomain wrapper for multiple problems.
Multidomain wrapper for multiple grid geometries.