Forked from
agnumpde / dune-matrix-vector
39 commits behind the upstream repository.
Code owners
Assign users and groups as approvers for specific file changes. Learn more.
algorithm.hh 1.48 KiB
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_MATRIX_VECTOR_ALGORITHM_HH
#define DUNE_MATRIX_VECTOR_ALGORITHM_HH
#include <dune/common/hybridutilities.hh>
#include <dune/matrix-vector/traitutilities.hh>
namespace Dune {
namespace MatrixVector {
//! \brief Hybrid for loop over sparse range (static/tuple-like candidate)
template <class Range, class F, std::enable_if_t<isTupleOrDerived<Range>(), int> = 0>
void sparseRangeFor(Range&& range, F&& f) {
using namespace Dune::Hybrid;
forEach(integralRange(size(range)), [&](auto&& i) {
f(range[i], i);
});
}
//! \brief Hybrid for loop over sparse range (dynamic/sparse candidate)
template<class Range, class F, std::enable_if_t<not isTupleOrDerived<Range>(), int> = 0>
void sparseRangeFor(Range&& range, F&& f)
{
auto it = range.begin();
auto end = range.end();
for (; it != end; ++it)
f(*it, it.index());
}
//! \brief Hybrid access to first sparse range element (static/tuple-like candiate)
template <class Range, class F, std::enable_if_t<isTupleOrDerived<Range>(), int> = 0>
void sparseRangeFirst(Range&& range, F&& f) {
f(range[Indices::_0]);
}
//! \brief Hybrid access to first sparse range element (dynamic/sparse candiate)
template<class Range, class F, std::enable_if_t<not isTupleOrDerived<Range>(), int> = 0>
void sparseRangeFirst(Range&& range, F&& f)
{
f(*range.begin());
}
} // namespace MatrixVector
} // namespace Dune
#endif // DUNE_MATRIX_VECTOR_ALGORITHM_HH