Kokkos Core Kernels Package Version of the Day
Loading...
Searching...
No Matches
Kokkos_Parallel_Reduce.hpp
1//@HEADER
2// ************************************************************************
3//
4// Kokkos v. 4.0
5// Copyright (2022) National Technology & Engineering
6// Solutions of Sandia, LLC (NTESS).
7//
8// Under the terms of Contract DE-NA0003525 with NTESS,
9// the U.S. Government retains certain rights in this software.
10//
11// Part of Kokkos, under the Apache License v2.0 with LLVM Exceptions.
12// See https://kokkos.org/LICENSE for license information.
13// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
14//
15//@HEADER
16
17#ifndef KOKKOS_IMPL_PUBLIC_INCLUDE
18#include <Kokkos_Macros.hpp>
19static_assert(false,
20 "Including non-public Kokkos header files is not allowed.");
21#endif
22#ifndef KOKKOS_PARALLEL_REDUCE_HPP
23#define KOKKOS_PARALLEL_REDUCE_HPP
24
25#include <Kokkos_ReductionIdentity.hpp>
26#include <Kokkos_View.hpp>
27#include <impl/Kokkos_FunctorAnalysis.hpp>
28#include <impl/Kokkos_Tools_Generic.hpp>
29#include <type_traits>
30#include <iostream>
31
32namespace Kokkos {
33
34template <class Scalar, class Space>
35struct Sum {
36 public:
37 // Required
38 using reducer = Sum<Scalar, Space>;
39 using value_type = std::remove_cv_t<Scalar>;
40 static_assert(!std::is_pointer_v<value_type> && !std::is_array_v<value_type>);
41
42 using result_view_type = Kokkos::View<value_type, Space>;
43
44 private:
45 result_view_type value;
46 bool references_scalar_v;
47
48 public:
49 KOKKOS_INLINE_FUNCTION
50 Sum(value_type& value_) : value(&value_), references_scalar_v(true) {}
51
52 KOKKOS_INLINE_FUNCTION
53 Sum(const result_view_type& value_)
54 : value(value_), references_scalar_v(false) {}
55
56 // Required
57 KOKKOS_INLINE_FUNCTION
58 void join(value_type& dest, const value_type& src) const { dest += src; }
59
60 KOKKOS_INLINE_FUNCTION
61 void init(value_type& val) const {
62 val = reduction_identity<value_type>::sum();
63 }
64
65 KOKKOS_INLINE_FUNCTION
66 value_type& reference() const { return *value.data(); }
67
68 KOKKOS_INLINE_FUNCTION
69 result_view_type view() const { return value; }
70
71 KOKKOS_INLINE_FUNCTION
72 bool references_scalar() const { return references_scalar_v; }
73};
74
75template <typename Scalar, typename... Properties>
77 ->Sum<Scalar, typename View<Scalar, Properties...>::memory_space>;
78
79template <class Scalar, class Space>
80struct Prod {
81 public:
82 // Required
83 using reducer = Prod<Scalar, Space>;
84 using value_type = std::remove_cv_t<Scalar>;
85 static_assert(!std::is_pointer_v<value_type> && !std::is_array_v<value_type>);
86
87 using result_view_type = Kokkos::View<value_type, Space>;
88
89 private:
90 result_view_type value;
91 bool references_scalar_v;
92
93 public:
94 KOKKOS_INLINE_FUNCTION
95 Prod(value_type& value_) : value(&value_), references_scalar_v(true) {}
96
97 KOKKOS_INLINE_FUNCTION
98 Prod(const result_view_type& value_)
99 : value(value_), references_scalar_v(false) {}
100
101 // Required
102 KOKKOS_INLINE_FUNCTION
103 void join(value_type& dest, const value_type& src) const { dest *= src; }
104
105 KOKKOS_INLINE_FUNCTION
106 void init(value_type& val) const {
107 val = reduction_identity<value_type>::prod();
108 }
109
110 KOKKOS_INLINE_FUNCTION
111 value_type& reference() const { return *value.data(); }
112
113 KOKKOS_INLINE_FUNCTION
114 result_view_type view() const { return value; }
115
116 KOKKOS_INLINE_FUNCTION
117 bool references_scalar() const { return references_scalar_v; }
118};
119
120template <typename Scalar, typename... Properties>
122 ->Prod<Scalar, typename View<Scalar, Properties...>::memory_space>;
123
124template <class Scalar, class Space>
125struct Min {
126 public:
127 // Required
128 using reducer = Min<Scalar, Space>;
129 using value_type = std::remove_cv_t<Scalar>;
130 static_assert(!std::is_pointer_v<value_type> && !std::is_array_v<value_type>);
131
132 using result_view_type = Kokkos::View<value_type, Space>;
133
134 private:
135 result_view_type value;
136 bool references_scalar_v;
137
138 public:
139 KOKKOS_INLINE_FUNCTION
140 Min(value_type& value_) : value(&value_), references_scalar_v(true) {}
141
142 KOKKOS_INLINE_FUNCTION
143 Min(const result_view_type& value_)
144 : value(value_), references_scalar_v(false) {}
145
146 // Required
147 KOKKOS_INLINE_FUNCTION
148 void join(value_type& dest, const value_type& src) const {
149 if (src < dest) dest = src;
150 }
151
152 KOKKOS_INLINE_FUNCTION
153 void init(value_type& val) const {
154 val = reduction_identity<value_type>::min();
155 }
156
157 KOKKOS_INLINE_FUNCTION
158 value_type& reference() const { return *value.data(); }
159
160 KOKKOS_INLINE_FUNCTION
161 result_view_type view() const { return value; }
162
163 KOKKOS_INLINE_FUNCTION
164 bool references_scalar() const { return references_scalar_v; }
165};
166
167template <typename Scalar, typename... Properties>
169 ->Min<Scalar, typename View<Scalar, Properties...>::memory_space>;
170
171template <class Scalar, class Space>
172struct Max {
173 public:
174 // Required
175 using reducer = Max<Scalar, Space>;
176 using value_type = std::remove_cv_t<Scalar>;
177 static_assert(!std::is_pointer_v<value_type> && !std::is_array_v<value_type>);
178
179 using result_view_type = Kokkos::View<value_type, Space>;
180
181 private:
182 result_view_type value;
183 bool references_scalar_v;
184
185 public:
186 KOKKOS_INLINE_FUNCTION
187 Max(value_type& value_) : value(&value_), references_scalar_v(true) {}
188
189 KOKKOS_INLINE_FUNCTION
190 Max(const result_view_type& value_)
191 : value(value_), references_scalar_v(false) {}
192
193 // Required
194 KOKKOS_INLINE_FUNCTION
195 void join(value_type& dest, const value_type& src) const {
196 if (src > dest) dest = src;
197 }
198
199 // Required
200 KOKKOS_INLINE_FUNCTION
201 void init(value_type& val) const {
202 val = reduction_identity<value_type>::max();
203 }
204
205 KOKKOS_INLINE_FUNCTION
206 value_type& reference() const { return *value.data(); }
207
208 KOKKOS_INLINE_FUNCTION
209 result_view_type view() const { return value; }
210
211 KOKKOS_INLINE_FUNCTION
212 bool references_scalar() const { return references_scalar_v; }
213};
214
215template <typename Scalar, typename... Properties>
217 ->Max<Scalar, typename View<Scalar, Properties...>::memory_space>;
218
219template <class Scalar, class Space>
220struct LAnd {
221 public:
222 // Required
223 using reducer = LAnd<Scalar, Space>;
224 using value_type = std::remove_cv_t<Scalar>;
225 static_assert(!std::is_pointer_v<value_type> && !std::is_array_v<value_type>);
226
227 using result_view_type = Kokkos::View<value_type, Space>;
228
229 private:
230 result_view_type value;
231 bool references_scalar_v;
232
233 public:
234 KOKKOS_INLINE_FUNCTION
235 LAnd(value_type& value_) : value(&value_), references_scalar_v(true) {}
236
237 KOKKOS_INLINE_FUNCTION
238 LAnd(const result_view_type& value_)
239 : value(value_), references_scalar_v(false) {}
240
241 KOKKOS_INLINE_FUNCTION
242 void join(value_type& dest, const value_type& src) const {
243 dest = dest && src;
244 }
245
246 KOKKOS_INLINE_FUNCTION
247 void init(value_type& val) const {
248 val = reduction_identity<value_type>::land();
249 }
250
251 KOKKOS_INLINE_FUNCTION
252 value_type& reference() const { return *value.data(); }
253
254 KOKKOS_INLINE_FUNCTION
255 result_view_type view() const { return value; }
256
257 KOKKOS_INLINE_FUNCTION
258 bool references_scalar() const { return references_scalar_v; }
259};
260
261template <typename Scalar, typename... Properties>
263 ->LAnd<Scalar, typename View<Scalar, Properties...>::memory_space>;
264
265template <class Scalar, class Space>
266struct LOr {
267 public:
268 // Required
269 using reducer = LOr<Scalar, Space>;
270 using value_type = std::remove_cv_t<Scalar>;
271 static_assert(!std::is_pointer_v<value_type> && !std::is_array_v<value_type>);
272
273 using result_view_type = Kokkos::View<value_type, Space>;
274
275 private:
276 result_view_type value;
277 bool references_scalar_v;
278
279 public:
280 KOKKOS_INLINE_FUNCTION
281 LOr(value_type& value_) : value(&value_), references_scalar_v(true) {}
282
283 KOKKOS_INLINE_FUNCTION
284 LOr(const result_view_type& value_)
285 : value(value_), references_scalar_v(false) {}
286
287 // Required
288 KOKKOS_INLINE_FUNCTION
289 void join(value_type& dest, const value_type& src) const {
290 dest = dest || src;
291 }
292
293 KOKKOS_INLINE_FUNCTION
294 void init(value_type& val) const {
295 val = reduction_identity<value_type>::lor();
296 }
297
298 KOKKOS_INLINE_FUNCTION
299 value_type& reference() const { return *value.data(); }
300
301 KOKKOS_INLINE_FUNCTION
302 result_view_type view() const { return value; }
303
304 KOKKOS_INLINE_FUNCTION
305 bool references_scalar() const { return references_scalar_v; }
306};
307
308template <typename Scalar, typename... Properties>
310 ->LOr<Scalar, typename View<Scalar, Properties...>::memory_space>;
311
312template <class Scalar, class Space>
313struct BAnd {
314 public:
315 // Required
316 using reducer = BAnd<Scalar, Space>;
317 using value_type = std::remove_cv_t<Scalar>;
318 static_assert(!std::is_pointer_v<value_type> && !std::is_array_v<value_type>);
319
320 using result_view_type = Kokkos::View<value_type, Space>;
321
322 private:
323 result_view_type value;
324 bool references_scalar_v;
325
326 public:
327 KOKKOS_INLINE_FUNCTION
328 BAnd(value_type& value_) : value(&value_), references_scalar_v(true) {}
329
330 KOKKOS_INLINE_FUNCTION
331 BAnd(const result_view_type& value_)
332 : value(value_), references_scalar_v(false) {}
333
334 // Required
335 KOKKOS_INLINE_FUNCTION
336 void join(value_type& dest, const value_type& src) const {
337 dest = dest & src;
338 }
339
340 KOKKOS_INLINE_FUNCTION
341 void init(value_type& val) const {
342 val = reduction_identity<value_type>::band();
343 }
344
345 KOKKOS_INLINE_FUNCTION
346 value_type& reference() const { return *value.data(); }
347
348 KOKKOS_INLINE_FUNCTION
349 result_view_type view() const { return value; }
350
351 KOKKOS_INLINE_FUNCTION
352 bool references_scalar() const { return references_scalar_v; }
353};
354
355template <typename Scalar, typename... Properties>
357 ->BAnd<Scalar, typename View<Scalar, Properties...>::memory_space>;
358
359template <class Scalar, class Space>
360struct BOr {
361 public:
362 // Required
363 using reducer = BOr<Scalar, Space>;
364 using value_type = std::remove_cv_t<Scalar>;
365 static_assert(!std::is_pointer_v<value_type> && !std::is_array_v<value_type>);
366
367 using result_view_type = Kokkos::View<value_type, Space>;
368
369 private:
370 result_view_type value;
371 bool references_scalar_v;
372
373 public:
374 KOKKOS_INLINE_FUNCTION
375 BOr(value_type& value_) : value(&value_), references_scalar_v(true) {}
376
377 KOKKOS_INLINE_FUNCTION
378 BOr(const result_view_type& value_)
379 : value(value_), references_scalar_v(false) {}
380
381 // Required
382 KOKKOS_INLINE_FUNCTION
383 void join(value_type& dest, const value_type& src) const {
384 dest = dest | src;
385 }
386
387 KOKKOS_INLINE_FUNCTION
388 void init(value_type& val) const {
389 val = reduction_identity<value_type>::bor();
390 }
391
392 KOKKOS_INLINE_FUNCTION
393 value_type& reference() const { return *value.data(); }
394
395 KOKKOS_INLINE_FUNCTION
396 result_view_type view() const { return value; }
397
398 KOKKOS_INLINE_FUNCTION
399 bool references_scalar() const { return references_scalar_v; }
400};
401
402template <typename Scalar, typename... Properties>
404 ->BOr<Scalar, typename View<Scalar, Properties...>::memory_space>;
405
406template <class Scalar, class Index>
407struct ValLocScalar {
408 Scalar val;
409 Index loc;
410};
411
412template <class Scalar, class Index, class Space>
413struct MinLoc {
414 private:
415 using scalar_type = std::remove_cv_t<Scalar>;
416 using index_type = std::remove_cv_t<Index>;
417 static_assert(!std::is_pointer_v<scalar_type> &&
418 !std::is_array_v<scalar_type>);
419
420 public:
421 // Required
422 using reducer = MinLoc<Scalar, Index, Space>;
423 using value_type = ValLocScalar<scalar_type, index_type>;
424
425 using result_view_type = Kokkos::View<value_type, Space>;
426
427 private:
428 result_view_type value;
429 bool references_scalar_v;
430
431 public:
432 KOKKOS_INLINE_FUNCTION
433 MinLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
434
435 KOKKOS_INLINE_FUNCTION
436 MinLoc(const result_view_type& value_)
437 : value(value_), references_scalar_v(false) {}
438
439 // Required
440 KOKKOS_INLINE_FUNCTION
441 void join(value_type& dest, const value_type& src) const {
442 if (src.val < dest.val) dest = src;
443 }
444
445 KOKKOS_INLINE_FUNCTION
446 void init(value_type& val) const {
447 val.val = reduction_identity<scalar_type>::min();
448 val.loc = reduction_identity<index_type>::min();
449 }
450
451 KOKKOS_INLINE_FUNCTION
452 value_type& reference() const { return *value.data(); }
453
454 KOKKOS_INLINE_FUNCTION
455 result_view_type view() const { return value; }
456
457 KOKKOS_INLINE_FUNCTION
458 bool references_scalar() const { return references_scalar_v; }
459};
460
461template <typename Scalar, typename Index, typename... Properties>
462MinLoc(View<ValLocScalar<Scalar, Index>, Properties...> const&)
463 ->MinLoc<Scalar, Index,
465 Properties...>::memory_space>;
466
467template <class Scalar, class Index, class Space>
468struct MaxLoc {
469 private:
470 using scalar_type = std::remove_cv_t<Scalar>;
471 using index_type = std::remove_cv_t<Index>;
472 static_assert(!std::is_pointer_v<scalar_type> &&
473 !std::is_array_v<scalar_type>);
474
475 public:
476 // Required
477 using reducer = MaxLoc<Scalar, Index, Space>;
478 using value_type = ValLocScalar<scalar_type, index_type>;
479
480 using result_view_type = Kokkos::View<value_type, Space>;
481
482 private:
483 result_view_type value;
484 bool references_scalar_v;
485
486 public:
487 KOKKOS_INLINE_FUNCTION
488 MaxLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
489
490 KOKKOS_INLINE_FUNCTION
491 MaxLoc(const result_view_type& value_)
492 : value(value_), references_scalar_v(false) {}
493
494 // Required
495 KOKKOS_INLINE_FUNCTION
496 void join(value_type& dest, const value_type& src) const {
497 if (src.val > dest.val) dest = src;
498 }
499
500 KOKKOS_INLINE_FUNCTION
501 void init(value_type& val) const {
502 val.val = reduction_identity<scalar_type>::max();
503 val.loc = reduction_identity<index_type>::min();
504 }
505
506 KOKKOS_INLINE_FUNCTION
507 value_type& reference() const { return *value.data(); }
508
509 KOKKOS_INLINE_FUNCTION
510 result_view_type view() const { return value; }
511
512 KOKKOS_INLINE_FUNCTION
513 bool references_scalar() const { return references_scalar_v; }
514};
515
516template <typename Scalar, typename Index, typename... Properties>
517MaxLoc(View<ValLocScalar<Scalar, Index>, Properties...> const&)
518 ->MaxLoc<Scalar, Index,
520 Properties...>::memory_space>;
521
522template <class Scalar>
523struct MinMaxScalar {
524 Scalar min_val, max_val;
525};
526
527template <class Scalar, class Space>
528struct MinMax {
529 private:
530 using scalar_type = std::remove_cv_t<Scalar>;
531 static_assert(!std::is_pointer_v<scalar_type> &&
532 !std::is_array_v<scalar_type>);
533
534 public:
535 // Required
536 using reducer = MinMax<Scalar, Space>;
537 using value_type = MinMaxScalar<scalar_type>;
538
539 using result_view_type = Kokkos::View<value_type, Space>;
540
541 private:
542 result_view_type value;
543 bool references_scalar_v;
544
545 public:
546 KOKKOS_INLINE_FUNCTION
547 MinMax(value_type& value_) : value(&value_), references_scalar_v(true) {}
548
549 KOKKOS_INLINE_FUNCTION
550 MinMax(const result_view_type& value_)
551 : value(value_), references_scalar_v(false) {}
552
553 // Required
554 KOKKOS_INLINE_FUNCTION
555 void join(value_type& dest, const value_type& src) const {
556 if (src.min_val < dest.min_val) {
557 dest.min_val = src.min_val;
558 }
559 if (src.max_val > dest.max_val) {
560 dest.max_val = src.max_val;
561 }
562 }
563
564 KOKKOS_INLINE_FUNCTION
565 void init(value_type& val) const {
566 val.max_val = reduction_identity<scalar_type>::max();
567 val.min_val = reduction_identity<scalar_type>::min();
568 }
569
570 KOKKOS_INLINE_FUNCTION
571 value_type& reference() const { return *value.data(); }
572
573 KOKKOS_INLINE_FUNCTION
574 result_view_type view() const { return value; }
575
576 KOKKOS_INLINE_FUNCTION
577 bool references_scalar() const { return references_scalar_v; }
578};
579
580template <typename Scalar, typename... Properties>
581MinMax(View<MinMaxScalar<Scalar>, Properties...> const&)
582 ->MinMax<Scalar,
583 typename View<MinMaxScalar<Scalar>, Properties...>::memory_space>;
584
585template <class Scalar, class Index>
586struct MinMaxLocScalar {
587 Scalar min_val, max_val;
588 Index min_loc, max_loc;
589};
590
591template <class Scalar, class Index, class Space>
592struct MinMaxLoc {
593 private:
594 using scalar_type = std::remove_cv_t<Scalar>;
595 using index_type = std::remove_cv_t<Index>;
596 static_assert(!std::is_pointer_v<scalar_type> &&
597 !std::is_array_v<scalar_type>);
598
599 public:
600 // Required
601 using reducer = MinMaxLoc<Scalar, Index, Space>;
602 using value_type = MinMaxLocScalar<scalar_type, index_type>;
603
604 using result_view_type = Kokkos::View<value_type, Space>;
605
606 private:
607 result_view_type value;
608 bool references_scalar_v;
609
610 public:
611 KOKKOS_INLINE_FUNCTION
612 MinMaxLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
613
614 KOKKOS_INLINE_FUNCTION
615 MinMaxLoc(const result_view_type& value_)
616 : value(value_), references_scalar_v(false) {}
617
618 // Required
619 KOKKOS_INLINE_FUNCTION
620 void join(value_type& dest, const value_type& src) const {
621 if (src.min_val < dest.min_val) {
622 dest.min_val = src.min_val;
623 dest.min_loc = src.min_loc;
624 }
625 if (src.max_val > dest.max_val) {
626 dest.max_val = src.max_val;
627 dest.max_loc = src.max_loc;
628 }
629 }
630
631 KOKKOS_INLINE_FUNCTION
632 void init(value_type& val) const {
633 val.max_val = reduction_identity<scalar_type>::max();
634 val.min_val = reduction_identity<scalar_type>::min();
635 val.max_loc = reduction_identity<index_type>::min();
636 val.min_loc = reduction_identity<index_type>::min();
637 }
638
639 KOKKOS_INLINE_FUNCTION
640 value_type& reference() const { return *value.data(); }
641
642 KOKKOS_INLINE_FUNCTION
643 result_view_type view() const { return value; }
644
645 KOKKOS_INLINE_FUNCTION
646 bool references_scalar() const { return references_scalar_v; }
647};
648
649template <typename Scalar, typename Index, typename... Properties>
650MinMaxLoc(View<MinMaxLocScalar<Scalar, Index>, Properties...> const&)
651 ->MinMaxLoc<Scalar, Index,
653 Properties...>::memory_space>;
654
655// --------------------------------------------------
656// reducers added to support std algorithms
657// --------------------------------------------------
658
659//
660// MaxFirstLoc
661//
662template <class Scalar, class Index, class Space>
663struct MaxFirstLoc {
664 private:
665 using scalar_type = std::remove_cv_t<Scalar>;
666 using index_type = std::remove_cv_t<Index>;
667 static_assert(!std::is_pointer_v<scalar_type> &&
668 !std::is_array_v<scalar_type>);
669 static_assert(std::is_integral_v<index_type>);
670
671 public:
672 // Required
673 using reducer = MaxFirstLoc<Scalar, Index, Space>;
674 using value_type = ::Kokkos::ValLocScalar<scalar_type, index_type>;
675
676 using result_view_type = ::Kokkos::View<value_type, Space>;
677
678 private:
679 result_view_type value;
680 bool references_scalar_v;
681
682 public:
683 KOKKOS_INLINE_FUNCTION
684 MaxFirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
685
686 KOKKOS_INLINE_FUNCTION
687 MaxFirstLoc(const result_view_type& value_)
688 : value(value_), references_scalar_v(false) {}
689
690 // Required
691 KOKKOS_INLINE_FUNCTION
692 void join(value_type& dest, const value_type& src) const {
693 if (dest.val < src.val) {
694 dest = src;
695 } else if (!(src.val < dest.val)) {
696 dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc;
697 }
698 }
699
700 KOKKOS_INLINE_FUNCTION
701 void init(value_type& val) const {
702 val.val = reduction_identity<scalar_type>::max();
703 val.loc = reduction_identity<index_type>::min();
704 }
705
706 KOKKOS_INLINE_FUNCTION
707 value_type& reference() const { return *value.data(); }
708
709 KOKKOS_INLINE_FUNCTION
710 result_view_type view() const { return value; }
711
712 KOKKOS_INLINE_FUNCTION
713 bool references_scalar() const { return references_scalar_v; }
714};
715
716template <typename Scalar, typename Index, typename... Properties>
717MaxFirstLoc(View<ValLocScalar<Scalar, Index>, Properties...> const&)
718 ->MaxFirstLoc<Scalar, Index,
720 Properties...>::memory_space>;
721
722//
723// MaxFirstLocCustomComparator
724// recall that comp(a,b) returns true is a < b
725//
726template <class Scalar, class Index, class ComparatorType, class Space>
727struct MaxFirstLocCustomComparator {
728 private:
729 using scalar_type = std::remove_cv_t<Scalar>;
730 using index_type = std::remove_cv_t<Index>;
731 static_assert(!std::is_pointer_v<scalar_type> &&
732 !std::is_array_v<scalar_type>);
733 static_assert(std::is_integral_v<index_type>);
734
735 public:
736 // Required
737 using reducer =
738 MaxFirstLocCustomComparator<Scalar, Index, ComparatorType, Space>;
739 using value_type = ::Kokkos::ValLocScalar<scalar_type, index_type>;
740
741 using result_view_type = ::Kokkos::View<value_type, Space>;
742
743 private:
744 result_view_type value;
745 bool references_scalar_v;
746 ComparatorType m_comp;
747
748 public:
749 KOKKOS_INLINE_FUNCTION
750 MaxFirstLocCustomComparator(value_type& value_, ComparatorType comp_)
751 : value(&value_), references_scalar_v(true), m_comp(comp_) {}
752
753 KOKKOS_INLINE_FUNCTION
754 MaxFirstLocCustomComparator(const result_view_type& value_,
755 ComparatorType comp_)
756 : value(value_), references_scalar_v(false), m_comp(comp_) {}
757
758 // Required
759 KOKKOS_INLINE_FUNCTION
760 void join(value_type& dest, const value_type& src) const {
761 if (m_comp(dest.val, src.val)) {
762 dest = src;
763 } else if (!m_comp(src.val, dest.val)) {
764 dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc;
765 }
766 }
767
768 KOKKOS_INLINE_FUNCTION
769 void init(value_type& val) const {
770 val.val = reduction_identity<scalar_type>::max();
771 val.loc = reduction_identity<index_type>::min();
772 }
773
774 KOKKOS_INLINE_FUNCTION
775 value_type& reference() const { return *value.data(); }
776
777 KOKKOS_INLINE_FUNCTION
778 result_view_type view() const { return value; }
779
780 KOKKOS_INLINE_FUNCTION
781 bool references_scalar() const { return references_scalar_v; }
782};
783
784template <typename Scalar, typename Index, typename ComparatorType,
785 typename... Properties>
786MaxFirstLocCustomComparator(
787 View<ValLocScalar<Scalar, Index>, Properties...> const&, ComparatorType)
788 ->MaxFirstLocCustomComparator<Scalar, Index, ComparatorType,
790 Properties...>::memory_space>;
791
792//
793// MinFirstLoc
794//
795template <class Scalar, class Index, class Space>
796struct MinFirstLoc {
797 private:
798 using scalar_type = std::remove_cv_t<Scalar>;
799 using index_type = std::remove_cv_t<Index>;
800 static_assert(!std::is_pointer_v<scalar_type> &&
801 !std::is_array_v<scalar_type>);
802 static_assert(std::is_integral_v<index_type>);
803
804 public:
805 // Required
806 using reducer = MinFirstLoc<Scalar, Index, Space>;
807 using value_type = ::Kokkos::ValLocScalar<scalar_type, index_type>;
808
809 using result_view_type = ::Kokkos::View<value_type, Space>;
810
811 private:
812 result_view_type value;
813 bool references_scalar_v;
814
815 public:
816 KOKKOS_INLINE_FUNCTION
817 MinFirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
818
819 KOKKOS_INLINE_FUNCTION
820 MinFirstLoc(const result_view_type& value_)
821 : value(value_), references_scalar_v(false) {}
822
823 // Required
824 KOKKOS_INLINE_FUNCTION
825 void join(value_type& dest, const value_type& src) const {
826 if (src.val < dest.val) {
827 dest = src;
828 } else if (!(dest.val < src.val)) {
829 dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc;
830 }
831 }
832
833 KOKKOS_INLINE_FUNCTION
834 void init(value_type& val) const {
835 val.val = reduction_identity<scalar_type>::min();
836 val.loc = reduction_identity<index_type>::min();
837 }
838
839 KOKKOS_INLINE_FUNCTION
840 value_type& reference() const { return *value.data(); }
841
842 KOKKOS_INLINE_FUNCTION
843 result_view_type view() const { return value; }
844
845 KOKKOS_INLINE_FUNCTION
846 bool references_scalar() const { return references_scalar_v; }
847};
848
849template <typename Scalar, typename Index, typename... Properties>
850MinFirstLoc(View<ValLocScalar<Scalar, Index>, Properties...> const&)
851 ->MinFirstLoc<Scalar, Index,
853 Properties...>::memory_space>;
854
855//
856// MinFirstLocCustomComparator
857// recall that comp(a,b) returns true is a < b
858//
859template <class Scalar, class Index, class ComparatorType, class Space>
860struct MinFirstLocCustomComparator {
861 private:
862 using scalar_type = std::remove_cv_t<Scalar>;
863 using index_type = std::remove_cv_t<Index>;
864 static_assert(!std::is_pointer_v<scalar_type> &&
865 !std::is_array_v<scalar_type>);
866 static_assert(std::is_integral_v<index_type>);
867
868 public:
869 // Required
870 using reducer =
871 MinFirstLocCustomComparator<Scalar, Index, ComparatorType, Space>;
872 using value_type = ::Kokkos::ValLocScalar<scalar_type, index_type>;
873
874 using result_view_type = ::Kokkos::View<value_type, Space>;
875
876 private:
877 result_view_type value;
878 bool references_scalar_v;
879 ComparatorType m_comp;
880
881 public:
882 KOKKOS_INLINE_FUNCTION
883 MinFirstLocCustomComparator(value_type& value_, ComparatorType comp_)
884 : value(&value_), references_scalar_v(true), m_comp(comp_) {}
885
886 KOKKOS_INLINE_FUNCTION
887 MinFirstLocCustomComparator(const result_view_type& value_,
888 ComparatorType comp_)
889 : value(value_), references_scalar_v(false), m_comp(comp_) {}
890
891 // Required
892 KOKKOS_INLINE_FUNCTION
893 void join(value_type& dest, const value_type& src) const {
894 if (m_comp(src.val, dest.val)) {
895 dest = src;
896 } else if (!m_comp(dest.val, src.val)) {
897 dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc;
898 }
899 }
900
901 KOKKOS_INLINE_FUNCTION
902 void init(value_type& val) const {
903 val.val = reduction_identity<scalar_type>::min();
904 val.loc = reduction_identity<index_type>::min();
905 }
906
907 KOKKOS_INLINE_FUNCTION
908 value_type& reference() const { return *value.data(); }
909
910 KOKKOS_INLINE_FUNCTION
911 result_view_type view() const { return value; }
912
913 KOKKOS_INLINE_FUNCTION
914 bool references_scalar() const { return references_scalar_v; }
915};
916
917template <typename Scalar, typename Index, typename ComparatorType,
918 typename... Properties>
919MinFirstLocCustomComparator(
920 View<ValLocScalar<Scalar, Index>, Properties...> const&, ComparatorType)
921 ->MinFirstLocCustomComparator<Scalar, Index, ComparatorType,
923 Properties...>::memory_space>;
924
925//
926// MinMaxFirstLastLoc
927//
928template <class Scalar, class Index, class Space>
929struct MinMaxFirstLastLoc {
930 private:
931 using scalar_type = std::remove_cv_t<Scalar>;
932 using index_type = std::remove_cv_t<Index>;
933 static_assert(!std::is_pointer_v<scalar_type> &&
934 !std::is_array_v<scalar_type>);
935 static_assert(std::is_integral_v<index_type>);
936
937 public:
938 // Required
939 using reducer = MinMaxFirstLastLoc<Scalar, Index, Space>;
940 using value_type = ::Kokkos::MinMaxLocScalar<scalar_type, index_type>;
941
942 using result_view_type = ::Kokkos::View<value_type, Space>;
943
944 private:
945 result_view_type value;
946 bool references_scalar_v;
947
948 public:
949 KOKKOS_INLINE_FUNCTION
950 MinMaxFirstLastLoc(value_type& value_)
951 : value(&value_), references_scalar_v(true) {}
952
953 KOKKOS_INLINE_FUNCTION
954 MinMaxFirstLastLoc(const result_view_type& value_)
955 : value(value_), references_scalar_v(false) {}
956
957 // Required
958 KOKKOS_INLINE_FUNCTION
959 void join(value_type& dest, const value_type& src) const {
960 if (src.min_val < dest.min_val) {
961 dest.min_val = src.min_val;
962 dest.min_loc = src.min_loc;
963 } else if (!(dest.min_val < src.min_val)) {
964 dest.min_loc = (src.min_loc < dest.min_loc) ? src.min_loc : dest.min_loc;
965 }
966
967 if (dest.max_val < src.max_val) {
968 dest.max_val = src.max_val;
969 dest.max_loc = src.max_loc;
970 } else if (!(src.max_val < dest.max_val)) {
971 dest.max_loc = (src.max_loc > dest.max_loc) ? src.max_loc : dest.max_loc;
972 }
973 }
974
975 KOKKOS_INLINE_FUNCTION
976 void init(value_type& val) const {
977 val.max_val = ::Kokkos::reduction_identity<scalar_type>::max();
978 val.min_val = ::Kokkos::reduction_identity<scalar_type>::min();
979 val.max_loc = ::Kokkos::reduction_identity<index_type>::max();
980 val.min_loc = ::Kokkos::reduction_identity<index_type>::min();
981 }
982
983 KOKKOS_INLINE_FUNCTION
984 value_type& reference() const { return *value.data(); }
985
986 KOKKOS_INLINE_FUNCTION
987 result_view_type view() const { return value; }
988
989 KOKKOS_INLINE_FUNCTION
990 bool references_scalar() const { return references_scalar_v; }
991};
992
993template <typename Scalar, typename Index, typename... Properties>
994MinMaxFirstLastLoc(View<MinMaxLocScalar<Scalar, Index>, Properties...> const&)
995 ->MinMaxFirstLastLoc<Scalar, Index,
997 Properties...>::memory_space>;
998
999//
1000// MinMaxFirstLastLocCustomComparator
1001// recall that comp(a,b) returns true is a < b
1002//
1003template <class Scalar, class Index, class ComparatorType, class Space>
1004struct MinMaxFirstLastLocCustomComparator {
1005 private:
1006 using scalar_type = std::remove_cv_t<Scalar>;
1007 using index_type = std::remove_cv_t<Index>;
1008 static_assert(!std::is_pointer_v<scalar_type> &&
1009 !std::is_array_v<scalar_type>);
1010 static_assert(std::is_integral_v<index_type>);
1011
1012 public:
1013 // Required
1014 using reducer =
1015 MinMaxFirstLastLocCustomComparator<Scalar, Index, ComparatorType, Space>;
1016 using value_type = ::Kokkos::MinMaxLocScalar<scalar_type, index_type>;
1017
1018 using result_view_type = ::Kokkos::View<value_type, Space>;
1019
1020 private:
1021 result_view_type value;
1022 bool references_scalar_v;
1023 ComparatorType m_comp;
1024
1025 public:
1026 KOKKOS_INLINE_FUNCTION
1027 MinMaxFirstLastLocCustomComparator(value_type& value_, ComparatorType comp_)
1028 : value(&value_), references_scalar_v(true), m_comp(comp_) {}
1029
1030 KOKKOS_INLINE_FUNCTION
1031 MinMaxFirstLastLocCustomComparator(const result_view_type& value_,
1032 ComparatorType comp_)
1033 : value(value_), references_scalar_v(false), m_comp(comp_) {}
1034
1035 // Required
1036 KOKKOS_INLINE_FUNCTION
1037 void join(value_type& dest, const value_type& src) const {
1038 if (m_comp(src.min_val, dest.min_val)) {
1039 dest.min_val = src.min_val;
1040 dest.min_loc = src.min_loc;
1041 } else if (!m_comp(dest.min_val, src.min_val)) {
1042 dest.min_loc = (src.min_loc < dest.min_loc) ? src.min_loc : dest.min_loc;
1043 }
1044
1045 if (m_comp(dest.max_val, src.max_val)) {
1046 dest.max_val = src.max_val;
1047 dest.max_loc = src.max_loc;
1048 } else if (!m_comp(src.max_val, dest.max_val)) {
1049 dest.max_loc = (src.max_loc > dest.max_loc) ? src.max_loc : dest.max_loc;
1050 }
1051 }
1052
1053 KOKKOS_INLINE_FUNCTION
1054 void init(value_type& val) const {
1055 val.max_val = ::Kokkos::reduction_identity<scalar_type>::max();
1056 val.min_val = ::Kokkos::reduction_identity<scalar_type>::min();
1057 val.max_loc = ::Kokkos::reduction_identity<index_type>::max();
1058 val.min_loc = ::Kokkos::reduction_identity<index_type>::min();
1059 }
1060
1061 KOKKOS_INLINE_FUNCTION
1062 value_type& reference() const { return *value.data(); }
1063
1064 KOKKOS_INLINE_FUNCTION
1065 result_view_type view() const { return value; }
1066
1067 KOKKOS_INLINE_FUNCTION
1068 bool references_scalar() const { return references_scalar_v; }
1069};
1070
1071template <typename Scalar, typename Index, typename ComparatorType,
1072 typename... Properties>
1073MinMaxFirstLastLocCustomComparator(
1074 View<MinMaxLocScalar<Scalar, Index>, Properties...> const&, ComparatorType)
1075 ->MinMaxFirstLastLocCustomComparator<
1076 Scalar, Index, ComparatorType,
1078 Properties...>::memory_space>;
1079
1080//
1081// FirstLoc
1082//
1083template <class Index>
1084struct FirstLocScalar {
1085 Index min_loc_true;
1086};
1087
1088template <class Index, class Space>
1089struct FirstLoc {
1090 private:
1091 using index_type = std::remove_cv_t<Index>;
1092 static_assert(std::is_integral_v<index_type>);
1093
1094 public:
1095 // Required
1096 using reducer = FirstLoc<Index, Space>;
1097 using value_type = FirstLocScalar<index_type>;
1098
1099 using result_view_type = ::Kokkos::View<value_type, Space>;
1100
1101 private:
1102 result_view_type value;
1103 bool references_scalar_v;
1104
1105 public:
1106 KOKKOS_INLINE_FUNCTION
1107 FirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
1108
1109 KOKKOS_INLINE_FUNCTION
1110 FirstLoc(const result_view_type& value_)
1111 : value(value_), references_scalar_v(false) {}
1112
1113 // Required
1114 KOKKOS_INLINE_FUNCTION
1115 void join(value_type& dest, const value_type& src) const {
1116 dest.min_loc_true = (src.min_loc_true < dest.min_loc_true)
1117 ? src.min_loc_true
1118 : dest.min_loc_true;
1119 }
1120
1121 KOKKOS_INLINE_FUNCTION
1122 void init(value_type& val) const {
1123 val.min_loc_true = ::Kokkos::reduction_identity<index_type>::min();
1124 }
1125
1126 KOKKOS_INLINE_FUNCTION
1127 value_type& reference() const { return *value.data(); }
1128
1129 KOKKOS_INLINE_FUNCTION
1130 result_view_type view() const { return value; }
1131
1132 KOKKOS_INLINE_FUNCTION
1133 bool references_scalar() const { return references_scalar_v; }
1134};
1135
1136template <typename Index, typename... Properties>
1137FirstLoc(View<FirstLocScalar<Index>, Properties...> const&)
1138 ->FirstLoc<Index, typename View<FirstLocScalar<Index>,
1139 Properties...>::memory_space>;
1140
1141//
1142// LastLoc
1143//
1144template <class Index>
1145struct LastLocScalar {
1146 Index max_loc_true;
1147};
1148
1149template <class Index, class Space>
1150struct LastLoc {
1151 private:
1152 using index_type = std::remove_cv_t<Index>;
1153 static_assert(std::is_integral_v<index_type>);
1154
1155 public:
1156 // Required
1157 using reducer = LastLoc<Index, Space>;
1158 using value_type = LastLocScalar<index_type>;
1159
1160 using result_view_type = ::Kokkos::View<value_type, Space>;
1161
1162 private:
1163 result_view_type value;
1164 bool references_scalar_v;
1165
1166 public:
1167 KOKKOS_INLINE_FUNCTION
1168 LastLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
1169
1170 KOKKOS_INLINE_FUNCTION
1171 LastLoc(const result_view_type& value_)
1172 : value(value_), references_scalar_v(false) {}
1173
1174 // Required
1175 KOKKOS_INLINE_FUNCTION
1176 void join(value_type& dest, const value_type& src) const {
1177 dest.max_loc_true = (src.max_loc_true > dest.max_loc_true)
1178 ? src.max_loc_true
1179 : dest.max_loc_true;
1180 }
1181
1182 KOKKOS_INLINE_FUNCTION
1183 void init(value_type& val) const {
1184 val.max_loc_true = ::Kokkos::reduction_identity<index_type>::max();
1185 }
1186
1187 KOKKOS_INLINE_FUNCTION
1188 value_type& reference() const { return *value.data(); }
1189
1190 KOKKOS_INLINE_FUNCTION
1191 result_view_type view() const { return value; }
1192
1193 KOKKOS_INLINE_FUNCTION
1194 bool references_scalar() const { return references_scalar_v; }
1195};
1196
1197template <typename Index, typename... Properties>
1198LastLoc(View<LastLocScalar<Index>, Properties...> const&)
1199 ->LastLoc<Index,
1200 typename View<LastLocScalar<Index>, Properties...>::memory_space>;
1201
1202template <class Index>
1203struct StdIsPartScalar {
1204 Index max_loc_true, min_loc_false;
1205};
1206
1207//
1208// StdIsPartitioned
1209//
1210template <class Index, class Space>
1211struct StdIsPartitioned {
1212 private:
1213 using index_type = std::remove_cv_t<Index>;
1214 static_assert(std::is_integral_v<index_type>);
1215
1216 public:
1217 // Required
1218 using reducer = StdIsPartitioned<Index, Space>;
1219 using value_type = StdIsPartScalar<index_type>;
1220
1221 using result_view_type = ::Kokkos::View<value_type, Space>;
1222
1223 private:
1224 result_view_type value;
1225 bool references_scalar_v;
1226
1227 public:
1228 KOKKOS_INLINE_FUNCTION
1229 StdIsPartitioned(value_type& value_)
1230 : value(&value_), references_scalar_v(true) {}
1231
1232 KOKKOS_INLINE_FUNCTION
1233 StdIsPartitioned(const result_view_type& value_)
1234 : value(value_), references_scalar_v(false) {}
1235
1236 // Required
1237 KOKKOS_INLINE_FUNCTION
1238 void join(value_type& dest, const value_type& src) const {
1239 dest.max_loc_true = (dest.max_loc_true < src.max_loc_true)
1240 ? src.max_loc_true
1241 : dest.max_loc_true;
1242
1243 dest.min_loc_false = (dest.min_loc_false < src.min_loc_false)
1244 ? dest.min_loc_false
1245 : src.min_loc_false;
1246 }
1247
1248 KOKKOS_INLINE_FUNCTION
1249 void init(value_type& val) const {
1250 val.max_loc_true = ::Kokkos::reduction_identity<index_type>::max();
1251 val.min_loc_false = ::Kokkos::reduction_identity<index_type>::min();
1252 }
1253
1254 KOKKOS_INLINE_FUNCTION
1255 value_type& reference() const { return *value.data(); }
1256
1257 KOKKOS_INLINE_FUNCTION
1258 result_view_type view() const { return value; }
1259
1260 KOKKOS_INLINE_FUNCTION
1261 bool references_scalar() const { return references_scalar_v; }
1262};
1263
1264template <typename Index, typename... Properties>
1265StdIsPartitioned(View<StdIsPartScalar<Index>, Properties...> const&)
1266 ->StdIsPartitioned<Index, typename View<StdIsPartScalar<Index>,
1267 Properties...>::memory_space>;
1268
1269template <class Index>
1270struct StdPartPointScalar {
1271 Index min_loc_false;
1272};
1273
1274//
1275// StdPartitionPoint
1276//
1277template <class Index, class Space>
1278struct StdPartitionPoint {
1279 private:
1280 using index_type = std::remove_cv_t<Index>;
1281 static_assert(std::is_integral_v<index_type>);
1282
1283 public:
1284 // Required
1285 using reducer = StdPartitionPoint<Index, Space>;
1286 using value_type = StdPartPointScalar<index_type>;
1287
1288 using result_view_type = ::Kokkos::View<value_type, Space>;
1289
1290 private:
1291 result_view_type value;
1292 bool references_scalar_v;
1293
1294 public:
1295 KOKKOS_INLINE_FUNCTION
1296 StdPartitionPoint(value_type& value_)
1297 : value(&value_), references_scalar_v(true) {}
1298
1299 KOKKOS_INLINE_FUNCTION
1300 StdPartitionPoint(const result_view_type& value_)
1301 : value(value_), references_scalar_v(false) {}
1302
1303 // Required
1304 KOKKOS_INLINE_FUNCTION
1305 void join(value_type& dest, const value_type& src) const {
1306 dest.min_loc_false = (dest.min_loc_false < src.min_loc_false)
1307 ? dest.min_loc_false
1308 : src.min_loc_false;
1309 }
1310
1311 KOKKOS_INLINE_FUNCTION
1312 void init(value_type& val) const {
1313 val.min_loc_false = ::Kokkos::reduction_identity<index_type>::min();
1314 }
1315
1316 KOKKOS_INLINE_FUNCTION
1317 value_type& reference() const { return *value.data(); }
1318
1319 KOKKOS_INLINE_FUNCTION
1320 result_view_type view() const { return value; }
1321
1322 KOKKOS_INLINE_FUNCTION
1323 bool references_scalar() const { return references_scalar_v; }
1324};
1325
1326template <typename Index, typename... Properties>
1327StdPartitionPoint(View<StdPartPointScalar<Index>, Properties...> const&)
1328 ->StdPartitionPoint<Index, typename View<StdPartPointScalar<Index>,
1329 Properties...>::memory_space>;
1330
1331} // namespace Kokkos
1332namespace Kokkos {
1333namespace Impl {
1334
1335template <typename FunctorType, typename FunctorAnalysisReducerType,
1336 typename Enable>
1337class CombinedFunctorReducer {
1338 public:
1339 using functor_type = FunctorType;
1340 using reducer_type = FunctorAnalysisReducerType;
1341 CombinedFunctorReducer(const FunctorType& functor,
1342 const FunctorAnalysisReducerType& reducer)
1343 : m_functor(functor), m_reducer(reducer) {}
1344 KOKKOS_FUNCTION const FunctorType& get_functor() const { return m_functor; }
1345 KOKKOS_FUNCTION const FunctorAnalysisReducerType& get_reducer() const {
1346 return m_reducer;
1347 }
1348
1349 private:
1350 FunctorType m_functor;
1351 FunctorAnalysisReducerType m_reducer;
1352};
1353template <typename FunctorType, typename FunctorAnalysisReducerType>
1354class CombinedFunctorReducer<
1355 FunctorType, FunctorAnalysisReducerType,
1356 std::enable_if_t<std::is_same_v<
1357 FunctorType, typename FunctorAnalysisReducerType::functor_type>>> {
1358 public:
1359 using functor_type = FunctorType;
1360 using reducer_type = FunctorAnalysisReducerType;
1361 CombinedFunctorReducer(const FunctorType& functor,
1362 const FunctorAnalysisReducerType&)
1363 : m_reducer(functor) {}
1364 KOKKOS_FUNCTION const FunctorType& get_functor() const {
1365 return m_reducer.get_functor();
1366 }
1367 KOKKOS_FUNCTION const FunctorAnalysisReducerType& get_reducer() const {
1368 return m_reducer;
1369 }
1370
1371 private:
1372 FunctorAnalysisReducerType m_reducer;
1373};
1374
1375template <class T, class ReturnType, class ValueTraits>
1376struct ParallelReduceReturnValue;
1377
1378template <class ReturnType, class FunctorType>
1379struct ParallelReduceReturnValue<
1380 std::enable_if_t<Kokkos::is_view<ReturnType>::value>, ReturnType,
1381 FunctorType> {
1382 using return_type = ReturnType;
1383 using reducer_type = InvalidType;
1384
1385 using value_type_scalar = typename return_type::value_type;
1386 using value_type_array = typename return_type::value_type* const;
1387
1388 using value_type = std::conditional_t<return_type::rank == 0,
1389 value_type_scalar, value_type_array>;
1390
1391 static return_type& return_value(ReturnType& return_val, const FunctorType&) {
1392 return return_val;
1393 }
1394};
1395
1396template <class ReturnType, class FunctorType>
1397struct ParallelReduceReturnValue<
1398 std::enable_if_t<!Kokkos::is_view<ReturnType>::value &&
1399 (!std::is_array<ReturnType>::value &&
1400 !std::is_pointer<ReturnType>::value) &&
1401 !Kokkos::is_reducer<ReturnType>::value>,
1402 ReturnType, FunctorType> {
1403 using return_type =
1404 Kokkos::View<ReturnType, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>;
1405
1406 using reducer_type = InvalidType;
1407
1408 using value_type = typename return_type::value_type;
1409
1410 static return_type return_value(ReturnType& return_val, const FunctorType&) {
1411 return return_type(&return_val);
1412 }
1413};
1414
1415template <class ReturnType, class FunctorType>
1416struct ParallelReduceReturnValue<
1417 std::enable_if_t<(std::is_array<ReturnType>::value ||
1418 std::is_pointer<ReturnType>::value)>,
1419 ReturnType, FunctorType> {
1420 using return_type = Kokkos::View<std::remove_const_t<ReturnType>,
1421 Kokkos::HostSpace, Kokkos::MemoryUnmanaged>;
1422
1423 using reducer_type = InvalidType;
1424
1425 using value_type = typename return_type::value_type[];
1426
1427 static return_type return_value(ReturnType& return_val,
1428 const FunctorType& functor) {
1429 if (std::is_array<ReturnType>::value)
1430 return return_type(return_val);
1431 else
1432 return return_type(return_val, functor.value_count);
1433 }
1434};
1435
1436template <class ReturnType, class FunctorType>
1437struct ParallelReduceReturnValue<
1438 std::enable_if_t<Kokkos::is_reducer<ReturnType>::value>, ReturnType,
1439 FunctorType> {
1440 using return_type = typename ReturnType::result_view_type;
1441 using reducer_type = ReturnType;
1442 using value_type = typename return_type::value_type;
1443
1444 static auto return_value(ReturnType& return_val, const FunctorType&) {
1445 return return_val.view();
1446 }
1447};
1448
1449template <class T, class ReturnType, class FunctorType>
1450struct ParallelReducePolicyType;
1451
1452template <class PolicyType, class FunctorType>
1453struct ParallelReducePolicyType<
1454 std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value>,
1455 PolicyType, FunctorType> {
1456 using policy_type = PolicyType;
1457 static PolicyType policy(const PolicyType& policy_) { return policy_; }
1458};
1459
1460template <class PolicyType, class FunctorType>
1461struct ParallelReducePolicyType<
1462 std::enable_if_t<std::is_integral<PolicyType>::value>, PolicyType,
1463 FunctorType> {
1464 using execution_space =
1465 typename Impl::FunctorPolicyExecutionSpace<FunctorType,
1466 void>::execution_space;
1467
1468 using policy_type = Kokkos::RangePolicy<execution_space>;
1469
1470 static policy_type policy(const PolicyType& policy_) {
1471 return policy_type(0, policy_);
1472 }
1473};
1474
1475template <class FunctorType, class ExecPolicy, class ValueType,
1476 class ExecutionSpace>
1477struct ParallelReduceFunctorType {
1478 using functor_type = FunctorType;
1479 static const functor_type& functor(const functor_type& functor) {
1480 return functor;
1481 }
1482};
1483
1484template <class PolicyType, class FunctorType, class ReturnType>
1485struct ParallelReduceAdaptor {
1486 using return_value_adapter =
1487 Impl::ParallelReduceReturnValue<void, ReturnType, FunctorType>;
1488
1489 static inline void execute_impl(const std::string& label,
1490 const PolicyType& policy,
1491 const FunctorType& functor,
1492 ReturnType& return_value) {
1493 using PassedReducerType = typename return_value_adapter::reducer_type;
1494 uint64_t kpID = 0;
1495
1496 PolicyType inner_policy = policy;
1497 Kokkos::Tools::Impl::begin_parallel_reduce<PassedReducerType>(
1498 inner_policy, functor, label, kpID);
1499
1500 using ReducerSelector =
1501 Kokkos::Impl::if_c<std::is_same<InvalidType, PassedReducerType>::value,
1502 FunctorType, PassedReducerType>;
1503 using Analysis = FunctorAnalysis<FunctorPatternInterface::REDUCE,
1504 PolicyType, typename ReducerSelector::type,
1505 typename return_value_adapter::value_type>;
1506 Kokkos::Impl::shared_allocation_tracking_disable();
1507 CombinedFunctorReducer functor_reducer(
1508 functor, typename Analysis::Reducer(
1509 ReducerSelector::select(functor, return_value)));
1510
1511 // FIXME Remove "Wrapper" once all backends implement the new interface
1512 Impl::ParallelReduce<decltype(functor_reducer), PolicyType,
1513 typename Impl::FunctorPolicyExecutionSpace<
1514 FunctorType, PolicyType>::execution_space>
1515 closure(functor_reducer, inner_policy,
1516 return_value_adapter::return_value(return_value, functor));
1517 Kokkos::Impl::shared_allocation_tracking_enable();
1518 closure.execute();
1519
1520 Kokkos::Tools::Impl::end_parallel_reduce<PassedReducerType>(
1521 inner_policy, functor, label, kpID);
1522 }
1523
1524 static constexpr bool is_array_reduction =
1525 Impl::FunctorAnalysis<
1526 Impl::FunctorPatternInterface::REDUCE, PolicyType, FunctorType,
1527 typename return_value_adapter::value_type>::StaticValueSize == 0;
1528
1529 template <typename Dummy = ReturnType>
1530 static inline std::enable_if_t<!(is_array_reduction &&
1531 std::is_pointer<Dummy>::value)>
1532 execute(const std::string& label, const PolicyType& policy,
1533 const FunctorType& functor, ReturnType& return_value) {
1534 execute_impl(label, policy, functor, return_value);
1535 }
1536};
1537} // namespace Impl
1538
1539//----------------------------------------------------------------------------
1540
1551
1552// Parallel Reduce Blocking behavior
1553
1554namespace Impl {
1555template <typename T>
1556struct ReducerHasTestReferenceFunction {
1557 template <typename E>
1558 static std::true_type test_func(decltype(&E::references_scalar));
1559 template <typename E>
1560 static std::false_type test_func(...);
1561
1562 enum {
1563 value = std::is_same<std::true_type, decltype(test_func<T>(nullptr))>::value
1564 };
1565};
1566
1567template <class ExecutionSpace, class Arg>
1568constexpr std::enable_if_t<
1569 // constraints only necessary because SFINAE lacks subsumption
1570 !ReducerHasTestReferenceFunction<Arg>::value &&
1571 !Kokkos::is_view<Arg>::value,
1572 // return type:
1573 bool>
1574parallel_reduce_needs_fence(ExecutionSpace const&, Arg const&) {
1575 return true;
1576}
1577
1578template <class ExecutionSpace, class Reducer>
1579constexpr std::enable_if_t<
1580 // equivalent to:
1581 // (requires (Reducer const& r) {
1582 // { reducer.references_scalar() } -> std::convertible_to<bool>;
1583 // })
1584 ReducerHasTestReferenceFunction<Reducer>::value,
1585 // return type:
1586 bool>
1587parallel_reduce_needs_fence(ExecutionSpace const&, Reducer const& reducer) {
1588 return reducer.references_scalar();
1589}
1590
1591template <class ExecutionSpace, class ViewLike>
1592constexpr std::enable_if_t<
1593 // requires Kokkos::ViewLike<ViewLike>
1594 Kokkos::is_view<ViewLike>::value,
1595 // return type:
1596 bool>
1597parallel_reduce_needs_fence(ExecutionSpace const&, ViewLike const&) {
1598 return false;
1599}
1600
1601template <class ExecutionSpace, class... Args>
1602struct ParallelReduceFence {
1603 template <class... ArgsDeduced>
1604 static void fence(const ExecutionSpace& ex, const std::string& name,
1605 ArgsDeduced&&... args) {
1606 if (Impl::parallel_reduce_needs_fence(ex, (ArgsDeduced &&) args...)) {
1607 ex.fence(name);
1608 }
1609 }
1610};
1611
1612} // namespace Impl
1613
1651
1652// ReturnValue is scalar or array: take by reference
1653
1654template <class PolicyType, class FunctorType, class ReturnType>
1655inline std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value &&
1656 !(Kokkos::is_view<ReturnType>::value ||
1657 Kokkos::is_reducer<ReturnType>::value ||
1658 std::is_pointer<ReturnType>::value)>
1659parallel_reduce(const std::string& label, const PolicyType& policy,
1660 const FunctorType& functor, ReturnType& return_value) {
1661 static_assert(
1662 !std::is_const<ReturnType>::value,
1663 "A const reduction result type is only allowed for a View, pointer or "
1664 "reducer return type!");
1665
1666 Impl::ParallelReduceAdaptor<PolicyType, FunctorType, ReturnType>::execute(
1667 label, policy, functor, return_value);
1668 Impl::ParallelReduceFence<typename PolicyType::execution_space, ReturnType>::
1669 fence(
1670 policy.space(),
1671 "Kokkos::parallel_reduce: fence due to result being value, not view",
1672 return_value);
1673}
1674
1675template <class PolicyType, class FunctorType, class ReturnType>
1676inline std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value &&
1677 !(Kokkos::is_view<ReturnType>::value ||
1678 Kokkos::is_reducer<ReturnType>::value ||
1679 std::is_pointer<ReturnType>::value)>
1680parallel_reduce(const PolicyType& policy, const FunctorType& functor,
1681 ReturnType& return_value) {
1682 static_assert(
1683 !std::is_const<ReturnType>::value,
1684 "A const reduction result type is only allowed for a View, pointer or "
1685 "reducer return type!");
1686
1687 Impl::ParallelReduceAdaptor<PolicyType, FunctorType, ReturnType>::execute(
1688 "", policy, functor, return_value);
1689 Impl::ParallelReduceFence<typename PolicyType::execution_space, ReturnType>::
1690 fence(
1691 policy.space(),
1692 "Kokkos::parallel_reduce: fence due to result being value, not view",
1693 return_value);
1694}
1695
1696template <class FunctorType, class ReturnType>
1697inline std::enable_if_t<!(Kokkos::is_view<ReturnType>::value ||
1698 Kokkos::is_reducer<ReturnType>::value ||
1699 std::is_pointer<ReturnType>::value)>
1700parallel_reduce(const size_t& policy, const FunctorType& functor,
1701 ReturnType& return_value) {
1702 static_assert(
1703 !std::is_const<ReturnType>::value,
1704 "A const reduction result type is only allowed for a View, pointer or "
1705 "reducer return type!");
1706
1707 using policy_type =
1708 typename Impl::ParallelReducePolicyType<void, size_t,
1709 FunctorType>::policy_type;
1710
1711 Impl::ParallelReduceAdaptor<policy_type, FunctorType, ReturnType>::execute(
1712 "", policy_type(0, policy), functor, return_value);
1713 Impl::ParallelReduceFence<typename policy_type::execution_space, ReturnType>::
1714 fence(
1715 typename policy_type::execution_space(),
1716 "Kokkos::parallel_reduce: fence due to result being value, not view",
1717 return_value);
1718}
1719
1720template <class FunctorType, class ReturnType>
1721inline std::enable_if_t<!(Kokkos::is_view<ReturnType>::value ||
1722 Kokkos::is_reducer<ReturnType>::value ||
1723 std::is_pointer<ReturnType>::value)>
1724parallel_reduce(const std::string& label, const size_t& policy,
1725 const FunctorType& functor, ReturnType& return_value) {
1726 static_assert(
1727 !std::is_const<ReturnType>::value,
1728 "A const reduction result type is only allowed for a View, pointer or "
1729 "reducer return type!");
1730
1731 using policy_type =
1732 typename Impl::ParallelReducePolicyType<void, size_t,
1733 FunctorType>::policy_type;
1734 Impl::ParallelReduceAdaptor<policy_type, FunctorType, ReturnType>::execute(
1735 label, policy_type(0, policy), functor, return_value);
1736 Impl::ParallelReduceFence<typename policy_type::execution_space, ReturnType>::
1737 fence(
1738 typename policy_type::execution_space(),
1739 "Kokkos::parallel_reduce: fence due to result being value, not view",
1740 return_value);
1741}
1742
1743// ReturnValue as View or Reducer: take by copy to allow for inline construction
1744
1745template <class PolicyType, class FunctorType, class ReturnType>
1746inline std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value &&
1747 (Kokkos::is_view<ReturnType>::value ||
1748 Kokkos::is_reducer<ReturnType>::value ||
1749 std::is_pointer<ReturnType>::value)>
1750parallel_reduce(const std::string& label, const PolicyType& policy,
1751 const FunctorType& functor, const ReturnType& return_value) {
1752 ReturnType return_value_impl = return_value;
1753 Impl::ParallelReduceAdaptor<PolicyType, FunctorType, ReturnType>::execute(
1754 label, policy, functor, return_value_impl);
1755 Impl::ParallelReduceFence<typename PolicyType::execution_space, ReturnType>::
1756 fence(
1757 policy.space(),
1758 "Kokkos::parallel_reduce: fence due to result being value, not view",
1759 return_value);
1760}
1761
1762template <class PolicyType, class FunctorType, class ReturnType>
1763inline std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value &&
1764 (Kokkos::is_view<ReturnType>::value ||
1765 Kokkos::is_reducer<ReturnType>::value ||
1766 std::is_pointer<ReturnType>::value)>
1767parallel_reduce(const PolicyType& policy, const FunctorType& functor,
1768 const ReturnType& return_value) {
1769 ReturnType return_value_impl = return_value;
1770 Impl::ParallelReduceAdaptor<PolicyType, FunctorType, ReturnType>::execute(
1771 "", policy, functor, return_value_impl);
1772 Impl::ParallelReduceFence<typename PolicyType::execution_space, ReturnType>::
1773 fence(
1774 policy.space(),
1775 "Kokkos::parallel_reduce: fence due to result being value, not view",
1776 return_value);
1777}
1778
1779template <class FunctorType, class ReturnType>
1780inline std::enable_if_t<Kokkos::is_view<ReturnType>::value ||
1781 Kokkos::is_reducer<ReturnType>::value ||
1782 std::is_pointer<ReturnType>::value>
1783parallel_reduce(const size_t& policy, const FunctorType& functor,
1784 const ReturnType& return_value) {
1785 using policy_type =
1786 typename Impl::ParallelReducePolicyType<void, size_t,
1787 FunctorType>::policy_type;
1788 ReturnType return_value_impl = return_value;
1789 Impl::ParallelReduceAdaptor<policy_type, FunctorType, ReturnType>::execute(
1790 "", policy_type(0, policy), functor, return_value_impl);
1791 Impl::ParallelReduceFence<typename policy_type::execution_space, ReturnType>::
1792 fence(
1793 typename policy_type::execution_space(),
1794 "Kokkos::parallel_reduce: fence due to result being value, not view",
1795 return_value);
1796}
1797
1798template <class FunctorType, class ReturnType>
1799inline std::enable_if_t<Kokkos::is_view<ReturnType>::value ||
1800 Kokkos::is_reducer<ReturnType>::value ||
1801 std::is_pointer<ReturnType>::value>
1802parallel_reduce(const std::string& label, const size_t& policy,
1803 const FunctorType& functor, const ReturnType& return_value) {
1804 using policy_type =
1805 typename Impl::ParallelReducePolicyType<void, size_t,
1806 FunctorType>::policy_type;
1807 ReturnType return_value_impl = return_value;
1808 Impl::ParallelReduceAdaptor<policy_type, FunctorType, ReturnType>::execute(
1809 label, policy_type(0, policy), functor, return_value_impl);
1810 Impl::ParallelReduceFence<typename policy_type::execution_space, ReturnType>::
1811 fence(
1812 typename policy_type::execution_space(),
1813 "Kokkos::parallel_reduce: fence due to result being value, not view",
1814 return_value);
1815}
1816
1817// No Return Argument
1818
1819template <class PolicyType, class FunctorType>
1820inline void parallel_reduce(
1821 const std::string& label, const PolicyType& policy,
1822 const FunctorType& functor,
1823 std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value>* =
1824 nullptr) {
1825 using FunctorAnalysis =
1826 Impl::FunctorAnalysis<Impl::FunctorPatternInterface::REDUCE, PolicyType,
1827 FunctorType, void>;
1828 using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0),
1829 typename FunctorAnalysis::value_type,
1830 typename FunctorAnalysis::pointer_type>;
1831
1832 static_assert(
1833 FunctorAnalysis::has_final_member_function,
1834 "Calling parallel_reduce without either return value or final function.");
1835
1836 using result_view_type =
1837 Kokkos::View<value_type, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>;
1838 result_view_type result_view;
1839
1840 Impl::ParallelReduceAdaptor<PolicyType, FunctorType,
1841 result_view_type>::execute(label, policy, functor,
1842 result_view);
1843}
1844
1845template <class PolicyType, class FunctorType>
1846inline void parallel_reduce(
1847 const PolicyType& policy, const FunctorType& functor,
1848 std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value>* =
1849 nullptr) {
1850 using FunctorAnalysis =
1851 Impl::FunctorAnalysis<Impl::FunctorPatternInterface::REDUCE, PolicyType,
1852 FunctorType, void>;
1853 using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0),
1854 typename FunctorAnalysis::value_type,
1855 typename FunctorAnalysis::pointer_type>;
1856
1857 static_assert(
1858 FunctorAnalysis::has_final_member_function,
1859 "Calling parallel_reduce without either return value or final function.");
1860
1861 using result_view_type =
1862 Kokkos::View<value_type, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>;
1863 result_view_type result_view;
1864
1865 Impl::ParallelReduceAdaptor<PolicyType, FunctorType,
1866 result_view_type>::execute("", policy, functor,
1867 result_view);
1868}
1869
1870template <class FunctorType>
1871inline void parallel_reduce(const size_t& policy, const FunctorType& functor) {
1872 using policy_type =
1873 typename Impl::ParallelReducePolicyType<void, size_t,
1874 FunctorType>::policy_type;
1875 using FunctorAnalysis =
1876 Impl::FunctorAnalysis<Impl::FunctorPatternInterface::REDUCE, policy_type,
1877 FunctorType, void>;
1878 using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0),
1879 typename FunctorAnalysis::value_type,
1880 typename FunctorAnalysis::pointer_type>;
1881
1882 static_assert(
1883 FunctorAnalysis::has_final_member_function,
1884 "Calling parallel_reduce without either return value or final function.");
1885
1886 using result_view_type =
1887 Kokkos::View<value_type, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>;
1888 result_view_type result_view;
1889
1890 Impl::ParallelReduceAdaptor<policy_type, FunctorType,
1891 result_view_type>::execute("",
1892 policy_type(0, policy),
1893 functor, result_view);
1894}
1895
1896template <class FunctorType>
1897inline void parallel_reduce(const std::string& label, const size_t& policy,
1898 const FunctorType& functor) {
1899 using policy_type =
1900 typename Impl::ParallelReducePolicyType<void, size_t,
1901 FunctorType>::policy_type;
1902 using FunctorAnalysis =
1903 Impl::FunctorAnalysis<Impl::FunctorPatternInterface::REDUCE, policy_type,
1904 FunctorType, void>;
1905 using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0),
1906 typename FunctorAnalysis::value_type,
1907 typename FunctorAnalysis::pointer_type>;
1908
1909 static_assert(
1910 FunctorAnalysis::has_final_member_function,
1911 "Calling parallel_reduce without either return value or final function.");
1912
1913 using result_view_type =
1914 Kokkos::View<value_type, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>;
1915 result_view_type result_view;
1916
1917 Impl::ParallelReduceAdaptor<policy_type, FunctorType,
1918 result_view_type>::execute(label,
1919 policy_type(0, policy),
1920 functor, result_view);
1921}
1922
1923} // namespace Kokkos
1924
1925#endif // KOKKOS_PARALLEL_REDUCE_HPP
View to an array of data.
ReturnType
ScopeGuard Some user scope issues have been identified with some Kokkos::finalize calls; ScopeGuard a...