3.2-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
34namespace Dumux {
35
41template<class MDTraits>
43{
44 using SolutionVector = typename MDTraits::SolutionVector;
45 static constexpr std::size_t numSubDomains = MDTraits::numSubDomains;
46
47 template<std::size_t i>
48 using GridGeometry = typename MDTraits::template SubDomain<i>::GridGeometry;
49 using GridGeometries = typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
50
51 template<std::size_t i>
52 using Problem = typename MDTraits::template SubDomain<i>::Problem;
53 using Problems = typename MDTraits::template TupleOfSharedPtrConst<Problem>;
54
55public:
57 template<std::size_t i>
58 using Type = typename MDTraits::template SubDomain<i>::GridVariables;
59
61 template<std::size_t i>
62 using PtrType = std::shared_ptr<Type<i>>;
63
65 using TupleType = typename MDTraits::template Tuple<PtrType>;
66
71
77 MultiDomainFVGridVariables(GridGeometries&& gridGeometries, Problems&& problems)
78 {
79 using namespace Dune::Hybrid;
80 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
81 {
82 constexpr auto i = std::decay_t<decltype(id)>::value;
83 elementAt(gridVars_, id) = std::make_shared<Type<i>>( std::get<i>(problems), std::get<i>(gridGeometries));
84 });
85 }
86
88 void init(const SolutionVector& sol)
89 {
90 using namespace Dune::Hybrid;
91 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
92 {
93 elementAt(gridVars_, id)->init(sol[id]);
94 });
95 }
96
98 void update(const SolutionVector& sol, bool forceFluxCacheUpdate = false)
99 {
100 using namespace Dune::Hybrid;
101 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
102 {
103 elementAt(gridVars_, id)->update(sol[id], forceFluxCacheUpdate);
104 });
105 }
106
108 void updateAfterGridAdaption(const SolutionVector& sol)
109 {
110 using namespace Dune::Hybrid;
111 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
112 {
113 elementAt(gridVars_, id)->updateAfterGridAdaption(sol[id]);
114 });
115 }
116
122 {
123 using namespace Dune::Hybrid;
124 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
125 {
126 elementAt(gridVars_, id)->advanceTimeStep();
127 });
128 }
129
131 void resetTimeStep(const SolutionVector& sol)
132 {
133 using namespace Dune::Hybrid;
134 forEach(std::make_index_sequence<numSubDomains>{}, [&](auto&& id)
135 {
136 elementAt(gridVars_, id)->resetTimeStep(sol[id]);
137 });
138 }
139
141 template<std::size_t i>
142 const Type<i>& operator[] (Dune::index_constant<i> id) const
143 { return *Dune::Hybrid::elementAt(gridVars_, id); }
144
146 template<std::size_t i>
147 Type<i>& operator[] (Dune::index_constant<i> id)
148 { return *Dune::Hybrid::elementAt(gridVars_, id); }
149
151 template<std::size_t i>
152 PtrType<i> get(Dune::index_constant<i> id = Dune::index_constant<i>{})
153 { return Dune::Hybrid::elementAt(gridVars_, id); }
154
156 template<std::size_t i>
157 void set(PtrType<i> p, Dune::index_constant<i> id = Dune::index_constant<i>{})
158 { Dune::Hybrid::elementAt(gridVars_, id) = p; }
159
165 { return gridVars_; }
166
167private:
168
170 TupleType gridVars_;
171};
172
173} // end namespace Dumux
174
175#endif
Definition: adapt.hh:29
A multidomain wrapper for multiple grid variables.
Definition: multidomain/fvgridvariables.hh:43
typename MDTraits::template Tuple< PtrType > TupleType
export type of tuple of pointers
Definition: multidomain/fvgridvariables.hh:65
typename MDTraits::template SubDomain< i >::GridVariables Type
export base types of the stored type
Definition: multidomain/fvgridvariables.hh:58
void set(PtrType< i > p, Dune::index_constant< i > id=Dune::index_constant< i >{})
set the pointer for sub domain i
Definition: multidomain/fvgridvariables.hh:157
void update(const SolutionVector &sol, bool forceFluxCacheUpdate=false)
update all variables
Definition: multidomain/fvgridvariables.hh:98
PtrType< i > get(Dune::index_constant< i > id=Dune::index_constant< i >{})
return the grid variables tuple we are wrapping
Definition: multidomain/fvgridvariables.hh:152
const Type< i > & operator[](Dune::index_constant< i > id) const
return the grid variables for domain with index i
Definition: multidomain/fvgridvariables.hh:142
void updateAfterGridAdaption(const SolutionVector &sol)
update all variables after grid adaption
Definition: multidomain/fvgridvariables.hh:108
std::shared_ptr< Type< i > > PtrType
export pointer types the stored type
Definition: multidomain/fvgridvariables.hh:62
MultiDomainFVGridVariables()=default
The default constructor.
void init(const SolutionVector &sol)
initialize all variables
Definition: multidomain/fvgridvariables.hh:88
MultiDomainFVGridVariables(GridGeometries &&gridGeometries, Problems &&problems)
Contruct the grid variables.
Definition: multidomain/fvgridvariables.hh:77
void resetTimeStep(const SolutionVector &sol)
resets state to the one before time integration
Definition: multidomain/fvgridvariables.hh:131
void advanceTimeStep()
Sets the current state as the previous for next time step.
Definition: multidomain/fvgridvariables.hh:121
TupleType getTuple()
return the grid variables tuple we are wrapping
Definition: multidomain/fvgridvariables.hh:164