3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
partial.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_COMMON_PARTIAL_HH
25#define DUMUX_COMMON_PARTIAL_HH
26
27#include <tuple>
28#include <type_traits>
29
30#include <dune/istl/multitypeblockvector.hh>
31
32namespace Dumux {
33
39template<class ...Args, std::size_t ...i>
40auto partial(Dune::MultiTypeBlockVector<Args...>& v, Dune::index_constant<i>... indices)
41{
42 return Dune::MultiTypeBlockVector<std::add_lvalue_reference_t<std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(v[indices]...);
43}
44
50template<class ...Args, std::size_t ...i>
51auto partial(const Dune::MultiTypeBlockVector<Args...>& v, Dune::index_constant<i>... indices)
52{
53 return Dune::MultiTypeBlockVector<std::add_lvalue_reference_t<const std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(v[indices]...);
54}
55
61template<class ...Args, std::size_t ...i>
62auto partial(std::tuple<Args...>& v, Dune::index_constant<i>... indices)
63{
64 return std::tuple<std::add_lvalue_reference_t<std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(std::get<indices>(v)...);
65}
66
72template<class ...Args, std::size_t ...i>
73auto partial(const std::tuple<Args...>& v, Dune::index_constant<i>... indices)
74{
75 return std::tuple<std::add_lvalue_reference_t<const std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(std::get<indices>(v)...);
76}
77
83template<class T, std::size_t ...i>
84auto partial(T& t, std::tuple<Dune::index_constant<i>...> indices)
85{
86 return partial(t, Dune::index_constant<i>{}...);
87}
88
89} // end namespace Dumux
90
91#endif
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
auto partial(Dune::MultiTypeBlockVector< Args... > &v, Dune::index_constant< i >... indices)
a function to get a MultiTypeBlockVector with references to some entries of another MultiTypeBlockVec...
Definition: partial.hh:40
Definition: variablesbackend.hh:43