Dune Core Modules (2.10.0)

Hybrid Utilities

Hybrid utility functions that work on homogeneous as well as heterogeneous containers. More...

Classes

class  Dune::Hybrid::HybridFunctor< Functor >
 Adapter of a hybrid functor that maintains results hybrid. More...
 

Functions

template<class T >
constexpr auto Dune::Hybrid::size (const T &t)
 Size query. More...
 
template<class Container , class Index >
constexpr decltype(auto) Dune::Hybrid::elementAt (Container &&c, Index &&i)
 Get element at given position from container. More...
 
template<class Begin , class End >
constexpr auto Dune::Hybrid::integralRange (const Begin &begin, const End &end)
 Create an integral range. More...
 
template<class End >
constexpr auto Dune::Hybrid::integralRange (const End &end)
 Create an integral range starting from 0. More...
 
template<class Range , class F >
constexpr void Dune::Hybrid::forEach (Range &&range, F &&f)
 Range based for loop. More...
 
template<class Range , class T , class F >
constexpr T Dune::Hybrid::accumulate (Range &&range, T value, F &&f)
 Accumulate values. More...
 
template<class Condition , class IfFunc , class ElseFunc >
decltype(auto) Dune::Hybrid::ifElse (const Condition &condition, IfFunc &&ifFunc, ElseFunc &&elseFunc)
 A conditional expression. More...
 
template<class Condition , class IfFunc >
void Dune::Hybrid::ifElse (const Condition &condition, IfFunc &&ifFunc)
 A conditional expression. More...
 
template<class... Args>
constexpr decltype(auto) Dune::Hybrid::HybridFunctor< Functor >::operator() (const Args &... args) const
 Adapter of a hybrid functor that keeps results hybrid. More...
 
template<class T1 , class T2 >
constexpr auto Dune::Hybrid::equals (T1 &&t1, T2 &&t2)
 Equality comparison. More...
 
template<class Cases , class Value , class Branches , class ElseBranch >
constexpr decltype(auto) Dune::Hybrid::switchCases (const Cases &cases, const Value &value, Branches &&branches, ElseBranch &&elseBranch)
 Switch statement. More...
 
template<class Cases , class Value , class Branches >
constexpr void Dune::Hybrid::switchCases (const Cases &cases, const Value &value, Branches &&branches)
 Switch statement. More...
 
template<class T , class Value , class Branches >
constexpr void Dune::Hybrid::switchCases (IntegralRange< T > range, const Value &value, Branches &&branches)
 Switch statement. More...
 

Variables

constexpr auto Dune::Hybrid::max = hybridFunctor(Impl::Max{})
 Function object that returns the greater of the given values. More...
 
constexpr auto Dune::Hybrid::min = hybridFunctor(Impl::Min{})
 Function object that returns the smaller of the given values. More...
 
constexpr auto Dune::Hybrid::plus = hybridFunctor(std::plus<>{})
 Function object for performing addition. More...
 
constexpr auto Dune::Hybrid::minus = hybridFunctor(std::minus<>{})
 Function object for performing subtraction. More...
 
constexpr auto Dune::Hybrid::equal_to = hybridFunctor(std::equal_to<>{})
 Function object for performing equality comparison. More...
 

Detailed Description

Hybrid utility functions that work on homogeneous as well as heterogeneous containers.

Function Documentation

◆ accumulate()

template<class Range , class T , class F >
constexpr T Dune::Hybrid::accumulate ( Range &&  range,
value,
F &&  f 
)
constexpr

Accumulate values.

Template Parameters
RangeType of given range
TType of accumulated value
FType of binary accumulation operator
Parameters
rangeThe range of values to accumulate
valueInitial value for accumulation
fBinary operator for accumulation

This supports looping over the same ranges as Hybrid::forEach

References Dune::Hybrid::forEach().

