1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 
4 #ifndef DUNE_ISTL_BLOCKLEVEL_HH
5 #define DUNE_ISTL_BLOCKLEVEL_HH
6 
7 #include <algorithm>
8 #include <type_traits>
9 
10 #include <dune/common/indices.hh>
11 #include <dune/common/typetraits.hh>
12 #include <dune/common/hybridutilities.hh>
13 
14 /*!
15  * \file
16  * \brief Helper functions for determining the vector/matrix block level
17  */
18 
19 // forward declaration
20 namespace Dune {
21 template<typename... Args>
22 class MultiTypeBlockVector;
23 template<typename FirstRow, typename... Args>
24 class MultiTypeBlockMatrix;
25 } // end namespace Dune
26 
27 namespace Dune::Impl {
28 
29 // forward declaration
30 template<typename T> struct MaxBlockLevel;
31 template<typename T> struct MinBlockLevel;
32 
33 //! recursively determine block level of a MultiTypeBlockMatrix
34 template<typename M, template<typename B> typename BlockLevel, typename Op>
blockLevelMultiTypeBlockMatrix(const Op & op)35 constexpr std::size_t blockLevelMultiTypeBlockMatrix(const Op& op)
36 {
37   // inialize with zeroth diagonal block
38   using namespace Dune::Indices;
39   using Block00 = typename std::decay_t<decltype(std::declval<M>()[_0][_0])>;
40   std::size_t blockLevel = BlockLevel<Block00>::value() + 1;
41   // iterate over all blocks to determine min/max block level
42   using namespace Dune::Hybrid;
43   forEach(integralRange(index_constant<M::N()>()), [&](auto&& i) {
44     using namespace Dune::Hybrid; // needed for icc, see issue #31
45     forEach(integralRange(index_constant<M::M()>()), [&](auto&& j) {
46       using Block = typename std::decay_t<decltype(std::declval<M>()[i][j])>;
47       blockLevel = op(blockLevel, BlockLevel<Block>::value() + 1);
48     });
49   });
50   return blockLevel;
51 }
52 
53 //! recursively determine block level of a MultiTypeBlockVector
54 template<typename V, template<typename B> typename BlockLevel, typename Op>
blockLevelMultiTypeBlockVector(const Op & op)55 constexpr std::size_t blockLevelMultiTypeBlockVector(const Op& op)
56 {
57   // inialize with zeroth block
58   using namespace Dune::Indices;
59   using Block0 = typename std::decay_t<decltype(std::declval<V>()[_0])>;
60   std::size_t blockLevel = BlockLevel<Block0>::value() + 1;
61   // iterate over all blocks to determine min/max block level
62   using namespace Dune::Hybrid;
63   forEach(integralRange(index_constant<V::size()>()), [&](auto&& i) {
64     using Block = typename std::decay_t<decltype(std::declval<V>()[i])>;
65     blockLevel = op(blockLevel, BlockLevel<Block>::value() + 1);
66   });
67   return blockLevel;
68 }
69 
70 template<typename T>
71 struct MaxBlockLevel
72 {
valueDune::Impl::MaxBlockLevel73   static constexpr std::size_t value(){
74     if constexpr (IsNumber<T>::value)
75       return 0;
76     else
77       return MaxBlockLevel<typename T::block_type>::value() + 1;
78   }
79 };
80 
81 template<typename T>
82 struct MinBlockLevel
83 {
84   // the default implementation assumes minBlockLevel == maxBlockLevel
valueDune::Impl::MinBlockLevel85   static constexpr std::size_t value()
86   { return MaxBlockLevel<T>::value(); }
87 };
88 
89 // max block level for MultiTypeBlockMatrix
90 template<typename FirstRow, typename... Args>
91 struct MaxBlockLevel<MultiTypeBlockMatrix<FirstRow, Args...>>
92 {
valueDune::Impl::MaxBlockLevel93   static constexpr std::size_t value()
94   {
95     using M = MultiTypeBlockMatrix<FirstRow, Args...>;
96     constexpr auto max = [](const auto& a, const auto& b){ return std::max(a,b); };
97     return blockLevelMultiTypeBlockMatrix<M, MaxBlockLevel>(max);
98   }
99 };
100 
101 // min block level for MultiTypeBlockMatrix
102 template<typename FirstRow, typename... Args>
103 struct MinBlockLevel<MultiTypeBlockMatrix<FirstRow, Args...>>
104 {
valueDune::Impl::MinBlockLevel105   static constexpr std::size_t value()
106   {
107     using M = MultiTypeBlockMatrix<FirstRow, Args...>;
108     constexpr auto min = [](const auto& a, const auto& b){ return std::min(a,b); };
109     return blockLevelMultiTypeBlockMatrix<M, MinBlockLevel>(min);
110   }
111 };
112 
113 // max block level for MultiTypeBlockVector
114 template<typename... Args>
115 struct MaxBlockLevel<MultiTypeBlockVector<Args...>>
116 {
valueDune::Impl::MaxBlockLevel117   static constexpr std::size_t value()
118   {
119     using V = MultiTypeBlockVector<Args...>;
120     constexpr auto max = [](const auto& a, const auto& b){ return std::max(a,b); };
121     return blockLevelMultiTypeBlockVector<V, MaxBlockLevel>(max);
122   }
123 };
124 
125 // min block level for MultiTypeBlockVector
126 template<typename... Args>
127 struct MinBlockLevel<MultiTypeBlockVector<Args...>>
128 {
valueDune::Impl::MinBlockLevel129   static constexpr std::size_t value()
130   {
131     using V = MultiTypeBlockVector<Args...>;
132     constexpr auto min = [](const auto& a, const auto& b){ return std::min(a,b); };
133     return blockLevelMultiTypeBlockVector<V, MinBlockLevel>(min);
134   }
135 };
136 
137 // special case: empty MultiTypeBlockVector
138 template<>
139 struct MaxBlockLevel<MultiTypeBlockVector<>>
140 {
141   static constexpr std::size_t value()
142   { return 0; };
143 };
144 
145 // special case: empty MultiTypeBlockVector
146 template<>
147 struct MinBlockLevel<MultiTypeBlockVector<>>
148 {
149   static constexpr std::size_t value()
150   { return 0; };
151 };
152 
153 } // end namespace Dune::Impl
154 
155 namespace Dune {
156 
157 //! Determine the maximum block level of a possibly nested vector/matrix type
158 template<typename T>
maxBlockLevel()159 constexpr std::size_t maxBlockLevel()
160 { return Impl::MaxBlockLevel<T>::value(); }
161 
162 //! Determine the minimum block level of a possibly nested vector/matrix type
163 template<typename T>
minBlockLevel()164 constexpr std::size_t minBlockLevel()
165 { return Impl::MinBlockLevel<T>::value(); }
166 
167 //! Determine if a vector/matrix has a uniquely determinable block level
168 template<typename T>
hasUniqueBlockLevel()169 constexpr bool hasUniqueBlockLevel()
170 { return maxBlockLevel<T>() == minBlockLevel<T>(); }
171 
172 //! Determine the block level of a possibly nested vector/matrix type
173 template<typename T>
blockLevel()174 constexpr std::size_t blockLevel()
175 {
176   static_assert(hasUniqueBlockLevel<T>(), "Block level cannot be uniquely determined!");
177   return Impl::MaxBlockLevel<T>::value();
178 }
179 
180 } // end namespace Dune
181 
182 #endif
183