ScatterView
¶
Header File: <Kokkos_ScatterView.hpp>
Warning
ScatterView
is still in the namespace Kokkos::Experimental
A Kokkos ScatterView serves as an interface for a standard Kokkos::View
implementing a scatter-add pattern either via atomics or data replication.
Construction of a ScatterView can be expensive, so you should try to reuse the same one if possible, in which case, you should call reset()
between uses.
ScatterView can not be addressed directly: each thread inside a parallel region needs to make a call to access()
and access the underlying View through the return value of access()
.
Following the parallel region, a call to the free function Kokkos::Experimental::contribute()
should be made to perform the final reduction.
template <typename DataType [, typename Layout [, typename ExecSpace [, typename Operation [, typename Duplication [, typename Contribution]]]]]> class ScatterViewParameters¶
Template parameters other than DataType
are optional, but if one is specified, preceding ones must also be specified. That means for example that Operation
can be omitted but if it is specified, Layout
and ExecSpace
must also be specified.
DataType
, Layout
and ExecSpace
need to be the same types as the one from the Kokkos::View this ScatterView is interfacing.
Operation
: Can take the values:
Kokkos::Experimental::ScatterSum
: performs a Sum. It is the default value.
Kokkos::Experimental::ScatterProd
: performs a Multiplication.
Kokkos::Experimental::ScatterMin
: takes the min.
Kokkos::Experimental::ScatterMax
: takes the max.
Duplication
: Whether to duplicate the grid or not; defaults to Kokkos::Experimental::ScatterDuplicated
, other option is Kokkos::Experimental::ScatterNonDuplicated
.
Contribution
: Whether to contribute to use atomics; defaults to Kokkos::Experimental::ScatterAtomics
, other option is Kokkoss::Experimental::ScatterNonAtomic
.
Creating a ScatterView with non default Operation
, Duplication
or Contribution
using this interface can become complicated, because you need to specify the exact type for DataType
, Layout
and ExecSpace
. This is why it is advised that you instead use the function Kokkos::Experimental::create_scatter_view()
.
Public Member Variables
Type of View passed to ScatterView constructor.
Value type of the original_view_type.
Reference type of the original_view_type.
DuplicatedDataType, a newly created DataType that has a new runtime dimension which becomes the largest-stride dimension, from the given View DataType.
Value type of data_type_info.
A View type created from the internal_data_type.
Constructors
The default constructor. Default constructs members.
Constructor from a Kokkos::View
. internal_view
member is copy constructed from this input view.
Constructor from variadic pack of dimension arguments. Constructs internal_view
member.
Constructor from variadic pack of dimension arguments. Constructs internal_view
member. This constructor allows passing an object created by Kokkos::view_alloc
as first argument, e.g., for specifying an execution space via Kokkos::view_alloc(exec_space, "label")
.
Public Methods
true if the internal_view
points to a valid memory location. This function works for both managed and unmanaged views. With the unmanaged view, there is no guarantee that referenced address is valid, only that it is a non-null pointer.
use within a kernel to return a ScatterAccess
member; this member accumulates a given thread’s contribution to the reduction.
a subview of a ScatterView
contribute ScatterView
array’s results into the input View dest
performs reset on destination array
tbd
resize a view with copying old data to new data at the corresponding indices
resize a view with discarding old data
Private Members
typedef original_view_type internal_view_type;
internal_view_type internal_view;
Free Functions
create a new ScatterView interfacing the View view
. Default value for Operation
is Kokkos::Experimental::ScatterSum
, Duplication
and Contribution
are chosen to make the ScatterView as efficient as possible when running on its ExecSpace
.
convenience function to perform final reduction of ScatterView results into a resultant View; may be called following parallel_reduce()
.
#include <Kokkos_Core.hpp> #include <Kokkos_ScatterView.hpp> KOKKOS_INLINE_FUNCTION int foo(int i) { return i; } KOKKOS_INLINE_FUNCTION double bar(int i) { return i*i; } int main (int argc, char* argv[]) { Kokkos::ScopeGuard guard(argc, argv); Kokkos::View<double*> results("results", 1); auto scatter = Kokkos::Experimental::create_scatter_view(results); Kokkos::parallel_for(1, KOKKOS_LAMBDA(int input_i) { auto access = scatter.access(); auto result_i = foo(input_i); auto contribution = bar(input_i); access(result_i) += contribution; }); Kokkos::Experimental::contribute(results, scatter); }
RetroSearch is an open source project built by @garambo | Open a GitHub Issue
Search and Browse the WWW like it's 1997 | Search results from DuckDuckGo
HTML:
3.2
| Encoding:
UTF-8
| Version:
0.7.4