Referenced by Dune::Functions::BSplinePreBasis< GV >::BSplinePreBasis(), Dune::BCRSMatrix< B, A >::nonzeroes(), Dune::MultiTypeBlockVector< Args >::one_norm(), Dune::MultiTypeBlockVector< Args >::one_norm_real(), Dune::BDFMCubeLocalBasis< D, R, 2, 1 >::partial(), Dune::BDFMCubeLocalBasis< D, R, 2, 2 >::partial(), Dune::BDFMCubeLocalBasis< D, R, 2, 3 >::partial(), Dune::BDM1Cube2DLocalBasis< D, R >::partial(), Dune::BDM1Simplex2DLocalBasis< D, R >::partial(), Dune::BDM2Cube2DLocalBasis< D, R >::partial(), Dune::BDM2Simplex2DLocalBasis< D, R >::partial(), Dune::HierarchicalSimplexP2LocalBasis< D, R, 2 >::partial(), Dune::RT02DLocalBasis< D, R >::partial(), Dune::RT0Cube2DLocalBasis< D, R >::partial(), Dune::RT12DLocalBasis< D, R >::partial(), Dune::RT1Cube2DLocalBasis< D, R >::partial(), Dune::RT2Cube2DLocalBasis< D, R >::partial(), Dune::RT3Cube2DLocalBasis< D, R >::partial(), Dune::RT4Cube2DLocalBasis< D, R >::partial(), Dune::RefinedP1LocalBasis< D, R, 2 >::partial(), Dune::BDM1Cube3DLocalBasis< D, R >::partial(), Dune::HierarchicalSimplexP2LocalBasis< D, R, 3 >::partial(), Dune::RT03DLocalBasis< D, R >::partial(), Dune::RT0Cube3DLocalBasis< D, R >::partial(), Dune::RT0PrismLocalBasis< D, R >::partial(), Dune::RT0PyramidLocalBasis< D, R >::partial(), Dune::RT1Cube3DLocalBasis< D, R >::partial(), Dune::RefinedP1LocalBasis< D, R, 3 >::partial(), Dune::PowerBasis< Backend, dimR >::partial(), Dune::P0LocalBasis< D, R, d >::partial(), Dune::EdgeS0_5Basis< Geometry, RF >::partial(), Dune::DualP1LocalBasis< D, R, dim, faceDualT >::partial(), Dune::DualQ1LocalBasis< D, R, dim >::partial(), Dune::RefinedP0LocalBasis< D, R, dim >::partial(), Dune::PolynomialBasis< Eval, CM, D, R >::partial(), Dune::Functions::BSplineLocalCoefficients< dim >::size(), and Dune::MultiTypeBlockVector< Args >::two_norm2().

◆ elementAt()

template<class Container , class Index >
constexpr decltype(auto) Dune::Hybrid::elementAt ( Container &&  c,
Index &&  i 
)
constexpr

Get element at given position from container.

Template Parameters
ContainerType of given container
IndexType of index
Parameters
cGiven container
iIndex of element to obtain
Returns
The element at position i, i.e. c[i]

If this returns the i-th entry of c. It supports the following containers

  • Containers providing dynamic access via operator[]
  • Heterogeneous containers providing access via operator[](integral_constant<...>)
  • std::tuple<...>
  • std::integer_sequence

References Dune::Hybrid::elementAt().

Referenced by Dune::Hybrid::elementAt(), Dune::RandomAccessIteratorFacade< T, V, R, D >::operator[](), and Dune::ArrayList< T, N, A >::operator[]().

◆ equals()

template<class T1 , class T2 >
constexpr auto Dune::Hybrid::equals ( T1 &&  t1,
T2 &&  t2 
)
constexpr

Equality comparison.

If both types have a static member value, the result of comparing these is returned as std::integral_constant<bool, *>. Otherwise the result of a runtime comparison of t1 and t2 is directly returned.

Deprecated:

References Dune::Hybrid::equal_to.

Referenced by Dune::BitSetVectorConstReference< block_size, Alloc >::operator!=(), Dune::operator!=(), Dune::BitSetVectorConstReference< block_size, Alloc >::operator==(), and Dune::operator==().

◆ forEach()

template<class Range , class F >
constexpr void Dune::Hybrid::forEach ( Range &&  range,
F &&  f 
)
constexpr

Range based for loop.

Template Parameters
RangeType of given range
FType of given predicate
Parameters
rangeThe range to loop over
fA predicate that will be called with each entry of the range

This supports looping over the following ranges

This especially included instances of std::integer_sequence, std::tuple, Dune::TupleVector, and Dune::MultiTypeBlockVector.

References Dune::Hybrid::forEach().

