Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-istl
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Model registry
Operate
Environments
Monitor
Incidents
Service Desk
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Patrick Jaap
dune-istl
Commits
9d4726b2
Commit
9d4726b2
authored
4 years ago
by
Timo Koch
Browse files
Options
Downloads
Patches
Plain Diff
[blocklevel] Make blockLevel work for vectors too
parent
8cb65bf9
No related branches found
Branches containing commit
No related tags found
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
dune/istl/blocklevel.hh
+64
-5
64 additions, 5 deletions
dune/istl/blocklevel.hh
with
64 additions
and
5 deletions
dune/istl/blocklevel.hh
+
64
−
5
View file @
9d4726b2
...
...
@@ -13,11 +13,13 @@
/*!
* \file
* \brief Helper functions for determining the matrix block level
* \brief Helper functions for determining the
vector/
matrix block level
*/
// forward declaration
namespace
Dune
{
template
<
typename
...
Args
>
class
MultiTypeBlockVector
;
template
<
typename
FirstRow
,
typename
...
Args
>
class
MultiTypeBlockMatrix
;
}
// end namespace Dune
...
...
@@ -48,6 +50,23 @@ constexpr std::size_t blockLevelMultiTypeBlockMatrix(const Op& op)
return
blockLevel
;
}
//! recursively determine block level of a MultiTypeBlockVector
template
<
typename
V
,
template
<
typename
B
>
typename
BlockLevel
,
typename
Op
>
constexpr
std
::
size_t
blockLevelMultiTypeBlockVector
(
const
Op
&
op
)
{
// inialize with zeroth block
using
namespace
Dune
::
Indices
;
using
Block0
=
typename
std
::
decay_t
<
decltype
(
std
::
declval
<
V
>
()[
_0
])
>
;
std
::
size_t
blockLevel
=
BlockLevel
<
Block0
>::
value
()
+
1
;
// iterate over all blocks to determine min/max block level
using
namespace
Dune
::
Hybrid
;
forEach
(
integralRange
(
index_constant
<
V
::
size
()
>
()),
[
&
](
auto
&&
i
)
{
using
Block
=
typename
std
::
decay_t
<
decltype
(
std
::
declval
<
V
>
()[
i
])
>
;
blockLevel
=
op
(
blockLevel
,
BlockLevel
<
Block
>::
value
()
+
1
);
});
return
blockLevel
;
}
template
<
typename
T
>
struct
MaxBlockLevel
{
...
...
@@ -91,26 +110,66 @@ struct MinBlockLevel<MultiTypeBlockMatrix<FirstRow, Args...>>
}
};
// max block level for MultiTypeBlockVector
template
<
typename
...
Args
>
struct
MaxBlockLevel
<
MultiTypeBlockVector
<
Args
...
>>
{
static
constexpr
std
::
size_t
value
()
{
using
V
=
MultiTypeBlockVector
<
Args
...
>
;
constexpr
auto
max
=
[](
const
auto
&
a
,
const
auto
&
b
){
return
std
::
max
(
a
,
b
);
};
return
blockLevelMultiTypeBlockVector
<
V
,
MaxBlockLevel
>
(
max
);
}
};
// min block level for MultiTypeBlockVector
template
<
typename
...
Args
>
struct
MinBlockLevel
<
MultiTypeBlockVector
<
Args
...
>>
{
static
constexpr
std
::
size_t
value
()
{
using
V
=
MultiTypeBlockVector
<
Args
...
>
;
constexpr
auto
min
=
[](
const
auto
&
a
,
const
auto
&
b
){
return
std
::
min
(
a
,
b
);
};
return
blockLevelMultiTypeBlockVector
<
V
,
MinBlockLevel
>
(
min
);
}
};
// special case: empty MultiTypeBlockVector
template
<
>
struct
MaxBlockLevel
<
MultiTypeBlockVector
<>>
{
static
constexpr
std
::
size_t
value
()
{
return
0
;
};
};
// special case: empty MultiTypeBlockVector
template
<
>
struct
MinBlockLevel
<
MultiTypeBlockVector
<>>
{
static
constexpr
std
::
size_t
value
()
{
return
0
;
};
};
}
// end namespace Dune::Impl
namespace
Dune
{
//! Determine the maximum block level of a possibly nested matrix type
//! Determine the maximum block level of a possibly nested
vector/
matrix type
template
<
typename
T
>
constexpr
std
::
size_t
maxBlockLevel
()
{
return
Impl
::
MaxBlockLevel
<
T
>::
value
();
}
//! Determine the minimum block level of a possibly nested matrix type
//! Determine the minimum block level of a possibly nested
vector/
matrix type
template
<
typename
T
>
constexpr
std
::
size_t
minBlockLevel
()
{
return
Impl
::
MinBlockLevel
<
T
>::
value
();
}
//! Determine if a matrix has a uniquely determinable block level
//! Determine if a
vector/
matrix has a uniquely determinable block level
template
<
typename
T
>
constexpr
bool
hasUniqueBlockLevel
()
{
return
maxBlockLevel
<
T
>
()
==
minBlockLevel
<
T
>
();
}
//! Determine the block level of a possibly nested matrix type
//! Determine the block level of a possibly nested
vector/
matrix type
template
<
typename
T
>
constexpr
std
::
size_t
blockLevel
()
{
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment