![]() |
Block-Structured AMR Software Framework
|
Performance-portable reduction and scan operations in AMReX. More...
Classes | |
| struct | amrex::ReduceOpSum |
| struct | amrex::ReduceOpMin |
| struct | amrex::ReduceOpMax |
| struct | amrex::ReduceOpLogicalAnd |
| struct | amrex::ReduceOpLogicalOr |
| class | amrex::ReduceData< Ts > |
| class | amrex::ReduceOps< Ps > |
| class | amrex::Reducer< Ops, Ts > |
| Class for local reductions (e.g., sum, min and max). More... | |
Functions | |
| template<typename T , typename N , std::enable_if_t< std::is_integral_v< N >, int > = 0> | |
| T | amrex::Reduce::Sum (N n, T const *v, T init_val=0) |
| Compute the sum of an array of values. | |
| template<typename T , typename N , typename F , std::enable_if_t< std::is_integral_v< N > &&!std::is_same_v< T *, std::decay_t< F > >, int > = 0> | |
| T | amrex::Reduce::Sum (N n, F const &f, T init_val=0) |
| Compute the sum of values generated by a callable. | |
| template<typename T , typename N , std::enable_if_t< std::is_integral_v< N >, int > = 0> | |
| T | amrex::Reduce::Min (N n, T const *v, T init_val=std::numeric_limits< T >::max()) |
| Compute the minimum of an array of values. | |
| template<typename T , typename N , typename F , std::enable_if_t< std::is_integral_v< N > &&!std::is_same_v< T *, std::decay_t< F > >, int > = 0> | |
| T | amrex::Reduce::Min (N n, F const &f, T init_val=std::numeric_limits< T >::max()) |
| Compute the minimum of values generated by a callable. | |
| template<typename T , typename N , std::enable_if_t< std::is_integral_v< N >, int > = 0> | |
| T | amrex::Reduce::Max (N n, T const *v, T init_val=std::numeric_limits< T >::lowest()) |
| Compute the maximum of an array of values. | |
| template<typename T , typename N , typename F , std::enable_if_t< std::is_integral_v< N > &&!std::is_same_v< T *, std::decay_t< F > >, int > = 0> | |
| T | amrex::Reduce::Max (N n, F const &f, T init_val=std::numeric_limits< T >::lowest()) |
| Compute the maximum of values generated by a callable. | |
| template<typename T , typename N , std::enable_if_t< std::is_integral_v< N >, int > = 0> | |
| std::pair< T, T > | amrex::Reduce::MinMax (N n, T const *v) |
| Compute the minimum and maximum of an array of values. | |
| template<typename T , typename N , typename F , std::enable_if_t< std::is_integral_v< N > &&!std::is_same_v< T *, std::decay_t< F > >, int > = 0> | |
| std::pair< T, T > | amrex::Reduce::MinMax (N n, F const &f) |
| Compute the minimum and maximum of values generated by a callable. | |
| template<typename T , typename N , typename P , std::enable_if_t< std::is_integral_v< N >, int > = 0> | |
| bool | amrex::Reduce::AnyOf (N n, T const *v, P const &pred) |
| Test whether any element in an array satisfies a unary predicate. | |
| template<typename P , int dim> | |
| bool | amrex::Reduce::AnyOf (BoxND< dim > const &box, P const &pred) |
| Test whether the predicate is true for any index in a Box. | |
| template<typename T , typename N , typename FIN , typename FOUT , typename TYPE , typename M = std::enable_if_t<std::is_integral_v<N> && (std::is_same_v<std::decay_t<TYPE>,Type::Inclusive> || std::is_same_v<std::decay_t<TYPE>,Type::Exclusive>)>> | |
| T | amrex::Scan::PrefixSum (N n, FIN const &fin, FOUT const &fout, TYPE, RetSum a_ret_sum=retSum) |
| template<typename N , typename T , typename M = std::enable_if_t<std::is_integral_v<N>>> | |
| T | amrex::Scan::InclusiveSum (N n, T const *in, T *out, RetSum a_ret_sum=retSum) |
| Inclusive sum. | |
| template<typename N , typename T , typename M = std::enable_if_t<std::is_integral_v<N>>> | |
| T | amrex::Scan::ExclusiveSum (N n, T const *in, T *out, RetSum a_ret_sum=retSum) |
| Exclusive sum. | |
| template<class InIter , class OutIter > | |
| OutIter | amrex::Gpu::inclusive_scan (InIter begin, InIter end, OutIter result) |
| template<class InIter , class OutIter > | |
| OutIter | amrex::Gpu::exclusive_scan (InIter begin, InIter end, OutIter result) |
Performance-portable reduction and scan operations in AMReX.
These interfaces provide reductions (e.g., sum, min, max) and scan operations (prefix sums) over index spaces and ranges on both CPUs and GPUs. They do not use MPI and operate on data local to the process.
Key functions and classes include:
| bool amrex::Reduce::AnyOf | ( | BoxND< dim > const & | box, |
| P const & | pred | ||
| ) |
Test whether the predicate is true for any index in a Box.
Note that this is local test that does not involve MPI.
| P | callable type |
| dim | dimension of the box |
| box | index region for the loop |
| pred | callable object with the signature of bool(int,int,int). |
| bool amrex::Reduce::AnyOf | ( | N | n, |
| T const * | v, | ||
| P const & | pred | ||
| ) |
Test whether any element in an array satisfies a unary predicate.
Note that this is local test that does not involve MPI.
| T | value type of the elements. |
| N | integral index type. |
| P | predicate type |
| n | number of elements to test. |
| v | pointer to the input data array. |
| pred | unary predicate object. |
| OutIter amrex::Gpu::exclusive_scan | ( | InIter | begin, |
| InIter | end, | ||
| OutIter | result | ||
| ) |
| T amrex::Scan::ExclusiveSum | ( | N | n, |
| T const * | in, | ||
| T * | out, | ||
| RetSum | a_ret_sum = retSum |
||
| ) |
Exclusive sum.
| n | number of elements |
| in | input |
| out | output |
| a_ret_sum | control the return value |
a_ret_sum is true | OutIter amrex::Gpu::inclusive_scan | ( | InIter | begin, |
| InIter | end, | ||
| OutIter | result | ||
| ) |
| T amrex::Scan::InclusiveSum | ( | N | n, |
| T const * | in, | ||
| T * | out, | ||
| RetSum | a_ret_sum = retSum |
||
| ) |
Inclusive sum.
| n | number of elements |
| in | input |
| out | output |
| a_ret_sum | control the return value |
a_ret_sum is true | T amrex::Reduce::Max | ( | N | n, |
| F const & | f, | ||
| T | init_val = std::numeric_limits<T>::lowest() |
||
| ) |
Compute the maximum of values generated by a callable.
Note that this is local reduction that does not involve MPI.
| T | value type of the reduction. |
| N | integral index type. |
| F | callable type with the signature T(N). |
| n | number of elements to reduce. |
| f | callable object invoked as f(i). For GPU builds, this needs to be a device lambda function or functor. |
| init_val | initial value (default: std::numeric_limits<T>::lowest()). |
| T amrex::Reduce::Max | ( | N | n, |
| T const * | v, | ||
| T | init_val = std::numeric_limits<T>::lowest() |
||
| ) |
Compute the maximum of an array of values.
Note that this is local reduction that does not involve MPI.
| T | value type of the reduction. |
| N | integral index type. |
| n | number of elements to reduce. |
| v | pointer to the input data array. |
| init_val | initial value (default: std::numeric_limits<T>::lowest()). |
| T amrex::Reduce::Min | ( | N | n, |
| F const & | f, | ||
| T | init_val = std::numeric_limits<T>::max() |
||
| ) |
Compute the minimum of values generated by a callable.
Note that this is local reduction that does not involve MPI.
| T | value type of the reduction. |
| N | integral index type. |
| F | callable type with the signature T(N). |
| n | number of elements to reduce. |
| f | callable object invoked as f(i). For GPU builds, this needs to be a device lambda function or functor. |
| init_val | initial value (default: std::numeric_limits<T>::max()). |
| T amrex::Reduce::Min | ( | N | n, |
| T const * | v, | ||
| T | init_val = std::numeric_limits<T>::max() |
||
| ) |
Compute the minimum of an array of values.
Note that this is local reduction that does not involve MPI.
| T | value type of the reduction. |
| N | integral index type. |
| n | number of elements to reduce. |
| v | pointer to the input data array. |
| init_val | initial value (default: std::numeric_limits<T>::max()). |
| std::pair< T, T > amrex::Reduce::MinMax | ( | N | n, |
| F const & | f | ||
| ) |
Compute the minimum and maximum of values generated by a callable.
Note that this is local reduction that does not involve MPI.
| T | value type of the reduction. |
| N | integral index type. |
| F | callable type with the signature T(N). |
| n | number of elements to reduce. |
| f | callable object invoked as f(i). For GPU builds, this needs to be a device lambda function or functor. |
| std::pair< T, T > amrex::Reduce::MinMax | ( | N | n, |
| T const * | v | ||
| ) |
Compute the minimum and maximum of an array of values.
Note that this is local reduction that does not involve MPI.
| T | value type of the reduction. |
| N | integral index type. |
| n | number of elements to reduce. |
| v | pointer to the input data array. |
| T amrex::Scan::PrefixSum | ( | N | n, |
| FIN const & | fin, | ||
| FOUT const & | fout, | ||
| TYPE | , | ||
| RetSum | a_ret_sum = retSum |
||
| ) |
| T amrex::Reduce::Sum | ( | N | n, |
| F const & | f, | ||
| T | init_val = 0 |
||
| ) |
Compute the sum of values generated by a callable.
Note that this is local reduction that does not involve MPI.
| T | value type of the reduction. |
| N | integral index type. |
| F | callable type with the signature T(N). |
| n | number of elements to reduce. |
| f | callable object invoked as f(i). For GPU builds, this needs to be a device lambda function or functor. |
| init_val | initial value (default: 0). |
| T amrex::Reduce::Sum | ( | N | n, |
| T const * | v, | ||
| T | init_val = 0 |
||
| ) |
Compute the sum of an array of values.
Note that this is local reduction that does not involve MPI.
| T | value type of the reduction. |
| N | integral index type. |
| n | number of elements to reduce. |
| v | pointer to the input data array. |
| init_val | initial value (default: 0). |