Kokkos Core Kernels Package  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Kokkos_ScatterView.hpp
Go to the documentation of this file.
1 /*
2 //@HEADER
3 // ************************************************************************
4 //
5 // Kokkos v. 2.0
6 // Copyright (2014) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Christian R. Trott (crtrott@sandia.gov)
39 //
40 // ************************************************************************
41 //@HEADER
42 */
43 
44 
50 
51 #ifndef KOKKOS_SCATTER_VIEW_HPP
52 #define KOKKOS_SCATTER_VIEW_HPP
53 
54 #include <Kokkos_Core.hpp>
55 #include <utility>
56 
57 namespace Kokkos {
58 namespace Experimental {
59 
60 /*
61  * Reduction Type list
62  * - These corresponds to subset of the reducers in parallel_reduce
63  * - See Implementations of ScatterValue for details.
64  */
65 enum : int {
66  ScatterSum,
67  ScatterProd,
68  ScatterMax,
69  ScatterMin,
70 };
71 
72 enum : int {
73  ScatterNonDuplicated = 0,
74  ScatterDuplicated = 1
75 };
76 
77 enum : int {
78  ScatterNonAtomic = 0,
79  ScatterAtomic = 1
80 };
81 
82 }} // Kokkos::Experimental
83 
84 namespace Kokkos {
85 namespace Impl {
86 namespace Experimental {
87 
88 template <typename ExecSpace>
89 struct DefaultDuplication;
90 
91 template <typename ExecSpace, int duplication>
92 struct DefaultContribution;
93 
94 #ifdef KOKKOS_ENABLE_SERIAL
95 template <>
96 struct DefaultDuplication<Kokkos::Serial> {
97  enum : int { value = Kokkos::Experimental::ScatterNonDuplicated };
98 };
99 template <>
100 struct DefaultContribution<Kokkos::Serial, Kokkos::Experimental::ScatterNonDuplicated> {
101  enum : int { value = Kokkos::Experimental::ScatterNonAtomic };
102 };
103 template <>
104 struct DefaultContribution<Kokkos::Serial, Kokkos::Experimental::ScatterDuplicated> {
105  enum : int { value = Kokkos::Experimental::ScatterNonAtomic };
106 };
107 #endif
108 
109 #ifdef KOKKOS_ENABLE_OPENMP
110 template <>
111 struct DefaultDuplication<Kokkos::OpenMP> {
112  enum : int { value = Kokkos::Experimental::ScatterDuplicated };
113 };
114 template <>
115 struct DefaultContribution<Kokkos::OpenMP, Kokkos::Experimental::ScatterNonDuplicated> {
116  enum : int { value = Kokkos::Experimental::ScatterAtomic };
117 };
118 template <>
119 struct DefaultContribution<Kokkos::OpenMP, Kokkos::Experimental::ScatterDuplicated> {
120  enum : int { value = Kokkos::Experimental::ScatterNonAtomic };
121 };
122 #endif
123 
124 #ifdef KOKKOS_ENABLE_HPX
125 template <>
126 struct DefaultDuplication<Kokkos::Experimental::HPX> {
127  enum : int { value = Kokkos::Experimental::ScatterDuplicated };
128 };
129 template <>
130 struct DefaultContribution<Kokkos::Experimental::HPX, Kokkos::Experimental::ScatterNonDuplicated> {
131  enum : int { value = Kokkos::Experimental::ScatterAtomic };
132 };
133 template <>
134 struct DefaultContribution<Kokkos::Experimental::HPX, Kokkos::Experimental::ScatterDuplicated> {
135  enum : int { value = Kokkos::Experimental::ScatterNonAtomic };
136 };
137 #endif
138 
139 #ifdef KOKKOS_ENABLE_THREADS
140 template <>
141 struct DefaultDuplication<Kokkos::Threads> {
142  enum : int { value = Kokkos::Experimental::ScatterDuplicated };
143 };
144 template <>
145 struct DefaultContribution<Kokkos::Threads, Kokkos::Experimental::ScatterNonDuplicated> {
146  enum : int { value = Kokkos::Experimental::ScatterAtomic };
147 };
148 template <>
149 struct DefaultContribution<Kokkos::Threads, Kokkos::Experimental::ScatterDuplicated> {
150  enum : int { value = Kokkos::Experimental::ScatterNonAtomic };
151 };
152 #endif
153 
154 #ifdef KOKKOS_ENABLE_CUDA
155 template <>
156 struct DefaultDuplication<Kokkos::Cuda> {
157  enum : int { value = Kokkos::Experimental::ScatterNonDuplicated };
158 };
159 template <>
160 struct DefaultContribution<Kokkos::Cuda, Kokkos::Experimental::ScatterNonDuplicated> {
161  enum : int { value = Kokkos::Experimental::ScatterAtomic };
162 };
163 template <>
164 struct DefaultContribution<Kokkos::Cuda, Kokkos::Experimental::ScatterDuplicated> {
165  enum : int { value = Kokkos::Experimental::ScatterAtomic };
166 };
167 #endif
168 
169 /* ScatterValue <Op=ScatterSum, contribution=ScatterNonAtomic> is the object returned by the access operator() of ScatterAccess,
170  This class inherits from the Sum<> reducer and it wraps join(dest, src) with convenient operator+=, etc.
171  Note the addition of update(ValueType const& rhs) and reset() so that all reducers can have common functions
172  See ReduceDuplicates and ResetDuplicates ) */
173 template <typename ValueType, int Op, int contribution>
174 struct ScatterValue;
175 
176 template <typename ValueType>
177 struct ScatterValue<ValueType, Kokkos::Experimental::ScatterSum, Kokkos::Experimental::ScatterNonAtomic> :
178  Sum<ValueType,Kokkos::DefaultExecutionSpace> {
179  public:
180  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ValueType& value_in) :
181  Sum<ValueType,Kokkos::DefaultExecutionSpace>(value_in)
182  {}
183  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ScatterValue&& other) :
184  Sum<ValueType,Kokkos::DefaultExecutionSpace>(other.reference())
185  {}
186  KOKKOS_FORCEINLINE_FUNCTION void operator+=(ValueType const& rhs) {
187  this->join( this->reference(), rhs );
188  }
189  KOKKOS_FORCEINLINE_FUNCTION void operator-=(ValueType const& rhs) {
190  this->join( this->reference(), -rhs );
191  }
192  KOKKOS_FORCEINLINE_FUNCTION void update(ValueType const& rhs) {
193  this->join( this->reference(), rhs );
194  }
195  KOKKOS_FORCEINLINE_FUNCTION void reset() {
196  this->init( this->reference() );
197  }
198 };
199 
200 /* ScatterValue <Op=ScatterSum, contribution=ScatterAtomic> is the object returned by the access operator()
201  * of ScatterAccess, similar to that returned by an Atomic View, it wraps Kokkos::atomic_add with convenient
202  operator+=, etc. This version also has the update(rhs) and reset() functions. */
203 template <typename ValueType>
204 struct ScatterValue<ValueType, Kokkos::Experimental::ScatterSum, Kokkos::Experimental::ScatterAtomic> :
205  Sum<ValueType,Kokkos::DefaultExecutionSpace> {
206  public:
207  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ValueType& value_in) :
208  Sum<ValueType,Kokkos::DefaultExecutionSpace>(value_in)
209  {}
210 
211  KOKKOS_FORCEINLINE_FUNCTION void operator+=(ValueType const& rhs) {
212  this->join(this->reference(), rhs);
213  }
214  KOKKOS_FORCEINLINE_FUNCTION void operator-=(ValueType const& rhs) {
215  this->join(this->reference(), -rhs);
216  }
217 
218  KOKKOS_INLINE_FUNCTION
219  void join(ValueType& dest, const ValueType& src) const {
220  Kokkos::atomic_add(&dest, src);
221  }
222 
223  KOKKOS_INLINE_FUNCTION
224  void join(volatile ValueType& dest, const volatile ValueType& src) const {
225  Kokkos::atomic_add(&dest, src);
226  }
227 
228  KOKKOS_FORCEINLINE_FUNCTION void update(ValueType const& rhs) {
229  this->join( this->reference(), rhs );
230  }
231 
232  KOKKOS_FORCEINLINE_FUNCTION void reset() {
233  this->init( this->reference() );
234  }
235 };
236 
237 /* ScatterValue <Op=ScatterProd, contribution=ScatterNonAtomic> is the object returned by the access operator() of ScatterAccess,
238  This class inherits from the Prod<> reducer and it wraps join(dest, src) with convenient operator*=, etc.
239  Note the addition of update(ValueType const& rhs) and reset() so that all reducers can have common functions
240  See ReduceDuplicates and ResetDuplicates ) */
241 template <typename ValueType>
242 struct ScatterValue<ValueType, Kokkos::Experimental::ScatterProd, Kokkos::Experimental::ScatterNonAtomic> :
243  Prod<ValueType,Kokkos::DefaultExecutionSpace> {
244  public:
245  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ValueType& value_in) :
246  Prod<ValueType,Kokkos::DefaultExecutionSpace>(value_in)
247  {}
248  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ScatterValue&& other) :
249  Prod<ValueType,Kokkos::DefaultExecutionSpace>(other.reference())
250  {}
251  KOKKOS_FORCEINLINE_FUNCTION void operator*=(ValueType const& rhs) {
252  this->join( this->reference(), rhs );
253  }
254  KOKKOS_FORCEINLINE_FUNCTION void operator/=(ValueType const& rhs) {
255  this->join( this->reference(), static_cast<ValueType>(1)/rhs );
256  }
257  KOKKOS_FORCEINLINE_FUNCTION void update(ValueType const& rhs) {
258  this->join( this->reference(), rhs );
259  }
260  KOKKOS_FORCEINLINE_FUNCTION void reset() {
261  this->init( this->reference() );
262  }
263 };
264 
265 /* ScatterValue <Op=ScatterProd, contribution=ScatterAtomic> is the object returned by the access operator()
266  * of ScatterAccess, similar to that returned by an Atomic View, it wraps and atomic_prod with convenient
267  operator*=, etc. atomic_prod uses the atomic_compare_exchange. This version also has the update(rhs) and reset() functions. */
268 template <typename ValueType>
269 struct ScatterValue<ValueType, Kokkos::Experimental::ScatterProd, Kokkos::Experimental::ScatterAtomic> :
270  Prod<ValueType,Kokkos::DefaultExecutionSpace> {
271  public:
272  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ValueType& value_in) :
273  Prod<ValueType,Kokkos::DefaultExecutionSpace>(value_in)
274  {}
275 
276  KOKKOS_FORCEINLINE_FUNCTION void operator*=(ValueType const& rhs) {
277  this->join(this->reference(), rhs);
278  }
279  KOKKOS_FORCEINLINE_FUNCTION void operator/=(ValueType const& rhs) {
280  this->join(this->reference(), static_cast<ValueType>(1)/rhs);
281  }
282 
283  KOKKOS_FORCEINLINE_FUNCTION
284  void atomic_prod(ValueType & dest, const ValueType& src) const {
285 
286  bool success = false;
287  while(!success) {
288  ValueType dest_old = dest;
289  ValueType dest_new = dest_old * src;
290  dest_new = Kokkos::atomic_compare_exchange<ValueType>(&dest,dest_old,dest_new);
291  success = ( (dest_new - dest_old)/dest_old <= 1e-15 );
292  }
293  }
294 
295  KOKKOS_INLINE_FUNCTION
296  void join(ValueType& dest, const ValueType& src) const {
297  atomic_prod(dest, src);
298  }
299 
300  KOKKOS_INLINE_FUNCTION
301  void join(volatile ValueType& dest, const volatile ValueType& src) const {
302  atomic_prod(dest, src);
303  }
304 
305  KOKKOS_FORCEINLINE_FUNCTION void update(ValueType const& rhs) {
306  this->join( this->reference(), rhs );
307  }
308  KOKKOS_FORCEINLINE_FUNCTION void reset() {
309  this->init( this->reference() );
310  }
311 
312 };
313 
314 /* ScatterValue <Op=ScatterMin, contribution=ScatterNonAtomic> is the object returned by the access operator() of ScatterAccess,
315  This class inherits from the Min<> reducer and it wraps join(dest, src) with convenient update(rhs).
316  Note the addition of update(ValueType const& rhs) and reset() are so that all reducers can have a common update function
317  See ReduceDuplicates and ResetDuplicates ) */
318 template <typename ValueType>
319 struct ScatterValue<ValueType, Kokkos::Experimental::ScatterMin, Kokkos::Experimental::ScatterNonAtomic> :
320  Min<ValueType,Kokkos::DefaultExecutionSpace> {
321  public:
322  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ValueType& value_in) :
323  Min<ValueType,Kokkos::DefaultExecutionSpace>(value_in)
324  {}
325  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ScatterValue&& other) :
326  Min<ValueType,Kokkos::DefaultExecutionSpace>(other.reference())
327  {}
328  KOKKOS_FORCEINLINE_FUNCTION void update(ValueType const& rhs) {
329  this->join( this->reference(), rhs );
330  }
331  KOKKOS_FORCEINLINE_FUNCTION void reset() {
332  this->init( this->reference() );
333  }
334 };
335 
336 /* ScatterValue <Op=ScatterMin, contribution=ScatterAtomic> is the object returned by the access operator()
337  * of ScatterAccess, similar to that returned by an Atomic View, it wraps and atomic_min with the update(rhs)
338  function. atomic_min uses the atomic_compare_exchange. This version also has the reset() function */
339 template <typename ValueType>
340 struct ScatterValue<ValueType, Kokkos::Experimental::ScatterMin, Kokkos::Experimental::ScatterAtomic> :
341  Min<ValueType,Kokkos::DefaultExecutionSpace> {
342  public:
343  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ValueType& value_in) :
344  Min<ValueType,Kokkos::DefaultExecutionSpace>(value_in)
345  {}
346 
347  KOKKOS_FORCEINLINE_FUNCTION
348  void atomic_min(ValueType & dest, const ValueType& src) const {
349 
350  bool success = false;
351  while(!success) {
352  ValueType dest_old = dest;
353  ValueType dest_new = ( dest_old > src ) ? src : dest_old;
354  dest_new = Kokkos::atomic_compare_exchange<ValueType>(&dest,dest_old,dest_new);
355  success = ( (dest_new - dest_old)/dest_old <= 1e-15 );
356  }
357  }
358 
359  KOKKOS_INLINE_FUNCTION
360  void join(ValueType& dest, const ValueType& src) const {
361  atomic_min(dest, src);
362  }
363 
364  KOKKOS_INLINE_FUNCTION
365  void join(volatile ValueType& dest, const volatile ValueType& src) const {
366  atomic_min(dest, src);
367  }
368 
369  KOKKOS_FORCEINLINE_FUNCTION void update(ValueType const& rhs) {
370  this->join( this->reference(), rhs );
371  }
372  KOKKOS_FORCEINLINE_FUNCTION void reset() {
373  this->init( this->reference() );
374  }
375 
376 };
377 
378 /* ScatterValue <Op=ScatterMax, contribution=ScatterNonAtomic> is the object returned by the access operataor() of ScatterAccess,
379  This class inherits from the Max<> reducer and it wraps join(dest, src) with convenient update(rhs).
380  Note the addition of update(ValueType const& rhs) and reset() are so that all reducers can have a common update function
381  See ReduceDuplicates and ResetDuplicates ) */
382 template <typename ValueType>
383 struct ScatterValue<ValueType, Kokkos::Experimental::ScatterMax, Kokkos::Experimental::ScatterNonAtomic> :
384  Max<ValueType,Kokkos::DefaultExecutionSpace> {
385  public:
386  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ValueType& value_in) :
387  Max<ValueType,Kokkos::DefaultExecutionSpace>(value_in)
388  {}
389  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ScatterValue&& other) :
390  Max<ValueType,Kokkos::DefaultExecutionSpace>(other.reference())
391  {}
392  KOKKOS_FORCEINLINE_FUNCTION void update(ValueType const& rhs) {
393  this->join( this->reference(), rhs );
394  }
395  KOKKOS_FORCEINLINE_FUNCTION void reset() {
396  this->init( this->reference() );
397  }
398 };
399 
400 /* ScatterValue <Op=ScatterMax, contribution=ScatterAtomic> is the object returned by the access operator()
401  * of ScatterAccess, similar to that returned by an Atomic View, it wraps and atomic_max with the update(rhs)
402  function. atomic_max uses the atomic_compare_exchange. This version also has the reset() function */
403 template <typename ValueType>
404 struct ScatterValue<ValueType, Kokkos::Experimental::ScatterMax, Kokkos::Experimental::ScatterAtomic> :
405  Max<ValueType,Kokkos::DefaultExecutionSpace> {
406  public:
407  KOKKOS_FORCEINLINE_FUNCTION ScatterValue(ValueType& value_in) :
408  Max<ValueType,Kokkos::DefaultExecutionSpace>(value_in)
409  {}
410 
411  KOKKOS_FORCEINLINE_FUNCTION
412  void atomic_max(ValueType & dest, const ValueType& src) const {
413 
414  bool success = false;
415  while(!success) {
416  ValueType dest_old = dest;
417  ValueType dest_new = ( dest_old < src ) ? src : dest_old;
418  dest_new = Kokkos::atomic_compare_exchange<ValueType>(&dest,dest_old,dest_new);
419  success = ( (dest_new - dest_old)/dest_old <= 1e-15 );
420  }
421  }
422 
423  KOKKOS_INLINE_FUNCTION
424  void join(ValueType& dest, const ValueType& src) const {
425  atomic_max(dest, src);
426  }
427 
428  KOKKOS_INLINE_FUNCTION
429  void join(volatile ValueType& dest, const volatile ValueType& src) const {
430  atomic_max(dest, src);
431  }
432 
433  KOKKOS_FORCEINLINE_FUNCTION void update(ValueType const& rhs) {
434  this->join( this->reference(), rhs );
435  }
436  KOKKOS_FORCEINLINE_FUNCTION void reset() {
437  this->init( this->reference() );
438  }
439 
440 };
441 
442 /* DuplicatedDataType, given a View DataType, will create a new DataType
443  that has a new runtime dimension which becomes the largest-stride dimension.
444  In the case of LayoutLeft, due to the limitation induced by the design of DataType
445  itself, it must convert any existing compile-time dimensions into runtime dimensions. */
446 template <typename T, typename Layout>
447 struct DuplicatedDataType;
448 
449 template <typename T>
450 struct DuplicatedDataType<T, Kokkos::LayoutRight> {
451  typedef T* value_type; // For LayoutRight, add a star all the way on the left
452 };
453 
454 template <typename T, size_t N>
455 struct DuplicatedDataType<T[N], Kokkos::LayoutRight> {
456  typedef typename DuplicatedDataType<T, Kokkos::LayoutRight>::value_type value_type[N];
457 };
458 
459 template <typename T>
460 struct DuplicatedDataType<T[], Kokkos::LayoutRight> {
461  typedef typename DuplicatedDataType<T, Kokkos::LayoutRight>::value_type value_type[];
462 };
463 
464 template <typename T>
465 struct DuplicatedDataType<T*, Kokkos::LayoutRight> {
466  typedef typename DuplicatedDataType<T, Kokkos::LayoutRight>::value_type* value_type;
467 };
468 
469 template <typename T>
470 struct DuplicatedDataType<T, Kokkos::LayoutLeft> {
471  typedef T* value_type;
472 };
473 
474 template <typename T, size_t N>
475 struct DuplicatedDataType<T[N], Kokkos::LayoutLeft> {
476  typedef typename DuplicatedDataType<T, Kokkos::LayoutLeft>::value_type* value_type;
477 };
478 
479 template <typename T>
480 struct DuplicatedDataType<T[], Kokkos::LayoutLeft> {
481  typedef typename DuplicatedDataType<T, Kokkos::LayoutLeft>::value_type* value_type;
482 };
483 
484 template <typename T>
485 struct DuplicatedDataType<T*, Kokkos::LayoutLeft> {
486  typedef typename DuplicatedDataType<T, Kokkos::LayoutLeft>::value_type* value_type;
487 };
488 
489 /* Insert integer argument pack into array */
490 
491 template<class T>
492 void args_to_array(size_t* array, int pos, T dim0) {
493  array[pos] = dim0;
494 }
495 template<class T, class ... Dims>
496 void args_to_array(size_t* array, int pos, T dim0, Dims ... dims) {
497  array[pos] = dim0;
498  args_to_array(array,pos+1,dims...);
499 }
500 
501 /* Slice is just responsible for stuffing the correct number of Kokkos::ALL
502  arguments on the correct side of the index in a call to subview() to get a
503  subview where the index specified is the largest-stride one. */
504 template <typename Layout, int rank, typename V, typename ... Args>
505 struct Slice {
506  typedef Slice<Layout, rank - 1, V, Kokkos::Impl::ALL_t, Args...> next;
507  typedef typename next::value_type value_type;
508 
509  static
510  value_type get(V const& src, const size_t i, Args ... args) {
511  return next::get(src, i, Kokkos::ALL, args...);
512  }
513 };
514 
515 template <typename V, typename ... Args>
516 struct Slice<Kokkos::LayoutRight, 1, V, Args...> {
517  typedef typename Kokkos::Impl::ViewMapping
518  < void
519  , V
520  , const size_t
521  , Args ...
522  >::type value_type;
523  static
524  value_type get(V const& src, const size_t i, Args ... args) {
525  return Kokkos::subview(src, i, args...);
526  }
527 };
528 
529 template <typename V, typename ... Args>
530 struct Slice<Kokkos::LayoutLeft, 1, V, Args...> {
531  typedef typename Kokkos::Impl::ViewMapping
532  < void
533  , V
534  , Args ...
535  , const size_t
536  >::type value_type;
537  static
538  value_type get(V const& src, const size_t i, Args ... args) {
539  return Kokkos::subview(src, args..., i);
540  }
541 };
542 
543 template <typename ExecSpace, typename ValueType, int Op>
544 struct ReduceDuplicates;
545 
546 template <typename ExecSpace, typename ValueType, int Op>
547 struct ReduceDuplicatesBase {
548  typedef ReduceDuplicates<ExecSpace, ValueType, Op> Derived;
549  ValueType const* src;
550  ValueType* dst;
551  size_t stride;
552  size_t start;
553  size_t n;
554  ReduceDuplicatesBase(ValueType const* src_in, ValueType* dest_in, size_t stride_in, size_t start_in, size_t n_in, std::string const& name)
555  : src(src_in)
556  , dst(dest_in)
557  , stride(stride_in)
558  , start(start_in)
559  , n(n_in)
560  {
561 #if defined(KOKKOS_ENABLE_PROFILING)
562  uint64_t kpID = 0;
563  if(Kokkos::Profiling::profileLibraryLoaded()) {
564  Kokkos::Profiling::beginParallelFor(std::string("reduce_") + name, 0, &kpID);
565  }
566 #endif
567  typedef RangePolicy<ExecSpace, size_t> policy_type;
569  const closure_type closure(*(static_cast<Derived*>(this)), policy_type(0, stride));
570  closure.execute();
571 #if defined(KOKKOS_ENABLE_PROFILING)
572  if(Kokkos::Profiling::profileLibraryLoaded()) {
573  Kokkos::Profiling::endParallelFor(kpID);
574  }
575 #endif
576  }
577 };
578 
579 /* ReduceDuplicates -- Perform reduction on destination array using strided source
580  * Use ScatterValue<> specific to operation to wrap destination array so that
581  * the reduction operation can be accessed via the update(rhs) function */
582 template <typename ExecSpace, typename ValueType, int Op>
583 struct ReduceDuplicates :
584  public ReduceDuplicatesBase<ExecSpace, ValueType, Op>
585 {
586  typedef ReduceDuplicatesBase<ExecSpace, ValueType, Op> Base;
587  ReduceDuplicates(ValueType const* src_in, ValueType* dst_in, size_t stride_in, size_t start_in, size_t n_in, std::string const& name):
588  Base(src_in, dst_in, stride_in, start_in, n_in, name)
589  {}
590  KOKKOS_FORCEINLINE_FUNCTION void operator()(size_t i) const {
591  for (size_t j = Base::start; j < Base::n; ++j) {
592  ScatterValue<ValueType, Op, Kokkos::Experimental::ScatterNonAtomic> sv(Base::dst[i]);
593  sv.update( Base::src[i + Base::stride * j] );
594  }
595  }
596 };
597 
598 
599 template <typename ExecSpace, typename ValueType, int Op>
600 struct ResetDuplicates;
601 
602 template <typename ExecSpace, typename ValueType, int Op>
603 struct ResetDuplicatesBase {
604  typedef ResetDuplicates<ExecSpace, ValueType, Op> Derived;
605  ValueType* data;
606  ResetDuplicatesBase(ValueType* data_in, size_t size_in, std::string const& name)
607  : data(data_in)
608  {
609 #if defined(KOKKOS_ENABLE_PROFILING)
610  uint64_t kpID = 0;
611  if(Kokkos::Profiling::profileLibraryLoaded()) {
612  Kokkos::Profiling::beginParallelFor(std::string("reduce_") + name, 0, &kpID);
613  }
614 #endif
615  typedef RangePolicy<ExecSpace, size_t> policy_type;
617  const closure_type closure(*(static_cast<Derived*>(this)), policy_type(0, size_in));
618  closure.execute();
619 #if defined(KOKKOS_ENABLE_PROFILING)
620  if(Kokkos::Profiling::profileLibraryLoaded()) {
621  Kokkos::Profiling::endParallelFor(kpID);
622  }
623 #endif
624  }
625 };
626 
627 /* ResetDuplicates -- Perform reset on destination array
628  * Use ScatterValue<> specific to operation to wrap destination array so that
629  * the reset operation can be accessed via the reset() function */
630 template <typename ExecSpace, typename ValueType, int Op>
631 struct ResetDuplicates :
632  public ResetDuplicatesBase<ExecSpace, ValueType, Op>
633 {
634  typedef ResetDuplicatesBase<ExecSpace, ValueType, Op> Base;
635  ResetDuplicates(ValueType* data_in, size_t size_in, std::string const& name):
636  Base(data_in, size_in, name)
637  {}
638  KOKKOS_FORCEINLINE_FUNCTION void operator()(size_t i) const {
639  ScatterValue<ValueType, Op, Kokkos::Experimental::ScatterNonAtomic> sv(Base::data[i]);
640  sv.reset();
641  }
642 };
643 
644 
645 }}} // Kokkos::Impl::Experimental
646 
647 namespace Kokkos {
648 namespace Experimental {
649 
650 template <typename DataType
651  ,typename Layout = Kokkos::DefaultExecutionSpace::array_layout
652  ,typename ExecSpace = Kokkos::DefaultExecutionSpace
653  ,int Op = ScatterSum
654  ,int duplication = Kokkos::Impl::Experimental::DefaultDuplication<ExecSpace>::value
655  ,int contribution = Kokkos::Impl::Experimental::DefaultContribution<ExecSpace, duplication>::value
656  >
657 class ScatterView;
658 
659 template <typename DataType
660  ,int Op
661  ,typename ExecSpace
662  ,typename Layout
663  ,int duplication
664  ,int contribution
665  ,int override_contribution
666  >
667 class ScatterAccess;
668 
669 // non-duplicated implementation
670 template <typename DataType
671  ,int Op
672  ,typename ExecSpace
673  ,typename Layout
674  ,int contribution
675  >
676 class ScatterView<DataType
677  ,Layout
678  ,ExecSpace
679  ,Op
680  ,ScatterNonDuplicated
681  ,contribution>
682 {
683 public:
684  typedef Kokkos::View<DataType, Layout, ExecSpace> original_view_type;
685  typedef typename original_view_type::value_type original_value_type;
686  typedef typename original_view_type::reference_type original_reference_type;
687  friend class ScatterAccess<DataType, Op, ExecSpace, Layout, ScatterNonDuplicated, contribution, ScatterNonAtomic>;
688  friend class ScatterAccess<DataType, Op, ExecSpace, Layout, ScatterNonDuplicated, contribution, ScatterAtomic>;
689 
690  ScatterView()
691  {
692  }
693 
694  template <typename RT, typename ... RP>
695  ScatterView(View<RT, RP...> const& original_view)
696  : internal_view(original_view)
697  {
698  }
699 
700  template <typename ... Dims>
701  ScatterView(std::string const& name, Dims ... dims)
702  : internal_view(name, dims ...)
703  {
704  }
705 
706  template <int override_contrib = contribution>
707  KOKKOS_FORCEINLINE_FUNCTION
708  ScatterAccess<DataType, Op, ExecSpace, Layout, ScatterNonDuplicated, contribution, override_contrib>
709  access() const {
710  return ScatterAccess<DataType, Op, ExecSpace, Layout, ScatterNonDuplicated, contribution, override_contrib>{*this};
711  }
712 
713  original_view_type subview() const {
714  return internal_view;
715  }
716 
717  template <typename DT, typename ... RP>
718  void contribute_into(View<DT, RP...> const& dest) const
719  {
720  typedef View<DT, RP...> dest_type;
721  static_assert(std::is_same<
722  typename dest_type::array_layout,
723  Layout>::value,
724  "ScatterView contribute destination has different layout");
725  static_assert(Kokkos::Impl::VerifyExecutionCanAccessMemorySpace<
726  typename ExecSpace::memory_space,
727  typename dest_type::memory_space>::value,
728  "ScatterView contribute destination memory space not accessible");
729  if (dest.data() == internal_view.data()) return;
730  Kokkos::Impl::Experimental::ReduceDuplicates<ExecSpace, original_value_type, Op>(
731  internal_view.data(),
732  dest.data(),
733  0,
734  0,
735  1,
736  internal_view.label());
737  }
738 
739  void reset() {
740  Kokkos::Impl::Experimental::ResetDuplicates<ExecSpace, original_value_type, Op>(
741  internal_view.data(),
742  internal_view.size(),
743  internal_view.label());
744  }
745  template <typename DT, typename ... RP>
746  void reset_except(View<DT, RP...> const& view) {
747  if (view.data() != internal_view.data()) reset();
748  }
749 
750  void resize(const size_t n0 = 0,
751  const size_t n1 = 0,
752  const size_t n2 = 0,
753  const size_t n3 = 0,
754  const size_t n4 = 0,
755  const size_t n5 = 0,
756  const size_t n6 = 0,
757  const size_t n7 = 0) {
758  ::Kokkos::resize(internal_view,n0,n1,n2,n3,n4,n5,n6,n7);
759  }
760 
761  void realloc(const size_t n0 = 0,
762  const size_t n1 = 0,
763  const size_t n2 = 0,
764  const size_t n3 = 0,
765  const size_t n4 = 0,
766  const size_t n5 = 0,
767  const size_t n6 = 0,
768  const size_t n7 = 0) {
769  ::Kokkos::realloc(internal_view,n0,n1,n2,n3,n4,n5,n6,n7);
770  }
771 
772 protected:
773  template <typename ... Args>
774  KOKKOS_FORCEINLINE_FUNCTION
775  original_reference_type at(Args ... args) const {
776  return internal_view(args...);
777  }
778 private:
779  typedef original_view_type internal_view_type;
780  internal_view_type internal_view;
781 };
782 
783 template <typename DataType
784  ,int Op
785  ,typename ExecSpace
786  ,typename Layout
787  ,int contribution
788  ,int override_contribution
789  >
790 class ScatterAccess<DataType
791  ,Op
792  ,ExecSpace
793  ,Layout
794  ,ScatterNonDuplicated
795  ,contribution
796  ,override_contribution>
797 {
798 public:
799  typedef ScatterView<DataType, Layout, ExecSpace, Op, ScatterNonDuplicated, contribution> view_type;
800  typedef typename view_type::original_value_type original_value_type;
801  typedef Kokkos::Impl::Experimental::ScatterValue<
802  original_value_type, Op, override_contribution> value_type;
803 
804  KOKKOS_INLINE_FUNCTION
805  ScatterAccess() :
806  view(view_type()) {
807  }
808 
809  KOKKOS_INLINE_FUNCTION
810  ScatterAccess(view_type const& view_in)
811  : view(view_in)
812  {
813  }
814 
815  KOKKOS_INLINE_FUNCTION
816  ~ScatterAccess()
817  {
818  }
819 
820  template <typename ... Args>
821  KOKKOS_FORCEINLINE_FUNCTION
822  value_type operator()(Args ... args) const {
823  return view.at(args...);
824  }
825 
826  template <typename Arg>
827  KOKKOS_FORCEINLINE_FUNCTION
828  typename std::enable_if<view_type::original_view_type::rank == 1 &&
829  std::is_integral<Arg>::value, value_type>::type
830  operator[](Arg arg) const {
831  return view.at(arg);
832  }
833 
834 private:
835  view_type const& view;
836 };
837 
838 // duplicated implementation
839 // LayoutLeft and LayoutRight are different enough that we'll just specialize each
840 
841 template <typename DataType
842  ,int Op
843  ,typename ExecSpace
844  ,int contribution
845  >
846 class ScatterView<DataType
847  ,Kokkos::LayoutRight
848  ,ExecSpace
849  ,Op
850  ,ScatterDuplicated
851  ,contribution>
852 {
853 public:
855  typedef typename original_view_type::value_type original_value_type;
856  typedef typename original_view_type::reference_type original_reference_type;
857  friend class ScatterAccess<DataType, Op, ExecSpace, Kokkos::LayoutRight, ScatterDuplicated, contribution, ScatterNonAtomic>;
858  friend class ScatterAccess<DataType, Op, ExecSpace, Kokkos::LayoutRight, ScatterDuplicated, contribution, ScatterAtomic>;
859  typedef typename Kokkos::Impl::Experimental::DuplicatedDataType<DataType, Kokkos::LayoutRight> data_type_info;
860  typedef typename data_type_info::value_type internal_data_type;
861  typedef Kokkos::View<internal_data_type, Kokkos::LayoutRight, ExecSpace> internal_view_type;
862 
863  ScatterView()
864  {
865  }
866 
867  template <typename RT, typename ... RP >
868  ScatterView(View<RT, RP...> const& original_view)
869  : unique_token()
870  , internal_view(Kokkos::ViewAllocateWithoutInitializing(
871  std::string("duplicated_") + original_view.label()),
872  unique_token.size(),
873 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
874  original_view.extent(0),
875  original_view.extent(1),
876  original_view.extent(2),
877  original_view.extent(3),
878  original_view.extent(4),
879  original_view.extent(5),
880  original_view.extent(6) )
881 #else
882  original_view.rank_dynamic > 0 ? original_view.extent(0): KOKKOS_IMPL_CTOR_DEFAULT_ARG,
883  original_view.rank_dynamic > 1 ? original_view.extent(1): KOKKOS_IMPL_CTOR_DEFAULT_ARG,
884  original_view.rank_dynamic > 2 ? original_view.extent(2): KOKKOS_IMPL_CTOR_DEFAULT_ARG,
885  original_view.rank_dynamic > 3 ? original_view.extent(3): KOKKOS_IMPL_CTOR_DEFAULT_ARG,
886  original_view.rank_dynamic > 4 ? original_view.extent(4): KOKKOS_IMPL_CTOR_DEFAULT_ARG,
887  original_view.rank_dynamic > 5 ? original_view.extent(5): KOKKOS_IMPL_CTOR_DEFAULT_ARG,
888  original_view.rank_dynamic > 6 ? original_view.extent(6): KOKKOS_IMPL_CTOR_DEFAULT_ARG)
889 
890 #endif
891  {
892  reset();
893  }
894 
895  template <typename ... Dims>
896  ScatterView(std::string const& name, Dims ... dims)
897  : internal_view(Kokkos::ViewAllocateWithoutInitializing(name), unique_token.size(), dims ...)
898  {
899  reset();
900  }
901 
902  template <int override_contribution = contribution>
903  KOKKOS_FORCEINLINE_FUNCTION
904  ScatterAccess<DataType, Op, ExecSpace, Kokkos::LayoutRight, ScatterDuplicated, contribution, override_contribution>
905  access() const {
906  return ScatterAccess<DataType, Op, ExecSpace, Kokkos::LayoutRight, ScatterDuplicated, contribution, override_contribution>{*this};
907  }
908 
909  typename Kokkos::Impl::Experimental::Slice<
910  Kokkos::LayoutRight, internal_view_type::rank, internal_view_type>::value_type
911  subview() const
912  {
913  return Kokkos::Impl::Experimental::Slice<
914  Kokkos::LayoutRight, internal_view_type::Rank, internal_view_type>::get(internal_view, 0);
915  }
916 
917  template <typename DT, typename ... RP>
918  void contribute_into(View<DT, RP...> const& dest) const
919  {
920  typedef View<DT, RP...> dest_type;
921  static_assert(std::is_same<
922  typename dest_type::array_layout,
923  Kokkos::LayoutRight>::value,
924  "ScatterView deep_copy destination has different layout");
925  static_assert(Kokkos::Impl::VerifyExecutionCanAccessMemorySpace<
926  typename ExecSpace::memory_space,
927  typename dest_type::memory_space>::value,
928  "ScatterView deep_copy destination memory space not accessible");
929  bool is_equal = (dest.data() == internal_view.data());
930  size_t start = is_equal ? 1 : 0;
931  Kokkos::Impl::Experimental::ReduceDuplicates<ExecSpace, original_value_type, Op>(
932  internal_view.data(),
933  dest.data(),
934  internal_view.stride(0),
935  start,
936  internal_view.extent(0),
937  internal_view.label());
938  }
939 
940  void reset() {
941  Kokkos::Impl::Experimental::ResetDuplicates<ExecSpace, original_value_type, Op>(
942  internal_view.data(),
943  internal_view.size(),
944  internal_view.label());
945  }
946  template <typename DT, typename ... RP>
947  void reset_except(View<DT, RP...> const& view) {
948  if (view.data() != internal_view.data()) {
949  reset();
950  return;
951  }
952  Kokkos::Impl::Experimental::ResetDuplicates<ExecSpace, original_value_type, Op>(
953  internal_view.data() + view.size(),
954  internal_view.size() - view.size(),
955  internal_view.label());
956  }
957 
958  void resize(const size_t n0 = 0,
959  const size_t n1 = 0,
960  const size_t n2 = 0,
961  const size_t n3 = 0,
962  const size_t n4 = 0,
963  const size_t n5 = 0,
964  const size_t n6 = 0) {
965  ::Kokkos::resize(internal_view,unique_token.size(),n0,n1,n2,n3,n4,n5,n6);
966  }
967 
968  void realloc(const size_t n0 = 0,
969  const size_t n1 = 0,
970  const size_t n2 = 0,
971  const size_t n3 = 0,
972  const size_t n4 = 0,
973  const size_t n5 = 0,
974  const size_t n6 = 0) {
975  ::Kokkos::realloc(internal_view,unique_token.size(),n0,n1,n2,n3,n4,n5,n6);
976  }
977 
978 protected:
979  template <typename ... Args>
980  KOKKOS_FORCEINLINE_FUNCTION
981  original_reference_type at(int rank, Args ... args) const {
982  return internal_view(rank, args...);
983  }
984 
985 protected:
987  ExecSpace, Kokkos::Experimental::UniqueTokenScope::Global> unique_token_type;
988 
989  unique_token_type unique_token;
990  internal_view_type internal_view;
991 };
992 
993 template <typename DataType
994  ,int Op
995  ,typename ExecSpace
996  ,int contribution
997  >
998 class ScatterView<DataType
999  ,Kokkos::LayoutLeft
1000  ,ExecSpace
1001  ,Op
1002  ,ScatterDuplicated
1003  ,contribution>
1004 {
1005 public:
1006  typedef Kokkos::View<DataType, Kokkos::LayoutLeft, ExecSpace> original_view_type;
1007  typedef typename original_view_type::value_type original_value_type;
1008  typedef typename original_view_type::reference_type original_reference_type;
1009  friend class ScatterAccess<DataType, Op, ExecSpace, Kokkos::LayoutLeft, ScatterDuplicated, contribution, ScatterNonAtomic>;
1010  friend class ScatterAccess<DataType, Op, ExecSpace, Kokkos::LayoutLeft, ScatterDuplicated, contribution, ScatterAtomic>;
1011  typedef typename Kokkos::Impl::Experimental::DuplicatedDataType<DataType, Kokkos::LayoutLeft> data_type_info;
1012  typedef typename data_type_info::value_type internal_data_type;
1013  typedef Kokkos::View<internal_data_type, Kokkos::LayoutLeft, ExecSpace> internal_view_type;
1014 
1015  ScatterView()
1016  {
1017  }
1018 
1019  template <typename RT, typename ... RP >
1020  ScatterView(View<RT, RP...> const& original_view)
1021  : unique_token()
1022  {
1023  size_t arg_N[8] = {
1024  original_view.rank>0?original_view.extent(0):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1025  original_view.rank>1?original_view.extent(1):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1026  original_view.rank>2?original_view.extent(2):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1027  original_view.rank>3?original_view.extent(3):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1028  original_view.rank>4?original_view.extent(4):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1029  original_view.rank>5?original_view.extent(5):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1030  original_view.rank>6?original_view.extent(6):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1031  KOKKOS_IMPL_CTOR_DEFAULT_ARG
1032  };
1033  arg_N[internal_view_type::rank - 1] = unique_token.size();
1034  internal_view = internal_view_type(
1035  Kokkos::ViewAllocateWithoutInitializing(
1036  std::string("duplicated_") + original_view.label()),
1037  arg_N[0], arg_N[1], arg_N[2], arg_N[3],
1038  arg_N[4], arg_N[5], arg_N[6], arg_N[7]);
1039  reset();
1040  }
1041 
1042  template <typename ... Dims>
1043  ScatterView(std::string const& name, Dims ... dims) {
1044  original_view_type original_view;
1045  size_t arg_N[8] = {
1046  original_view.rank>0?original_view.static_extent(0):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1047  original_view.rank>1?original_view.static_extent(1):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1048  original_view.rank>2?original_view.static_extent(2):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1049  original_view.rank>3?original_view.static_extent(3):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1050  original_view.rank>4?original_view.static_extent(4):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1051  original_view.rank>5?original_view.static_extent(5):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1052  original_view.rank>6?original_view.static_extent(6):KOKKOS_IMPL_CTOR_DEFAULT_ARG,
1053  KOKKOS_IMPL_CTOR_DEFAULT_ARG
1054  };
1055  Kokkos::Impl::Experimental::args_to_array(arg_N,0,dims ...);
1056  arg_N[internal_view_type::rank - 1] = unique_token.size();
1057  internal_view = internal_view_type(Kokkos::ViewAllocateWithoutInitializing(name),
1058  arg_N[0], arg_N[1], arg_N[2], arg_N[3],
1059  arg_N[4], arg_N[5], arg_N[6], arg_N[7]);
1060  reset();
1061  }
1062 
1063  template <int override_contribution = contribution>
1064  KOKKOS_FORCEINLINE_FUNCTION
1065  ScatterAccess<DataType, Op, ExecSpace, Kokkos::LayoutLeft, ScatterDuplicated, contribution, override_contribution>
1066  access() const {
1067  return ScatterAccess<DataType, Op, ExecSpace, Kokkos::LayoutLeft, ScatterDuplicated, contribution, override_contribution>{*this};
1068  }
1069 
1070  typename Kokkos::Impl::Experimental::Slice<
1071  Kokkos::LayoutLeft, internal_view_type::rank, internal_view_type>::value_type
1072  subview() const
1073  {
1074  return Kokkos::Impl::Experimental::Slice<
1075  Kokkos::LayoutLeft, internal_view_type::rank, internal_view_type>::get(internal_view, 0);
1076  }
1077 
1078  template <typename ... RP>
1079  void contribute_into(View<RP...> const& dest) const
1080  {
1081  typedef View<RP...> dest_type;
1082  static_assert(std::is_same<
1083  typename dest_type::value_type,
1084  typename original_view_type::non_const_value_type>::value,
1085  "ScatterView deep_copy destination has wrong value_type");
1086  static_assert(std::is_same<
1087  typename dest_type::array_layout,
1088  Kokkos::LayoutLeft>::value,
1089  "ScatterView deep_copy destination has different layout");
1090  static_assert(Kokkos::Impl::VerifyExecutionCanAccessMemorySpace<
1091  typename ExecSpace::memory_space,
1092  typename dest_type::memory_space>::value,
1093  "ScatterView deep_copy destination memory space not accessible");
1094  auto extent = internal_view.extent(
1095  internal_view_type::rank - 1);
1096  bool is_equal = (dest.data() == internal_view.data());
1097  size_t start = is_equal ? 1 : 0;
1098  Kokkos::Impl::Experimental::ReduceDuplicates<ExecSpace, original_value_type, Op>(
1099  internal_view.data(),
1100  dest.data(),
1101  internal_view.stride(internal_view_type::rank - 1),
1102  start,
1103  extent,
1104  internal_view.label());
1105  }
1106 
1107  void reset() {
1108  Kokkos::Impl::Experimental::ResetDuplicates<ExecSpace, original_value_type, Op>(
1109  internal_view.data(),
1110  internal_view.size(),
1111  internal_view.label());
1112  }
1113  template <typename DT, typename ... RP>
1114  void reset_except(View<DT, RP...> const& view) {
1115  if (view.data() != internal_view.data()) {
1116  reset();
1117  return;
1118  }
1119  Kokkos::Impl::Experimental::ResetDuplicates<ExecSpace, original_value_type, Op>(
1120  internal_view.data() + view.size(),
1121  internal_view.size() - view.size(),
1122  internal_view.label());
1123  }
1124 
1125  void resize(const size_t n0 = 0,
1126  const size_t n1 = 0,
1127  const size_t n2 = 0,
1128  const size_t n3 = 0,
1129  const size_t n4 = 0,
1130  const size_t n5 = 0,
1131  const size_t n6 = 0) {
1132 
1133  size_t arg_N[8] = {n0,n1,n2,n3,n4,n5,n6,0};
1134  const int i = internal_view.rank-1;
1135  arg_N[i] = unique_token.size();
1136 
1137  ::Kokkos::resize(internal_view,
1138  arg_N[0], arg_N[1], arg_N[2], arg_N[3],
1139  arg_N[4], arg_N[5], arg_N[6], arg_N[7]);
1140  }
1141 
1142  void realloc(const size_t n0 = 0,
1143  const size_t n1 = 0,
1144  const size_t n2 = 0,
1145  const size_t n3 = 0,
1146  const size_t n4 = 0,
1147  const size_t n5 = 0,
1148  const size_t n6 = 0) {
1149 
1150  size_t arg_N[8] = {n0,n1,n2,n3,n4,n5,n6,0};
1151  const int i = internal_view.rank-1;
1152  arg_N[i] = unique_token.size();
1153 
1154  ::Kokkos::realloc(internal_view,
1155  arg_N[0], arg_N[1], arg_N[2], arg_N[3],
1156  arg_N[4], arg_N[5], arg_N[6], arg_N[7]);
1157  }
1158 
1159 protected:
1160  template <typename ... Args>
1161  inline original_reference_type at(int thread_id, Args ... args) const {
1162  return internal_view(args..., thread_id);
1163  }
1164 
1165 protected:
1167  ExecSpace, Kokkos::Experimental::UniqueTokenScope::Global> unique_token_type;
1168 
1169  unique_token_type unique_token;
1170  internal_view_type internal_view;
1171 };
1172 
1173 
1174 /* This object has to be separate in order to store the thread ID, which cannot
1175  be obtained until one is inside a parallel construct, and may be relatively
1176  expensive to obtain at every contribution
1177  (calls a non-inlined function, looks up a thread-local variable).
1178  Due to the expense, it is sensible to query it at most once per parallel iterate
1179  (ideally once per thread, but parallel_for doesn't expose that)
1180  and then store it in a stack variable.
1181  ScatterAccess serves as a non-const object on the stack which can store the thread ID */
1182 
1183 template <typename DataType
1184  ,int Op
1185  ,typename ExecSpace
1186  ,typename Layout
1187  ,int contribution
1188  ,int override_contribution
1189  >
1190 class ScatterAccess<DataType
1191  ,Op
1192  ,ExecSpace
1193  ,Layout
1194  ,ScatterDuplicated
1195  ,contribution
1196  ,override_contribution>
1197 {
1198 public:
1199  typedef ScatterView<DataType, Layout, ExecSpace, Op, ScatterDuplicated, contribution> view_type;
1200  typedef typename view_type::original_value_type original_value_type;
1201  typedef Kokkos::Impl::Experimental::ScatterValue<
1202  original_value_type, Op, override_contribution> value_type;
1203 
1204  KOKKOS_FORCEINLINE_FUNCTION
1205  ScatterAccess(view_type const& view_in)
1206  : view(view_in)
1207  , thread_id(view_in.unique_token.acquire()) {
1208  }
1209 
1210  KOKKOS_FORCEINLINE_FUNCTION
1211  ~ScatterAccess() {
1212  if (thread_id != ~thread_id_type(0)) view.unique_token.release(thread_id);
1213  }
1214 
1215  template <typename ... Args>
1216  KOKKOS_FORCEINLINE_FUNCTION
1217  value_type operator()(Args ... args) const {
1218  return view.at(thread_id, args...);
1219  }
1220 
1221  template <typename Arg>
1222  KOKKOS_FORCEINLINE_FUNCTION
1223  typename std::enable_if<view_type::original_view_type::rank == 1 &&
1224  std::is_integral<Arg>::value, value_type>::type
1225  operator[](Arg arg) const {
1226  return view.at(thread_id, arg);
1227  }
1228 
1229 private:
1230 
1231  view_type const& view;
1232 
1233  // simplify RAII by disallowing copies
1234  ScatterAccess(ScatterAccess const& other) = delete;
1235  ScatterAccess& operator=(ScatterAccess const& other) = delete;
1236  ScatterAccess& operator=(ScatterAccess&& other) = delete;
1237 
1238 public:
1239  // do need to allow moves though, for the common
1240  // auto b = a.access();
1241  // that assignments turns into a move constructor call
1242  KOKKOS_FORCEINLINE_FUNCTION
1243  ScatterAccess(ScatterAccess&& other)
1244  : view(other.view)
1245  , thread_id(other.thread_id)
1246  {
1247  other.thread_id = ~thread_id_type(0);
1248  }
1249 
1250 private:
1251 
1252  typedef typename view_type::unique_token_type unique_token_type;
1253  typedef typename unique_token_type::size_type thread_id_type;
1254  thread_id_type thread_id;
1255 };
1256 
1257 template <int Op = Kokkos::Experimental::ScatterSum,
1258  int duplication = -1,
1259  int contribution = -1,
1260  typename RT, typename ... RP>
1261 ScatterView
1262  < RT
1263  , typename ViewTraits<RT, RP...>::array_layout
1264  , typename ViewTraits<RT, RP...>::execution_space
1265  , Op
1266  /* just setting defaults if not specified... things got messy because the view type
1267  does not come before the duplication/contribution settings in the
1268  template parameter list */
1269  , duplication == -1 ? Kokkos::Impl::Experimental::DefaultDuplication<typename ViewTraits<RT, RP...>::execution_space>::value : duplication
1270  , contribution == -1 ?
1271  Kokkos::Impl::Experimental::DefaultContribution<
1272  typename ViewTraits<RT, RP...>::execution_space,
1273  (duplication == -1 ?
1274  Kokkos::Impl::Experimental::DefaultDuplication<
1275  typename ViewTraits<RT, RP...>::execution_space
1276  >::value
1277  : duplication
1278  )
1279  >::value
1280  : contribution
1281  >
1282 create_scatter_view(View<RT, RP...> const& original_view) {
1283  return original_view; // implicit ScatterView constructor call
1284 }
1285 
1286 }} // namespace Kokkos::Experimental
1287 
1288 namespace Kokkos {
1289 namespace Experimental {
1290 
1291 template <typename DT1, typename DT2, typename LY, typename ES, int OP, int CT, int DP, typename ... VP>
1292 void
1293 contribute(View<DT1, VP...>& dest, Kokkos::Experimental::ScatterView<DT2, LY, ES, OP, CT, DP> const& src)
1294 {
1295  src.contribute_into(dest);
1296 }
1297 
1298 }} // namespace Kokkos::Experimental
1299 
1300 namespace Kokkos {
1301 
1302 template <typename DT, typename LY, typename ES, int OP, int CT, int DP, typename ... IS>
1303 void
1304 realloc(Kokkos::Experimental::ScatterView<DT, LY, ES, OP, CT, DP>& scatter_view, IS ... is)
1305 {
1306  scatter_view.realloc(is ...);
1307 }
1308 
1309 template <typename DT, typename LY, typename ES, int OP, int CT, int DP, typename ... IS>
1310 void
1311 resize(Kokkos::Experimental::ScatterView<DT, LY, ES, OP, CT, DP>& scatter_view, IS ... is)
1312 {
1313  scatter_view.resize(is ...);
1314 }
1315 
1316 } // namespace Kokkos
1317 
1318 #endif
Memory layout tag indicating left-to-right (Fortran scheme) striding of multi-indices.
class to generate unique ids base on the required amount of concurrency
std::enable_if< std::is_same< typename Kokkos::View< T, P...>::array_layout, Kokkos::LayoutLeft >::value||std::is_same< typename Kokkos::View< T, P...>::array_layout, Kokkos::LayoutRight >::value >::type resize(Kokkos::View< T, P...> &v, const size_t n0=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n1=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n2=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n3=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n4=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n5=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n6=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n7=KOKKOS_IMPL_CTOR_DEFAULT_ARG)
Resize a view with copying old data to new data at the corresponding indices.
View to an array of data.
Memory layout tag indicating right-to-left (C or lexigraphical scheme) striding of multi-indices...
void resize(DynRankView< T, P...> &v, const size_t n0=KOKKOS_INVALID_INDEX, const size_t n1=KOKKOS_INVALID_INDEX, const size_t n2=KOKKOS_INVALID_INDEX, const size_t n3=KOKKOS_INVALID_INDEX, const size_t n4=KOKKOS_INVALID_INDEX, const size_t n5=KOKKOS_INVALID_INDEX, const size_t n6=KOKKOS_INVALID_INDEX, const size_t n7=KOKKOS_INVALID_INDEX)
Resize a view with copying old data to new data at the corresponding indices.
Implementation of the ParallelFor operator that has a partial specialization for the device...
std::enable_if< std::is_same< typename Kokkos::View< T, P...>::array_layout, Kokkos::LayoutLeft >::value||std::is_same< typename Kokkos::View< T, P...>::array_layout, Kokkos::LayoutRight >::value >::type realloc(Kokkos::View< T, P...> &v, const size_t n0=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n1=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n2=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n3=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n4=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n5=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n6=KOKKOS_IMPL_CTOR_DEFAULT_ARG, const size_t n7=KOKKOS_IMPL_CTOR_DEFAULT_ARG)
Resize a view with discarding old data.
KOKKOS_INLINE_FUNCTION constexpr unsigned rank(const View< D, P...> &V)
Temporary free function rank() until rank() is implemented in the View.
void realloc(DynRankView< T, P...> &v, const size_t n0=KOKKOS_INVALID_INDEX, const size_t n1=KOKKOS_INVALID_INDEX, const size_t n2=KOKKOS_INVALID_INDEX, const size_t n3=KOKKOS_INVALID_INDEX, const size_t n4=KOKKOS_INVALID_INDEX, const size_t n5=KOKKOS_INVALID_INDEX, const size_t n6=KOKKOS_INVALID_INDEX, const size_t n7=KOKKOS_INVALID_INDEX)
Resize a view with copying old data to new data at the corresponding indices.