Referenced by Dune::Hybrid::accumulate(), Dune::MultiTypeBlockVector< Args >::axpy(), Dune::Functions::CompositePreBasis< IMS, SPB >::CompositePreBasis(), Dune::MultiTypeBlockVector< Args >::dim(), Dune::Functions::CompositePreBasis< IMS, SPB >::dimension(), Dune::TypeTree::HybridTreePath< T >::element(), Dune::fillGridViewInfoSerial(), Dune::flatVectorForEach(), Dune::Hybrid::forEach(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::frobenius_norm2(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::infinity_norm(), Dune::MultiTypeBlockVector< Args >::infinity_norm(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::infinity_norm_real(), Dune::MultiTypeBlockVector< Args >::infinity_norm_real(), Dune::Functions::CompositePreBasis< IMS, SPB >::initializeIndices(), Dune::LocalFiniteElementVariantCache< Base >::LocalFiniteElementVariantCache(), Dune::Functions::CompositePreBasis< IMS, SPB >::makeNode(), Dune::Functions::CompositePreBasis< IMS, SPB >::maxNodeSize(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator*=(), Dune::MultiTypeBlockVector< Args >::operator*=(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator+=(), Dune::MultiTypeBlockVector< Args >::operator+=(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator-=(), Dune::MultiTypeBlockVector< Args >::operator-=(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator/=(), Dune::MultiTypeBlockVector< Args >::operator/=(), Dune::operator<<(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::operator=(), Dune::MultiTypeBlockVector< Args >::operator=(), Dune::operator>>(), Dune::TypeTree::HybridTreePath< T >::operator[](), Dune::SizeCache< GridImp >::size(), and Dune::Functions::CompositePreBasis< IMS, SPB >::update().

◆ ifElse() [1/2]

template<class Condition , class IfFunc >
void Dune::Hybrid::ifElse ( const Condition &  condition,
IfFunc &&  ifFunc 
)

A conditional expression.

This provides an ifElse conditional with empty else clause.

References Dune::Hybrid::ifElse().

Referenced by Dune::Hybrid::ifElse().

◆ ifElse() [2/2]

template<class Condition , class IfFunc , class ElseFunc >
decltype(auto) Dune::Hybrid::ifElse ( const Condition &  condition,
IfFunc &&  ifFunc,
ElseFunc &&  elseFunc 
)

A conditional expression.

This will call either ifFunc or elseFunc depending on the condition. In any case a single argument will be passed to the called function. This will always be the identity function. Passing an expression through this function will lead to lazy evaluation. This way both 'branches' can contain expressions that are only valid within this branch if the condition is a std::integral_constant<bool,*>.

In order to do this, the passed functors must have a single argument of type auto.

Due to the lazy evaluation mechanism and support for std::integral_constant<bool,*> this allows to emulate a static if statement.

References Dune::Hybrid::ifElse().

◆ integralRange() [1/2]

template<class Begin , class End >
constexpr auto Dune::Hybrid::integralRange ( const Begin &  begin,
const End &  end 
)
constexpr

◆ integralRange() [2/2]

template<class End >
constexpr auto Dune::Hybrid::integralRange ( const End &  end)
constexpr

Create an integral range starting from 0.

Template Parameters
EndType of end entry of the range
Parameters
endOne past the last entry of the range
Returns
An object encoding the given range

This is a short cut for integralRange(_0, end).

References Dune::Indices::_0, and Dune::Hybrid::integralRange().

Referenced by Dune::Hybrid::integralRange().

◆ operator()()

template<class Functor >
template<class... Args>
constexpr decltype(auto) Dune::Hybrid::HybridFunctor< Functor >::operator() ( const Args &...  args) const
inlineconstexpr

Adapter of a hybrid functor that keeps results hybrid.

Implements an operator that promotes the results of the underlying functor to an integral constant if all the function arguments are integral constants, otherwise, usual promotion rules apply.

◆ size()

template<class T >
constexpr auto Dune::Hybrid::size ( const T &  t)
constexpr

Size query.

Template Parameters
TType of container whose size is queried
Parameters
tContainer whose size is queried
Returns
Size of t

If the size of t is known at compile type the size is returned as std::integral_constant<std::size_t, size>. Otherwise the result of t.size() is returned.

Supported types for deriving the size at compile time are:

  • instances of std::integer_sequence
  • all types std::tuple_size is implemented for
  • all types that have a static constexpr method size() The latter e.g. includes Dune::FieldVector

References Dune::Hybrid::size().

Referenced by Dune::Hybrid::size().

◆ switchCases() [1/3]

template<class Cases , class Value , class Branches >
constexpr void Dune::Hybrid::switchCases ( const Cases &  cases,
const Value &  value,
Branches &&  branches 
)
constexpr

Switch statement.

Template Parameters
CasesType of case range
ValueType of value to check against the cases
BranchesType of branch function
Parameters
casesA range of cases to check for
valueThe value to check against the cases
branchesA callback that will be executed with matching entry from case list

Value is checked against all entries of the given range. If one matches, then branches is executed with the matching value as single argument. If the range is an std::integer_sequence, or StaticIntegralRange, the value is passed as std::integral_constant. If non of the entries matches, the behavior is undefined.

◆ switchCases() [2/3]

template<class Cases , class Value , class Branches , class ElseBranch >
constexpr decltype(auto) Dune::Hybrid::switchCases ( const Cases &  cases,
const Value &  value,
Branches &&  branches,
ElseBranch &&  elseBranch 
)
constexpr

Switch statement.

Template Parameters
CasesType of case range
ValueType of value to check against the cases
BranchesType of branch function
ElseBranchType of branch function
Parameters
casesA range of cases to check for
valueThe value to check against the cases
branchesA callback that will be executed with matching entry from case list
elseBranchA callback that will be executed if no other entry matches

Value is checked against all entries of the given range. If one matches, then branches is executed with the matching value as single argument. If the range is an std::integer_sequence, or StaticIntegralRange, the value is passed as std::integral_constant. If non of the entries matches, then elseBranch is executed without any argument.

Notice that this short circuits, e.g., if one case matches, the others are no longer evaluated.

The return value will be deduced from the else branch.

◆ switchCases() [3/3]

template<class T , class Value , class Branches >
constexpr void Dune::Hybrid::switchCases ( IntegralRange< T >  range,
const Value &  value,
Branches &&  branches 
)
constexpr

Switch statement.

Template Parameters
TThe type of the cases
ValueType of value to check against the cases
BranchesType of branch function
Parameters
rangeA dynamic range of cases to check for
valueThe value to check against the cases
branchesA callback that will be executed with matching entry from case list

This overload of the switchCases utility is selected if the range of cases is passed as an IntegralRange. If the value is contained in that range, it is passed as single argument to the callback branches. If not, the behavior is undefined.

References Dune::IntegralRange< T >::contains().

Variable Documentation

◆ equal_to

constexpr auto Dune::Hybrid::equal_to = hybridFunctor(std::equal_to<>{})
inlineconstexpr

Function object for performing equality comparison.

See also
HybridFunctor

If both arguments have a static member value, the result of comparing these for equality is returned as std::integral_constant<bool, *>. Otherwise the result of a comparison of the two arguments is directly returned.

using namespace Dune::Indices;
{ // hybrid transformation!
auto j = Dune::Hybrid::equal_to( 2, 1); // -> false
auto j = Dune::Hybrid::equal_to( 2, _1); // -> false
auto k = Dune::Hybrid::equal_to(_2, _1); // -> std::false_type
// independent of the context, `k` encodes its value in the type system
}
constexpr index_constant< 1 > _1
Compile time index with value 1.
Definition: indices.hh:55
constexpr index_constant< 2 > _2
Compile time index with value 2.
Definition: indices.hh:58
constexpr auto equal_to
Function object for performing equality comparison.
Definition: hybridutilities.hh:572
Namespace with predefined compile time indices for the range [0,19].
Definition: indices.hh:50

Referenced by Dune::Hybrid::equals().

◆ max

constexpr auto Dune::Hybrid::max = hybridFunctor(Impl::Max{})
inlineconstexpr

Function object that returns the greater of the given values.

See also
HybridFunctor

If all arguments have a static member value, the maximum value of these is returned as std::integral_constant<*, *>. Otherwise the result of a direct max of the arguments is returned.

using namespace Dune::Indices;
{ // hybrid transformation!
auto j = Dune::Hybrid::max( 1, 2); // -> 2
auto j = Dune::Hybrid::max( 1, _2); // -> 2
auto k = Dune::Hybrid::max(_1, _2); // -> Dune::Indices::_2
// independent of the context, `k` encodes its value in the type system
}
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484

Referenced by Dune::Cholmod< Vector, Index >::apply(), checkElementDataMapper(), checkMixedDataMapper(), checkVertexDataMapper(), Dune::Functions::BSplinePreBasis< GV >::evaluate(), Dune::Functions::BSplinePreBasis< GV >::evaluateAll(), Dune::Functions::BSplinePreBasis< GV >::evaluateFunction(), Dune::Functions::BSplinePreBasis< GV >::evaluateJacobian(), Dune::fillGridViewInfoSerial(), Dune::Amg::MatrixHierarchy< M, PI, A >::getCoarsestAggregatesOnFinest(), Dune::Functions::BSplinePreBasis< GV >::indices(), Dune::DenseMatrix< MAT >::infinity_norm(), Dune::DenseVector< V >::infinity_norm(), Dune::BCRSMatrix< B, A >::infinity_norm(), Dune::Matrix< T, A >::infinity_norm(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::infinity_norm(), Dune::MultiTypeBlockVector< Args >::infinity_norm(), Dune::DenseMatrix< MAT >::infinity_norm_real(), Dune::DenseVector< V >::infinity_norm_real(), Dune::BCRSMatrix< B, A >::infinity_norm_real(), Dune::Matrix< T, A >::infinity_norm_real(), Dune::MultiTypeBlockMatrix< FirstRow, Args >::infinity_norm_real(), Dune::MultiTypeBlockVector< Args >::infinity_norm_real(), Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::interior(), Dune::YGridComponent< Coordinates >::intersection(), Dune::Functions::LagrangePreBasis< GV, k, R >::LagrangePreBasis(), Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::LocalCoordinate(), Dune::YaspGrid< dim, Coordinates >::mark(), Dune::IndicesSyncer< T >::DefaultNumberer::operator()(), and Dune::IndicesSyncer< T >::sync().

◆ min

constexpr auto Dune::Hybrid::min = hybridFunctor(Impl::Min{})
inlineconstexpr

Function object that returns the smaller of the given values.

See also
HybridFunctor

If all arguments have a static member value, the minimum value of these is returned as std::integral_constant<*, *>. Otherwise the result of a direct min of the arguments is returned.

using namespace Dune::Indices;
{ // hybrid transformation!
auto j = Dune::Hybrid::min( 1, 2); // -> 1
auto j = Dune::Hybrid::min( 1, _2); // -> 1
auto k = Dune::Hybrid::min(_1, _2); // -> Dune::Indices::_1
// independent of the context, `k` encodes its value in the type system
}
constexpr auto min
Function object that returns the smaller of the given values.
Definition: hybridutilities.hh:506

Referenced by Dune::GeneralizedPCGSolver< X >::apply(), Dune::RestartedFCGSolver< X >::apply(), checkElementDataMapper(), checkMixedDataMapper(), checkVertexDataMapper(), Dune::DiagonalMatrix< K, n >::DiagonalMatrix(), Dune::fillGridViewInfoSerial(), Dune::DenseVector< V >::find(), Dune::MatrixImp::DenseMatrixBase< B, A >::find(), Dune::VariableBlockVector< B, A >::find(), Dune::YGridComponent< Coordinates >::intersection(), Dune::ReservedVector< T, n >::operator<(), Dune::SubsamplingVTKWriter< GridView >::writeGridPoints(), Dune::VTKWriter< GridView >::writeGridPoints(), and Dune::YaspHierarchicIterator< GridImp >::YaspHierarchicIterator().

◆ minus

constexpr auto Dune::Hybrid::minus = hybridFunctor(std::minus<>{})
inlineconstexpr

Function object for performing subtraction.

See also
HybridFunctor

If all arguments have a static member value, the subtracted value of these is returned as std::integral_constant<*, *>. Otherwise the result of a direct subtraction of the arguments is returned.

using namespace Dune::Indices;
{ // hybrid transformation!
auto j = Dune::Hybrid::minus( 2, 1); // -> 1
auto j = Dune::Hybrid::minus( 2, _1); // -> 1
auto k = Dune::Hybrid::minus(_2, _1); // -> Dune::Indices::_1
// independent of the context, `k` encodes its value in the type system
}
constexpr auto minus
Function object for performing subtraction.
Definition: hybridutilities.hh:550

◆ plus

constexpr auto Dune::Hybrid::plus = hybridFunctor(std::plus<>{})
inlineconstexpr

Function object for performing addition.

See also
HybridFunctor

If all arguments have a static member value, the added value of these is returned as std::integral_constant<*, *>. Otherwise the result of a direct addition of the arguments is returned.

using namespace Dune::Indices;
{ // hybrid transformation!
auto j = Dune::Hybrid::plus( 1, 2); // -> 3
auto j = Dune::Hybrid::plus( 1, _2); // -> 3
auto k = Dune::Hybrid::plus(_1, _2); // -> Dune::Indices::_3
// independent of the context, `k` encodes its value in the type system
}
constexpr auto plus
Function object for performing addition.
Definition: hybridutilities.hh:528

Referenced by Dune::TypeTree::accumulate_back(), and Dune::TypeTree::accumulate_front().

Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 26, 23:30, 2024)