Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::Reducer< Ops, Ts > Class Template Reference

Class for local reductions (e.g., sum, min and max). More...

#include <AMReX_Reduce.H>

Inheritance diagram for amrex::Reducer< Ops, Ts >:

Public Types

using Result_t = typename Base::Result_t
 Reduction result type, GpuTuple<U...>, where U... are the types in Ts.
 

Public Member Functions

 Reducer ()=default
 
 ~Reducer ()=default
 
template<typename F , int dim>
std::enable_if_t< IsCallable< F, int, int, int >::value||IsCallable< F, IntVectND< dim > >::value > eval (BoxND< dim > const &box, F &&f)
 Reduction over a Box.
 
template<typename F , int dim>
std::enable_if_t< IsCallable< F, int, int, int, int >::value||IsCallable< F, IntVectND< dim >, int >::value > eval (BoxND< dim > const &box, int ncomp, F &&f)
 Reduction over a Box plus component index.
 
template<typename MF , typename F >
std::enable_if_t< IsFabArray< MF >::value &&IsCallable< F, int, int, int, int >::value > eval (MF const &mf, IntVect const &nghost, F &&f)
 Reduction over a MultiFab-like object.
 
template<typename MF , typename F >
std::enable_if_t< IsFabArray< MF >::value &&IsCallable< F, int, int, int, int, int >::value > eval (MF const &mf, IntVect const &nghost, int ncomp, F &&f)
 Reduction over a MultiFab-like object.
 
template<typename N , typename F >
std::enable_if_t< IsCallable< F, N >::value > eval (N n, F &&f)
 Reduction over a 1D index range.
 
Result_t getResult ()
 Get the final reduction result.
 

Static Public Attributes

static constexpr int size = GpuTupleSize<Result_t>::value
 Number of reduction operators.
 

Detailed Description

template<typename Ops, typename Ts>
class amrex::Reducer< Ops, Ts >

Class for local reductions (e.g., sum, min and max).

This class wraps ReduceOps and ReduceData to provide support for local reductions. A typical usage is:

using Result_t = typename decltype(reducer)::Result_t;
// In the code below, mf is a MulitFab and imf is an iMultiFab.
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (MFIter mfi(mf,TilingIfNotGpu()); mfi.isValid(); ++mfi) {
Box const& box = mfi.tilebox();
Array4<Real const> const& ar = mf.const_array(mfi);
Array4<int const> const& ai = imf.const_array(mfi);
reducer.eval(box, [=] AMREX_GPU_DEVICE (int i, int j, int k) -> Result_t
{
auto vr = ar(i,j,k);
auto vi = ai(i,j,k);
return {vi,vr,std::abs(vr),vr*vr};
});
}
auto result = reducer.getResult();
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
Class for local reductions (e.g., sum, min and max).
Definition AMReX_Reduce.H:1607
typename Base::Result_t Result_t
Reduction result type, GpuTuple<U...>, where U... are the types in Ts.
Definition AMReX_Reduce.H:1611
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Struct for holding types.
Definition AMReX_TypeList.H:13

The code above performs a sum reduction on an iMultiFab, and on a MultiFab, it computes the maximum, the maximum of the absolute value, and the sum of squares. The data type of result is GpuTuple<int,Real,Real,Real>.

If the type list is long, amrex::TypeMultiplier can be used to reduce verbosity. For example, the equvalient definition of Reducer above is

TypeAt< 0, decltype(detail::TApply< TParam >((TypeList<>{}+...+detail::SingleTypeMultiplier(std::declval< Types >()))))> TypeMultiplier
Return the first template argument with the later arguments applied to it. Types of the form T[N] are...
Definition AMReX_TypeList.H:221

If there is only a single reduction operator, the use of TypeList may be omitted. For example,

Multiple variants of the eval function are provided for looping over Box, Box plus component, 1D index ranges, and fused loops over MultiFab like objects.

User defined types are supported as long as they provide the operations required by the reduction operators. amrex::ValLocPair is such an example of a custom type that supports ReduceOpMin and ReduceOpMax.

Note
This is a local reduction and does not involve MPI communication. For a global reduction across ranks, reduce the local result with MPI.
A Reducer instance is intended to be used as an ephemeral object. Create a new instance for a new reduction rather than reusing one across unrelated loops.
Template Parameters
OpsReduction operator type(s). May be a single operator (e.g., ReduceOpSum) or a TypeList of operators.
TsValue type(s) corresponding to Ops. May be a single type or a TypeList of types.

Member Typedef Documentation

◆ Result_t

template<typename Ops , typename Ts >
using amrex::Reducer< Ops, Ts >::Result_t = typename Base::Result_t

Reduction result type, GpuTuple<U...>, where U... are the types in Ts.

Constructor & Destructor Documentation

◆ Reducer()

template<typename Ops , typename Ts >
amrex::Reducer< Ops, Ts >::Reducer ( )
default

