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