12#ifndef DUMUX_MULTIDOMAIN_FV_GRID_VARIABLES_HH
13#define DUMUX_MULTIDOMAIN_FV_GRID_VARIABLES_HH
19#include <dune/common/hybridutilities.hh>
20#include <dune/common/indices.hh>
32template<
class MDTraits>
35 using SolutionVector =
typename MDTraits::SolutionVector;
36 static constexpr std::size_t numSubDomains = MDTraits::numSubDomains;
38 template<std::
size_t i>
39 using GridGeometry =
typename MDTraits::template SubDomain<i>::GridGeometry;
40 using GridGeometries =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
42 template<std::
size_t i>
43 using Problem =
typename MDTraits::template SubDomain<i>::Problem;
44 using Problems =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
48 template<std::
size_t i>
49 using Type =
typename MDTraits::template SubDomain<i>::GridVariables;
52 template<std::
size_t i>
56 using TupleType =
typename MDTraits::template Tuple<PtrType>;
65 using namespace Dune::Hybrid;
66 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
68 constexpr auto i = std::decay_t<
decltype(id)>::value;
69 std::get<i>(gridVars_) = std::make_shared<Type<i>>(
70 problems.template get<i>(), gridGeometries.template get<i>()
80 : gridVars_(std::move(ggTuple))
84 void init(
const SolutionVector& sol)
86 using namespace Dune::Hybrid;
87 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
89 elementAt(gridVars_,
id)->init(sol[
id]);
94 void update(
const SolutionVector& sol,
bool forceFluxCacheUpdate =
false)
96 using namespace Dune::Hybrid;
97 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
99 elementAt(gridVars_,
id)->update(sol[
id], forceFluxCacheUpdate);
106 using namespace Dune::Hybrid;
107 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
109 elementAt(gridVars_,
id)->updateAfterGridAdaption(sol[
id]);
119 using namespace Dune::Hybrid;
120 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
122 elementAt(gridVars_,
id)->advanceTimeStep();
129 using namespace Dune::Hybrid;
130 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
132 elementAt(gridVars_,
id)->resetTimeStep(sol[
id]);
137 template<std::
size_t i>
139 {
return *std::get<i>(gridVars_); }
142 template<std::
size_t i>
144 {
return *std::get<i>(gridVars_); }
147 template<std::
size_t i>
148 const PtrType<i>&
get(Dune::index_constant<i>
id = Dune::index_constant<i>{})
const
149 {
return std::get<i>(gridVars_); }
152 template<std::
size_t i>
154 {
return std::get<i>(gridVars_); }
160 {
return gridVars_; }
166 {
return gridVars_; }
A multidomain wrapper for multiple grid geometries.
Definition: multidomain/fvgridgeometry.hh:46
A multidomain wrapper for multiple grid variables.
Definition: multidomain/fvgridvariables.hh:34
typename MDTraits::template Tuple< PtrType > TupleType
export type of tuple of pointers
Definition: multidomain/fvgridvariables.hh:56
MultiDomainFVGridVariables(TupleType ggTuple)
Construct wrapper from a tuple of grid variables.
Definition: multidomain/fvgridvariables.hh:79
typename MDTraits::template SubDomain< i >::GridVariables Type
export base types of the stored type
Definition: multidomain/fvgridvariables.hh:49
const TupleType & asTuple() const
Access the underlying tuple representation.
Definition: multidomain/fvgridvariables.hh:165
void update(const SolutionVector &sol, bool forceFluxCacheUpdate=false)
update all variables
Definition: multidomain/fvgridvariables.hh:94
MultiDomainFVGridVariables(MultiDomainFVGridGeometry< MDTraits > gridGeometries, MultiDomainFVProblem< MDTraits > problems)
Construct the grid variables.
Definition: multidomain/fvgridvariables.hh:63
const Type< i > & operator[](Dune::index_constant< i > id) const
return the grid variables for domain with index i
Definition: multidomain/fvgridvariables.hh:138
void updateAfterGridAdaption(const SolutionVector &sol)
update all variables after grid adaption
Definition: multidomain/fvgridvariables.hh:104
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:148
std::shared_ptr< Type< i > > PtrType
export pointer types the stored type
Definition: multidomain/fvgridvariables.hh:53
void init(const SolutionVector &sol)
initialize all variables
Definition: multidomain/fvgridvariables.hh:84
void resetTimeStep(const SolutionVector &sol)
resets state to the one before time integration
Definition: multidomain/fvgridvariables.hh:127
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:153
void advanceTimeStep()
Sets the current state as the previous for next time step.
Definition: multidomain/fvgridvariables.hh:117
TupleType & asTuple()
Access the underlying tuple representation.
Definition: multidomain/fvgridvariables.hh:159
A multidomain wrapper for multiple problems.
Definition: multidomain/fvproblem.hh:33
Multidomain wrapper for multiple grid geometries.
Multidomain wrapper for multiple problems.