|
| | 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.
|
| |
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:
#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
Box const& box = mfi.tilebox();
{
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
-
| Ops | Reduction operator type(s). May be a single operator (e.g., ReduceOpSum) or a TypeList of operators. |
| Ts | Value type(s) corresponding to Ops. May be a single type or a TypeList of types. |