3.2-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#include <dune/common/version.hh>
32
33namespace Dumux {
34
35#if DUNE_VERSION_GT_REV(DUNE_ISTL,2,6,0)
41template<class ...Args, std::size_t ...i>
42auto partial(Dune::MultiTypeBlockVector<Args...>& v, Dune::index_constant<i>... indices)
43{
44 return Dune::MultiTypeBlockVector<std::add_lvalue_reference_t<std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(v[indices]...);
45}
46
47#else
48
49// for backwards-compatibility of partial with v2.6.0
50template<typename ... Args>
51class MultiTypeBlockVectorProxy : public std::tuple<Args...>
52{
53 typedef std::tuple<Args...> ParentType;
54public:
55 using ParentType::ParentType;
56
58 typedef double field_type;
59
62 template< std::size_t index >
63 typename std::tuple_element<index, ParentType>::type&
64 operator[] ( const std::integral_constant< std::size_t, index > indexVariable )
65 {
66 DUNE_UNUSED_PARAMETER(indexVariable);
67 return std::get<index>(*this);
68 }
69
72 template< std::size_t index >
73 const typename std::tuple_element<index, ParentType>::type&
74 operator[] ( const std::integral_constant< std::size_t, index > indexVariable ) const
75 {
76 DUNE_UNUSED_PARAMETER(indexVariable);
77 return std::get<index>(*this);
78 }
79};
80
81template<class T> struct isMultiTypeBlockVector;
82
83template<class... T>
84struct isMultiTypeBlockVector<MultiTypeBlockVectorProxy<T...> > : public std::true_type {};
85
91template<class ...Args, std::size_t ...i>
92auto partial(Dune::MultiTypeBlockVector<Args...>& v, Dune::index_constant<i>... indices)
93{
94 return MultiTypeBlockVectorProxy<std::add_lvalue_reference_t<std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(v[indices]...);
95}
96
97#endif
98
104template<class ...Args, std::size_t ...i>
105auto partial(std::tuple<Args...>& v, Dune::index_constant<i>... indices)
106{
107 return std::tuple<std::add_lvalue_reference_t<std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(std::get<indices>(v)...);
108}
109
115template<class T, std::size_t ...i>
116auto partial(T& t, std::tuple<Dune::index_constant<i>...> indices)
117{
118 return partial(t, Dune::index_constant<i>{}...);
119}
120
121} // end namespace Dumux
122
123#endif
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:92
Definition: partial.hh:52
std::tuple_element< index, ParentType >::type & operator[](const std::integral_constant< std::size_t, index > indexVariable)
Random-access operator.
Definition: partial.hh:64
double field_type
Definition: partial.hh:58
MultiTypeBlockVectorProxy< Args... > type
Definition: partial.hh:57
Helper type to determine whether a given type is a Dune::MultiTypeBlockVector.
Definition: vector.hh:34