◆ ~Reducer()

template<typename Ops , typename Ts >
amrex::Reducer< Ops, Ts >::~Reducer ( )
default

Member Function Documentation

◆ eval() [1/5]

template<typename Ops , typename Ts >
template<typename F , int dim>
std::enable_if_t< IsCallable< F, int, int, int >::value|| IsCallable< F, IntVectND< dim > >::value > amrex::Reducer< Ops, Ts >::eval ( BoxND< dim > const &  box,
F &&  f 
)
inline

Reduction over a Box.

This loops over all integer indices in box, i.e., a dim-dimenstional iteration space.

Template Parameters
Fcallable type.
dimdimension of the Box.
Parameters
boxBox
fcallable object returning Result_t for reduction. Its signature is either Result_t(int,int,int) or Result_t(IntVectND<dim>).

◆ eval() [2/5]

template<typename Ops , typename Ts >
template<typename F , int dim>
std::enable_if_t< IsCallable< F, int, int, int, int >::value|| IsCallable< F, IntVectND< dim >, int >::value > amrex::Reducer< Ops, Ts >::eval ( BoxND< dim > const &  box,
int  ncomp,
F &&  f 
)
inline

Reduction over a Box plus component index.

This loops over all integer indices in box and over components n in [0,ncomp), i.e., a (dim+1)-dimenstional iteration space.

Template Parameters
Fcallable type.
dimdimension of the Box.
Parameters
boxBox
ncompnumber of components
fcallable object returning Result_t for reduction. Its signature is either Result_t(int,int,int,int) or Result_t(IntVectND<dim>,int), with the last argument the component index in both versions.

◆ eval() [3/5]

template<typename Ops , typename Ts >
template<typename MF , typename F >
std::enable_if_t< IsFabArray< MF >::value && IsCallable< F, int, int, int, int >::value > amrex::Reducer< Ops, Ts >::eval ( MF const &  mf,
IntVect const &  nghost,
F &&  f 
)
inline

Reduction over a MultiFab-like object.

This loops over the valid region of a MultiFab-like object plus nghost layers of ghost cells.

The MultiFab-like object mf is used only to define the iteration space. The values participating in the reduction are provided by the user-supplied callable f.

Template Parameters
MFMultiFab-like type.
Fcallable type.
Parameters
mfMultiFab-like object for providing the iteration space.
nghostnumber of ghost-cell layers included in the iteration space.
fcallable object returning Result_t for reduction. Its signature is Result_t(int,int,int,int), where the arguments correspond to the local box index in mf and the $(i,j,k)$ cell indices.

◆ eval() [4/5]

template<typename Ops , typename Ts >
template<typename MF , typename F >
std::enable_if_t< IsFabArray< MF >::value && IsCallable< F, int, int, int, int, int >::value > amrex::Reducer< Ops, Ts >::eval ( MF const &  mf,
IntVect const &  nghost,
int  ncomp,
F &&  f 
)
inline

Reduction over a MultiFab-like object.

This loops over the valid region of a MultiFab-like object plus nghost layers of ghost cells and over components, i.e., a (AMREX_SPACEDIM+1)-dimensional iteration space.

The MultiFab-like object mf is used only to define the iteration space. The values participating in the reduction are provided by the user-supplied callable f.

Template Parameters
MFMultiFab-like type.
Fcallable type.
Parameters
mfMultiFab-like object for providing the iteration space.
nghostnumber of ghost-cell layers included in the iteration space.
ncompnumber of components.
fcallable object returning Result_t for reduction. Its signature is Result_t(int,int,int,int,int), where the arguments correspond to the local box index in mf, the $(i,j,k)$ cell indices, and the component index.

◆ eval() [5/5]

template<typename Ops , typename Ts >
template<typename N , typename F >
std::enable_if_t< IsCallable< F, N >::value > amrex::Reducer< Ops, Ts >::eval ( n,
F &&  f 
)
inline

Reduction over a 1D index range.

This loops over [0,n).

Template Parameters
Nintegral index type
Fcallable type.
Parameters
nexclusive upper bound of the 1d loop.
fcallable object returning Result_t for reduction. It must be invocable as Result_t(N).

◆ getResult()

template<typename Ops , typename Ts >
Result_t amrex::Reducer< Ops, Ts >::getResult ( )
inline

Get the final reduction result.

The return data type is Result_t (i.e., amrex::GpuTuple). Individual elements can be accessed using amrex::get<I>(). If all underlying value types are the same, the result can also be converted to an amrex::GpuArray via amrex::tupleToArray() and then accessed using operator[].

Returns
reduction result as Result_t.

Member Data Documentation

◆ size

template<typename Ops , typename Ts >
constexpr int amrex::Reducer< Ops, Ts >::size = GpuTupleSize<Result_t>::value
staticconstexpr

Number of reduction operators.


The documentation for this class was generated from the following file: