24#ifndef DUMUX_MULTIDOMAIN_FV_GRID_VARIABLES_HH
25#define DUMUX_MULTIDOMAIN_FV_GRID_VARIABLES_HH
31#include <dune/common/hybridutilities.hh>
32#include <dune/common/indices.hh>
44template<
class MDTraits>
47 using SolutionVector =
typename MDTraits::SolutionVector;
48 static constexpr std::size_t numSubDomains = MDTraits::numSubDomains;
50 template<std::
size_t i>
51 using GridGeometry =
typename MDTraits::template SubDomain<i>::GridGeometry;
52 using GridGeometries =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
54 template<std::
size_t i>
55 using Problem =
typename MDTraits::template SubDomain<i>::Problem;
56 using Problems =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
60 template<std::
size_t i>
61 using Type =
typename MDTraits::template SubDomain<i>::GridVariables;
64 template<std::
size_t i>
68 using TupleType =
typename MDTraits::template Tuple<PtrType>;
77 using namespace Dune::Hybrid;
78 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
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>()
92 : gridVars_(std::move(ggTuple))
96 void init(
const SolutionVector& sol)
98 using namespace Dune::Hybrid;
99 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
101 elementAt(gridVars_,
id)->init(sol[
id]);
106 void update(
const SolutionVector& sol,
bool forceFluxCacheUpdate =
false)
108 using namespace Dune::Hybrid;
109 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
111 elementAt(gridVars_,
id)->update(sol[
id], forceFluxCacheUpdate);
118 using namespace Dune::Hybrid;
119 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
121 elementAt(gridVars_,
id)->updateAfterGridAdaption(sol[
id]);
131 using namespace Dune::Hybrid;
132 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
134 elementAt(gridVars_,
id)->advanceTimeStep();
141 using namespace Dune::Hybrid;
142 forEach(std::make_index_sequence<numSubDomains>{}, [&](
auto&& id)
144 elementAt(gridVars_,
id)->resetTimeStep(sol[
id]);
149 template<std::
size_t i>
151 {
return *std::get<i>(gridVars_); }
154 template<std::
size_t i>
156 {
return *std::get<i>(gridVars_); }
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_); }
164 template<std::
size_t i>
166 {
return std::get<i>(gridVars_); }
172 {
return gridVars_; }
178 {
return gridVars_; }
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.