Kokkos Core Kernels Package  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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>
19 static_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 
31 namespace Kokkos {
32 
33 template <class Scalar, class Space>
34 struct 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 
74 template <typename Scalar, typename... Properties>
75 KOKKOS_DEDUCTION_GUIDE Sum(View<Scalar, Properties...> const&)
76  ->Sum<Scalar, typename View<Scalar, Properties...>::memory_space>;
77 
78 template <class Scalar, class Space>
79 struct 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 
119 template <typename Scalar, typename... Properties>
120 KOKKOS_DEDUCTION_GUIDE Prod(View<Scalar, Properties...> const&)
121  ->Prod<Scalar, typename View<Scalar, Properties...>::memory_space>;
122 
123 template <class Scalar, class Space>
124 struct 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 
166 template <typename Scalar, typename... Properties>
167 KOKKOS_DEDUCTION_GUIDE Min(View<Scalar, Properties...> const&)
168  ->Min<Scalar, typename View<Scalar, Properties...>::memory_space>;
169 
170 template <class Scalar, class Space>
171 struct 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 
214 template <typename Scalar, typename... Properties>
215 KOKKOS_DEDUCTION_GUIDE Max(View<Scalar, Properties...> const&)
216  ->Max<Scalar, typename View<Scalar, Properties...>::memory_space>;
217 
218 template <class Scalar, class Space>
219 struct 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 
260 template <typename Scalar, typename... Properties>
261 KOKKOS_DEDUCTION_GUIDE LAnd(View<Scalar, Properties...> const&)
262  ->LAnd<Scalar, typename View<Scalar, Properties...>::memory_space>;
263 
264 template <class Scalar, class Space>
265 struct 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 
307 template <typename Scalar, typename... Properties>
308 KOKKOS_DEDUCTION_GUIDE LOr(View<Scalar, Properties...> const&)
309  ->LOr<Scalar, typename View<Scalar, Properties...>::memory_space>;
310 
311 template <class Scalar, class Space>
312 struct 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 
354 template <typename Scalar, typename... Properties>
355 KOKKOS_DEDUCTION_GUIDE BAnd(View<Scalar, Properties...> const&)
356  ->BAnd<Scalar, typename View<Scalar, Properties...>::memory_space>;
357 
358 template <class Scalar, class Space>
359 struct 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 
401 template <typename Scalar, typename... Properties>
402 KOKKOS_DEDUCTION_GUIDE BOr(View<Scalar, Properties...> const&)
403  ->BOr<Scalar, typename View<Scalar, Properties...>::memory_space>;
404 
405 template <class Scalar, class Index>
406 struct ValLocScalar {
407  Scalar val;
408  Index loc;
409 };
410 
411 template <class Scalar, class Index, class Space>
412 struct 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) dest = src;
442  }
443 
444  KOKKOS_INLINE_FUNCTION
445  void init(value_type& val) const {
446  val.val = reduction_identity<scalar_type>::min();
447  val.loc = reduction_identity<index_type>::min();
448  }
449 
450  KOKKOS_INLINE_FUNCTION
451  value_type& reference() const { return *value.data(); }
452 
453  KOKKOS_INLINE_FUNCTION
454  result_view_type view() const { return value; }
455 
456  KOKKOS_INLINE_FUNCTION
457  bool references_scalar() const { return references_scalar_v; }
458 };
459 
460 template <typename Scalar, typename Index, typename... Properties>
461 KOKKOS_DEDUCTION_GUIDE MinLoc(
462  View<ValLocScalar<Scalar, Index>, Properties...> const&)
463  ->MinLoc<Scalar, Index,
464  typename View<ValLocScalar<Scalar, Index>,
465  Properties...>::memory_space>;
466 
467 template <class Scalar, class Index, class Space>
468 struct 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 
516 template <typename Scalar, typename Index, typename... Properties>
517 KOKKOS_DEDUCTION_GUIDE MaxLoc(
518  View<ValLocScalar<Scalar, Index>, Properties...> const&)
519  ->MaxLoc<Scalar, Index,
520  typename View<ValLocScalar<Scalar, Index>,
521  Properties...>::memory_space>;
522 
523 template <class Scalar>
524 struct MinMaxScalar {
525  Scalar min_val, max_val;
526 };
527 
528 template <class Scalar, class Space>
529 struct MinMax {
530  private:
531  using scalar_type = std::remove_cv_t<Scalar>;
532  static_assert(!std::is_pointer_v<scalar_type> &&
533  !std::is_array_v<scalar_type>);
534 
535  public:
536  // Required
537  using reducer = MinMax<Scalar, Space>;
538  using value_type = MinMaxScalar<scalar_type>;
539 
540  using result_view_type = Kokkos::View<value_type, Space>;
541 
542  private:
543  result_view_type value;
544  bool references_scalar_v;
545 
546  public:
547  KOKKOS_INLINE_FUNCTION
548  MinMax(value_type& value_) : value(&value_), references_scalar_v(true) {}
549 
550  KOKKOS_INLINE_FUNCTION
551  MinMax(const result_view_type& value_)
552  : value(value_), references_scalar_v(false) {}
553 
554  // Required
555  KOKKOS_INLINE_FUNCTION
556  void join(value_type& dest, const value_type& src) const {
557  if (src.min_val < dest.min_val) {
558  dest.min_val = src.min_val;
559  }
560  if (src.max_val > dest.max_val) {
561  dest.max_val = src.max_val;
562  }
563  }
564 
565  KOKKOS_INLINE_FUNCTION
566  void init(value_type& val) const {
567  val.max_val = reduction_identity<scalar_type>::max();
568  val.min_val = reduction_identity<scalar_type>::min();
569  }
570 
571  KOKKOS_INLINE_FUNCTION
572  value_type& reference() const { return *value.data(); }
573 
574  KOKKOS_INLINE_FUNCTION
575  result_view_type view() const { return value; }
576 
577  KOKKOS_INLINE_FUNCTION
578  bool references_scalar() const { return references_scalar_v; }
579 };
580 
581 template <typename Scalar, typename... Properties>
582 KOKKOS_DEDUCTION_GUIDE MinMax(View<MinMaxScalar<Scalar>, Properties...> const&)
583  ->MinMax<Scalar,
584  typename View<MinMaxScalar<Scalar>, Properties...>::memory_space>;
585 
586 template <class Scalar, class Index>
587 struct MinMaxLocScalar {
588  Scalar min_val, max_val;
589  Index min_loc, max_loc;
590 };
591 
592 template <class Scalar, class Index, class Space>
593 struct MinMaxLoc {
594  private:
595  using scalar_type = std::remove_cv_t<Scalar>;
596  using index_type = std::remove_cv_t<Index>;
597  static_assert(!std::is_pointer_v<scalar_type> &&
598  !std::is_array_v<scalar_type>);
599 
600  public:
601  // Required
602  using reducer = MinMaxLoc<Scalar, Index, Space>;
603  using value_type = MinMaxLocScalar<scalar_type, index_type>;
604 
605  using result_view_type = Kokkos::View<value_type, Space>;
606 
607  private:
608  result_view_type value;
609  bool references_scalar_v;
610 
611  public:
612  KOKKOS_INLINE_FUNCTION
613  MinMaxLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
614 
615  KOKKOS_INLINE_FUNCTION
616  MinMaxLoc(const result_view_type& value_)
617  : value(value_), references_scalar_v(false) {}
618 
619  // Required
620  KOKKOS_INLINE_FUNCTION
621  void join(value_type& dest, const value_type& src) const {
622  if (src.min_val < dest.min_val) {
623  dest.min_val = src.min_val;
624  dest.min_loc = src.min_loc;
625  }
626  if (src.max_val > dest.max_val) {
627  dest.max_val = src.max_val;
628  dest.max_loc = src.max_loc;
629  }
630  }
631 
632  KOKKOS_INLINE_FUNCTION
633  void init(value_type& val) const {
634  val.max_val = reduction_identity<scalar_type>::max();
635  val.min_val = reduction_identity<scalar_type>::min();
636  val.max_loc = reduction_identity<index_type>::min();
637  val.min_loc = reduction_identity<index_type>::min();
638  }
639 
640  KOKKOS_INLINE_FUNCTION
641  value_type& reference() const { return *value.data(); }
642 
643  KOKKOS_INLINE_FUNCTION
644  result_view_type view() const { return value; }
645 
646  KOKKOS_INLINE_FUNCTION
647  bool references_scalar() const { return references_scalar_v; }
648 };
649 
650 template <typename Scalar, typename Index, typename... Properties>
651 KOKKOS_DEDUCTION_GUIDE MinMaxLoc(
652  View<MinMaxLocScalar<Scalar, Index>, Properties...> const&)
653  ->MinMaxLoc<Scalar, Index,
654  typename View<MinMaxLocScalar<Scalar, Index>,
655  Properties...>::memory_space>;
656 
657 // --------------------------------------------------
658 // reducers added to support std algorithms
659 // --------------------------------------------------
660 
661 //
662 // MaxFirstLoc
663 //
664 template <class Scalar, class Index, class Space>
665 struct MaxFirstLoc {
666  private:
667  using scalar_type = std::remove_cv_t<Scalar>;
668  using index_type = std::remove_cv_t<Index>;
669  static_assert(!std::is_pointer_v<scalar_type> &&
670  !std::is_array_v<scalar_type>);
671  static_assert(std::is_integral_v<index_type>);
672 
673  public:
674  // Required
675  using reducer = MaxFirstLoc<Scalar, Index, Space>;
676  using value_type = ::Kokkos::ValLocScalar<scalar_type, index_type>;
677 
678  using result_view_type = ::Kokkos::View<value_type, Space>;
679 
680  private:
681  result_view_type value;
682  bool references_scalar_v;
683 
684  public:
685  KOKKOS_INLINE_FUNCTION
686  MaxFirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
687 
688  KOKKOS_INLINE_FUNCTION
689  MaxFirstLoc(const result_view_type& value_)
690  : value(value_), references_scalar_v(false) {}
691 
692  // Required
693  KOKKOS_INLINE_FUNCTION
694  void join(value_type& dest, const value_type& src) const {
695  if (dest.val < src.val) {
696  dest = src;
697  } else if (!(src.val < dest.val)) {
698  dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc;
699  }
700  }
701 
702  KOKKOS_INLINE_FUNCTION
703  void init(value_type& val) const {
704  val.val = reduction_identity<scalar_type>::max();
705  val.loc = reduction_identity<index_type>::min();
706  }
707 
708  KOKKOS_INLINE_FUNCTION
709  value_type& reference() const { return *value.data(); }
710 
711  KOKKOS_INLINE_FUNCTION
712  result_view_type view() const { return value; }
713 
714  KOKKOS_INLINE_FUNCTION
715  bool references_scalar() const { return references_scalar_v; }
716 };
717 
718 template <typename Scalar, typename Index, typename... Properties>
719 KOKKOS_DEDUCTION_GUIDE MaxFirstLoc(
720  View<ValLocScalar<Scalar, Index>, Properties...> const&)
721  ->MaxFirstLoc<Scalar, Index,
722  typename View<ValLocScalar<Scalar, Index>,
723  Properties...>::memory_space>;
724 
725 //
726 // MaxFirstLocCustomComparator
727 // recall that comp(a,b) returns true is a < b
728 //
729 template <class Scalar, class Index, class ComparatorType, class Space>
730 struct MaxFirstLocCustomComparator {
731  private:
732  using scalar_type = std::remove_cv_t<Scalar>;
733  using index_type = std::remove_cv_t<Index>;
734  static_assert(!std::is_pointer_v<scalar_type> &&
735  !std::is_array_v<scalar_type>);
736  static_assert(std::is_integral_v<index_type>);
737 
738  public:
739  // Required
740  using reducer =
741  MaxFirstLocCustomComparator<Scalar, Index, ComparatorType, Space>;
742  using value_type = ::Kokkos::ValLocScalar<scalar_type, index_type>;
743 
744  using result_view_type = ::Kokkos::View<value_type, Space>;
745 
746  private:
747  result_view_type value;
748  bool references_scalar_v;
749  ComparatorType m_comp;
750 
751  public:
752  KOKKOS_INLINE_FUNCTION
753  MaxFirstLocCustomComparator(value_type& value_, ComparatorType comp_)
754  : value(&value_), references_scalar_v(true), m_comp(comp_) {}
755 
756  KOKKOS_INLINE_FUNCTION
757  MaxFirstLocCustomComparator(const result_view_type& value_,
758  ComparatorType comp_)
759  : value(value_), references_scalar_v(false), m_comp(comp_) {}
760 
761  // Required
762  KOKKOS_INLINE_FUNCTION
763  void join(value_type& dest, const value_type& src) const {
764  if (m_comp(dest.val, src.val)) {
765  dest = src;
766  } else if (!m_comp(src.val, dest.val)) {
767  dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc;
768  }
769  }
770 
771  KOKKOS_INLINE_FUNCTION
772  void init(value_type& val) const {
773  val.val = reduction_identity<scalar_type>::max();
774  val.loc = reduction_identity<index_type>::min();
775  }
776 
777  KOKKOS_INLINE_FUNCTION
778  value_type& reference() const { return *value.data(); }
779 
780  KOKKOS_INLINE_FUNCTION
781  result_view_type view() const { return value; }
782 
783  KOKKOS_INLINE_FUNCTION
784  bool references_scalar() const { return references_scalar_v; }
785 };
786 
787 template <typename Scalar, typename Index, typename ComparatorType,
788  typename... Properties>
789 KOKKOS_DEDUCTION_GUIDE MaxFirstLocCustomComparator(
790  View<ValLocScalar<Scalar, Index>, Properties...> const&, ComparatorType)
791  ->MaxFirstLocCustomComparator<Scalar, Index, ComparatorType,
792  typename View<ValLocScalar<Scalar, Index>,
793  Properties...>::memory_space>;
794 
795 //
796 // MinFirstLoc
797 //
798 template <class Scalar, class Index, class Space>
799 struct MinFirstLoc {
800  private:
801  using scalar_type = std::remove_cv_t<Scalar>;
802  using index_type = std::remove_cv_t<Index>;
803  static_assert(!std::is_pointer_v<scalar_type> &&
804  !std::is_array_v<scalar_type>);
805  static_assert(std::is_integral_v<index_type>);
806 
807  public:
808  // Required
809  using reducer = MinFirstLoc<Scalar, Index, Space>;
810  using value_type = ::Kokkos::ValLocScalar<scalar_type, index_type>;
811 
812  using result_view_type = ::Kokkos::View<value_type, Space>;
813 
814  private:
815  result_view_type value;
816  bool references_scalar_v;
817 
818  public:
819  KOKKOS_INLINE_FUNCTION
820  MinFirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
821 
822  KOKKOS_INLINE_FUNCTION
823  MinFirstLoc(const result_view_type& value_)
824  : value(value_), references_scalar_v(false) {}
825 
826  // Required
827  KOKKOS_INLINE_FUNCTION
828  void join(value_type& dest, const value_type& src) const {
829  if (src.val < dest.val) {
830  dest = src;
831  } else if (!(dest.val < src.val)) {
832  dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc;
833  }
834  }
835 
836  KOKKOS_INLINE_FUNCTION
837  void init(value_type& val) const {
838  val.val = reduction_identity<scalar_type>::min();
839  val.loc = reduction_identity<index_type>::min();
840  }
841 
842  KOKKOS_INLINE_FUNCTION
843  value_type& reference() const { return *value.data(); }
844 
845  KOKKOS_INLINE_FUNCTION
846  result_view_type view() const { return value; }
847 
848  KOKKOS_INLINE_FUNCTION
849  bool references_scalar() const { return references_scalar_v; }
850 };
851 
852 template <typename Scalar, typename Index, typename... Properties>
853 KOKKOS_DEDUCTION_GUIDE MinFirstLoc(
854  View<ValLocScalar<Scalar, Index>, Properties...> const&)
855  ->MinFirstLoc<Scalar, Index,
856  typename View<ValLocScalar<Scalar, Index>,
857  Properties...>::memory_space>;
858 
859 //
860 // MinFirstLocCustomComparator
861 // recall that comp(a,b) returns true is a < b
862 //
863 template <class Scalar, class Index, class ComparatorType, class Space>
864 struct MinFirstLocCustomComparator {
865  private:
866  using scalar_type = std::remove_cv_t<Scalar>;
867  using index_type = std::remove_cv_t<Index>;
868  static_assert(!std::is_pointer_v<scalar_type> &&
869  !std::is_array_v<scalar_type>);
870  static_assert(std::is_integral_v<index_type>);
871 
872  public:
873  // Required
874  using reducer =
875  MinFirstLocCustomComparator<Scalar, Index, ComparatorType, Space>;
876  using value_type = ::Kokkos::ValLocScalar<scalar_type, index_type>;
877 
878  using result_view_type = ::Kokkos::View<value_type, Space>;
879 
880  private:
881  result_view_type value;
882  bool references_scalar_v;
883  ComparatorType m_comp;
884 
885  public:
886  KOKKOS_INLINE_FUNCTION
887  MinFirstLocCustomComparator(value_type& value_, ComparatorType comp_)
888  : value(&value_), references_scalar_v(true), m_comp(comp_) {}
889 
890  KOKKOS_INLINE_FUNCTION
891  MinFirstLocCustomComparator(const result_view_type& value_,
892  ComparatorType comp_)
893  : value(value_), references_scalar_v(false), m_comp(comp_) {}
894 
895  // Required
896  KOKKOS_INLINE_FUNCTION
897  void join(value_type& dest, const value_type& src) const {
898  if (m_comp(src.val, dest.val)) {
899  dest = src;
900  } else if (!m_comp(dest.val, src.val)) {
901  dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc;
902  }
903  }
904 
905  KOKKOS_INLINE_FUNCTION
906  void init(value_type& val) const {
907  val.val = reduction_identity<scalar_type>::min();
908  val.loc = reduction_identity<index_type>::min();
909  }
910 
911  KOKKOS_INLINE_FUNCTION
912  value_type& reference() const { return *value.data(); }
913 
914  KOKKOS_INLINE_FUNCTION
915  result_view_type view() const { return value; }
916 
917  KOKKOS_INLINE_FUNCTION
918  bool references_scalar() const { return references_scalar_v; }
919 };
920 
921 template <typename Scalar, typename Index, typename ComparatorType,
922  typename... Properties>
923 KOKKOS_DEDUCTION_GUIDE MinFirstLocCustomComparator(
924  View<ValLocScalar<Scalar, Index>, Properties...> const&, ComparatorType)
925  ->MinFirstLocCustomComparator<Scalar, Index, ComparatorType,
926  typename View<ValLocScalar<Scalar, Index>,
927  Properties...>::memory_space>;
928 
929 //
930 // MinMaxFirstLastLoc
931 //
932 template <class Scalar, class Index, class Space>
933 struct MinMaxFirstLastLoc {
934  private:
935  using scalar_type = std::remove_cv_t<Scalar>;
936  using index_type = std::remove_cv_t<Index>;
937  static_assert(!std::is_pointer_v<scalar_type> &&
938  !std::is_array_v<scalar_type>);
939  static_assert(std::is_integral_v<index_type>);
940 
941  public:
942  // Required
943  using reducer = MinMaxFirstLastLoc<Scalar, Index, Space>;
944  using value_type = ::Kokkos::MinMaxLocScalar<scalar_type, index_type>;
945 
946  using result_view_type = ::Kokkos::View<value_type, Space>;
947 
948  private:
949  result_view_type value;
950  bool references_scalar_v;
951 
952  public:
953  KOKKOS_INLINE_FUNCTION
954  MinMaxFirstLastLoc(value_type& value_)
955  : value(&value_), references_scalar_v(true) {}
956 
957  KOKKOS_INLINE_FUNCTION
958  MinMaxFirstLastLoc(const result_view_type& value_)
959  : value(value_), references_scalar_v(false) {}
960 
961  // Required
962  KOKKOS_INLINE_FUNCTION
963  void join(value_type& dest, const value_type& src) const {
964  if (src.min_val < dest.min_val) {
965  dest.min_val = src.min_val;
966  dest.min_loc = src.min_loc;
967  } else if (!(dest.min_val < src.min_val)) {
968  dest.min_loc = (src.min_loc < dest.min_loc) ? src.min_loc : dest.min_loc;
969  }
970 
971  if (dest.max_val < src.max_val) {
972  dest.max_val = src.max_val;
973  dest.max_loc = src.max_loc;
974  } else if (!(src.max_val < dest.max_val)) {
975  dest.max_loc = (src.max_loc > dest.max_loc) ? src.max_loc : dest.max_loc;
976  }
977  }
978 
979  KOKKOS_INLINE_FUNCTION
980  void init(value_type& val) const {
981  val.max_val = ::Kokkos::reduction_identity<scalar_type>::max();
982  val.min_val = ::Kokkos::reduction_identity<scalar_type>::min();
983  val.max_loc = ::Kokkos::reduction_identity<index_type>::max();
984  val.min_loc = ::Kokkos::reduction_identity<index_type>::min();
985  }
986 
987  KOKKOS_INLINE_FUNCTION
988  value_type& reference() const { return *value.data(); }
989 
990  KOKKOS_INLINE_FUNCTION
991  result_view_type view() const { return value; }
992 
993  KOKKOS_INLINE_FUNCTION
994  bool references_scalar() const { return references_scalar_v; }
995 };
996 
997 template <typename Scalar, typename Index, typename... Properties>
998 KOKKOS_DEDUCTION_GUIDE MinMaxFirstLastLoc(
999  View<MinMaxLocScalar<Scalar, Index>, Properties...> const&)
1000  ->MinMaxFirstLastLoc<Scalar, Index,
1001  typename View<MinMaxLocScalar<Scalar, Index>,
1002  Properties...>::memory_space>;
1003 
1004 //
1005 // MinMaxFirstLastLocCustomComparator
1006 // recall that comp(a,b) returns true is a < b
1007 //
1008 template <class Scalar, class Index, class ComparatorType, class Space>
1009 struct MinMaxFirstLastLocCustomComparator {
1010  private:
1011  using scalar_type = std::remove_cv_t<Scalar>;
1012  using index_type = std::remove_cv_t<Index>;
1013  static_assert(!std::is_pointer_v<scalar_type> &&
1014  !std::is_array_v<scalar_type>);
1015  static_assert(std::is_integral_v<index_type>);
1016 
1017  public:
1018  // Required
1019  using reducer =
1020  MinMaxFirstLastLocCustomComparator<Scalar, Index, ComparatorType, Space>;
1021  using value_type = ::Kokkos::MinMaxLocScalar<scalar_type, index_type>;
1022 
1023  using result_view_type = ::Kokkos::View<value_type, Space>;
1024 
1025  private:
1026  result_view_type value;
1027  bool references_scalar_v;
1028  ComparatorType m_comp;
1029 
1030  public:
1031  KOKKOS_INLINE_FUNCTION
1032  MinMaxFirstLastLocCustomComparator(value_type& value_, ComparatorType comp_)
1033  : value(&value_), references_scalar_v(true), m_comp(comp_) {}
1034 
1035  KOKKOS_INLINE_FUNCTION
1036  MinMaxFirstLastLocCustomComparator(const result_view_type& value_,
1037  ComparatorType comp_)
1038  : value(value_), references_scalar_v(false), m_comp(comp_) {}
1039 
1040  // Required
1041  KOKKOS_INLINE_FUNCTION
1042  void join(value_type& dest, const value_type& src) const {
1043  if (m_comp(src.min_val, dest.min_val)) {
1044  dest.min_val = src.min_val;
1045  dest.min_loc = src.min_loc;
1046  } else if (!m_comp(dest.min_val, src.min_val)) {
1047  dest.min_loc = (src.min_loc < dest.min_loc) ? src.min_loc : dest.min_loc;
1048  }
1049 
1050  if (m_comp(dest.max_val, src.max_val)) {
1051  dest.max_val = src.max_val;
1052  dest.max_loc = src.max_loc;
1053  } else if (!m_comp(src.max_val, dest.max_val)) {
1054  dest.max_loc = (src.max_loc > dest.max_loc) ? src.max_loc : dest.max_loc;
1055  }
1056  }
1057 
1058  KOKKOS_INLINE_FUNCTION
1059  void init(value_type& val) const {
1060  val.max_val = ::Kokkos::reduction_identity<scalar_type>::max();
1061  val.min_val = ::Kokkos::reduction_identity<scalar_type>::min();
1062  val.max_loc = ::Kokkos::reduction_identity<index_type>::max();
1063  val.min_loc = ::Kokkos::reduction_identity<index_type>::min();
1064  }
1065 
1066  KOKKOS_INLINE_FUNCTION
1067  value_type& reference() const { return *value.data(); }
1068 
1069  KOKKOS_INLINE_FUNCTION
1070  result_view_type view() const { return value; }
1071 
1072  KOKKOS_INLINE_FUNCTION
1073  bool references_scalar() const { return references_scalar_v; }
1074 };
1075 
1076 template <typename Scalar, typename Index, typename ComparatorType,
1077  typename... Properties>
1078 KOKKOS_DEDUCTION_GUIDE MinMaxFirstLastLocCustomComparator(
1079  View<MinMaxLocScalar<Scalar, Index>, Properties...> const&, ComparatorType)
1080  ->MinMaxFirstLastLocCustomComparator<
1081  Scalar, Index, ComparatorType,
1082  typename View<MinMaxLocScalar<Scalar, Index>,
1083  Properties...>::memory_space>;
1084 
1085 //
1086 // FirstLoc
1087 //
1088 template <class Index>
1089 struct FirstLocScalar {
1090  Index min_loc_true;
1091 };
1092 
1093 template <class Index, class Space>
1094 struct FirstLoc {
1095  private:
1096  using index_type = std::remove_cv_t<Index>;
1097  static_assert(std::is_integral_v<index_type>);
1098 
1099  public:
1100  // Required
1101  using reducer = FirstLoc<Index, Space>;
1102  using value_type = FirstLocScalar<index_type>;
1103 
1104  using result_view_type = ::Kokkos::View<value_type, Space>;
1105 
1106  private:
1107  result_view_type value;
1108  bool references_scalar_v;
1109 
1110  public:
1111  KOKKOS_INLINE_FUNCTION
1112  FirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
1113 
1114  KOKKOS_INLINE_FUNCTION
1115  FirstLoc(const result_view_type& value_)
1116  : value(value_), references_scalar_v(false) {}
1117 
1118  // Required
1119  KOKKOS_INLINE_FUNCTION
1120  void join(value_type& dest, const value_type& src) const {
1121  dest.min_loc_true = (src.min_loc_true < dest.min_loc_true)
1122  ? src.min_loc_true
1123  : dest.min_loc_true;
1124  }
1125 
1126  KOKKOS_INLINE_FUNCTION
1127  void init(value_type& val) const {
1128  val.min_loc_true = ::Kokkos::reduction_identity<index_type>::min();
1129  }
1130 
1131  KOKKOS_INLINE_FUNCTION
1132  value_type& reference() const { return *value.data(); }
1133 
1134  KOKKOS_INLINE_FUNCTION
1135  result_view_type view() const { return value; }
1136 
1137  KOKKOS_INLINE_FUNCTION
1138  bool references_scalar() const { return references_scalar_v; }
1139 };
1140 
1141 template <typename Index, typename... Properties>
1142 KOKKOS_DEDUCTION_GUIDE FirstLoc(
1143  View<FirstLocScalar<Index>, Properties...> const&)
1144  ->FirstLoc<Index, typename View<FirstLocScalar<Index>,
1145  Properties...>::memory_space>;
1146 
1147 //
1148 // LastLoc
1149 //
1150 template <class Index>
1151 struct LastLocScalar {
1152  Index max_loc_true;
1153 };
1154 
1155 template <class Index, class Space>
1156 struct LastLoc {
1157  private:
1158  using index_type = std::remove_cv_t<Index>;
1159  static_assert(std::is_integral_v<index_type>);
1160 
1161  public:
1162  // Required
1163  using reducer = LastLoc<Index, Space>;
1164  using value_type = LastLocScalar<index_type>;
1165 
1166  using result_view_type = ::Kokkos::View<value_type, Space>;
1167 
1168  private:
1169  result_view_type value;
1170  bool references_scalar_v;
1171 
1172  public:
1173  KOKKOS_INLINE_FUNCTION
1174  LastLoc(value_type& value_) : value(&value_), references_scalar_v(true) {}
1175 
1176  KOKKOS_INLINE_FUNCTION
1177  LastLoc(const result_view_type& value_)
1178  : value(value_), references_scalar_v(false) {}
1179 
1180  // Required
1181  KOKKOS_INLINE_FUNCTION
1182  void join(value_type& dest, const value_type& src) const {
1183  dest.max_loc_true = (src.max_loc_true > dest.max_loc_true)
1184  ? src.max_loc_true
1185  : dest.max_loc_true;
1186  }
1187 
1188  KOKKOS_INLINE_FUNCTION
1189  void init(value_type& val) const {
1190  val.max_loc_true = ::Kokkos::reduction_identity<index_type>::max();
1191  }
1192 
1193  KOKKOS_INLINE_FUNCTION
1194  value_type& reference() const { return *value.data(); }
1195 
1196  KOKKOS_INLINE_FUNCTION
1197  result_view_type view() const { return value; }
1198 
1199  KOKKOS_INLINE_FUNCTION
1200  bool references_scalar() const { return references_scalar_v; }
1201 };
1202 
1203 template <typename Index, typename... Properties>
1204 KOKKOS_DEDUCTION_GUIDE LastLoc(View<LastLocScalar<Index>, Properties...> const&)
1205  ->LastLoc<Index,
1206  typename View<LastLocScalar<Index>, Properties...>::memory_space>;
1207 
1208 template <class Index>
1209 struct StdIsPartScalar {
1210  Index max_loc_true, min_loc_false;
1211 };
1212 
1213 //
1214 // StdIsPartitioned
1215 //
1216 template <class Index, class Space>
1217 struct StdIsPartitioned {
1218  private:
1219  using index_type = std::remove_cv_t<Index>;
1220  static_assert(std::is_integral_v<index_type>);
1221 
1222  public:
1223  // Required
1224  using reducer = StdIsPartitioned<Index, Space>;
1225  using value_type = StdIsPartScalar<index_type>;
1226 
1227  using result_view_type = ::Kokkos::View<value_type, Space>;
1228 
1229  private:
1230  result_view_type value;
1231  bool references_scalar_v;
1232 
1233  public:
1234  KOKKOS_INLINE_FUNCTION
1235  StdIsPartitioned(value_type& value_)
1236  : value(&value_), references_scalar_v(true) {}
1237 
1238  KOKKOS_INLINE_FUNCTION
1239  StdIsPartitioned(const result_view_type& value_)
1240  : value(value_), references_scalar_v(false) {}
1241 
1242  // Required
1243  KOKKOS_INLINE_FUNCTION
1244  void join(value_type& dest, const value_type& src) const {
1245  dest.max_loc_true = (dest.max_loc_true < src.max_loc_true)
1246  ? src.max_loc_true
1247  : dest.max_loc_true;
1248 
1249  dest.min_loc_false = (dest.min_loc_false < src.min_loc_false)
1250  ? dest.min_loc_false
1251  : src.min_loc_false;
1252  }
1253 
1254  KOKKOS_INLINE_FUNCTION
1255  void init(value_type& val) const {
1256  val.max_loc_true = ::Kokkos::reduction_identity<index_type>::max();
1257  val.min_loc_false = ::Kokkos::reduction_identity<index_type>::min();
1258  }
1259 
1260  KOKKOS_INLINE_FUNCTION
1261  value_type& reference() const { return *value.data(); }
1262 
1263  KOKKOS_INLINE_FUNCTION
1264  result_view_type view() const { return value; }
1265 
1266  KOKKOS_INLINE_FUNCTION
1267  bool references_scalar() const { return references_scalar_v; }
1268 };
1269 
1270 template <typename Index, typename... Properties>
1271 KOKKOS_DEDUCTION_GUIDE StdIsPartitioned(
1272  View<StdIsPartScalar<Index>, Properties...> const&)
1273  ->StdIsPartitioned<Index, typename View<StdIsPartScalar<Index>,
1274  Properties...>::memory_space>;
1275 
1276 template <class Index>
1277 struct StdPartPointScalar {
1278  Index min_loc_false;
1279 };
1280 
1281 //
1282 // StdPartitionPoint
1283 //
1284 template <class Index, class Space>
1285 struct StdPartitionPoint {
1286  private:
1287  using index_type = std::remove_cv_t<Index>;
1288  static_assert(std::is_integral_v<index_type>);
1289 
1290  public:
1291  // Required
1292  using reducer = StdPartitionPoint<Index, Space>;
1293  using value_type = StdPartPointScalar<index_type>;
1294 
1295  using result_view_type = ::Kokkos::View<value_type, Space>;
1296 
1297  private:
1298  result_view_type value;
1299  bool references_scalar_v;
1300 
1301  public:
1302  KOKKOS_INLINE_FUNCTION
1303  StdPartitionPoint(value_type& value_)
1304  : value(&value_), references_scalar_v(true) {}
1305 
1306  KOKKOS_INLINE_FUNCTION
1307  StdPartitionPoint(const result_view_type& value_)
1308  : value(value_), references_scalar_v(false) {}
1309 
1310  // Required
1311  KOKKOS_INLINE_FUNCTION
1312  void join(value_type& dest, const value_type& src) const {
1313  dest.min_loc_false = (dest.min_loc_false < src.min_loc_false)
1314  ? dest.min_loc_false
1315  : src.min_loc_false;
1316  }
1317 
1318  KOKKOS_INLINE_FUNCTION
1319  void init(value_type& val) const {
1320  val.min_loc_false = ::Kokkos::reduction_identity<index_type>::min();
1321  }
1322 
1323  KOKKOS_INLINE_FUNCTION
1324  value_type& reference() const { return *value.data(); }
1325 
1326  KOKKOS_INLINE_FUNCTION
1327  result_view_type view() const { return value; }
1328 
1329  KOKKOS_INLINE_FUNCTION
1330  bool references_scalar() const { return references_scalar_v; }
1331 };
1332 
1333 template <typename Index, typename... Properties>
1334 KOKKOS_DEDUCTION_GUIDE StdPartitionPoint(
1335  View<StdPartPointScalar<Index>, Properties...> const&)
1336  ->StdPartitionPoint<Index, typename View<StdPartPointScalar<Index>,
1337  Properties...>::memory_space>;
1338 
1339 } // namespace Kokkos
1340 namespace Kokkos {
1341 namespace Impl {
1342 
1343 template <typename FunctorType, typename FunctorAnalysisReducerType,
1344  typename Enable>
1345 class CombinedFunctorReducer {
1346  public:
1347  using functor_type = FunctorType;
1348  using reducer_type = FunctorAnalysisReducerType;
1349  CombinedFunctorReducer(const FunctorType& functor,
1350  const FunctorAnalysisReducerType& reducer)
1351  : m_functor(functor), m_reducer(reducer) {}
1352  KOKKOS_FUNCTION const FunctorType& get_functor() const { return m_functor; }
1353  KOKKOS_FUNCTION const FunctorAnalysisReducerType& get_reducer() const {
1354  return m_reducer;
1355  }
1356 
1357  private:
1358  FunctorType m_functor;
1359  FunctorAnalysisReducerType m_reducer;
1360 };
1361 template <typename FunctorType, typename FunctorAnalysisReducerType>
1362 class CombinedFunctorReducer<
1363  FunctorType, FunctorAnalysisReducerType,
1364  std::enable_if_t<std::is_same_v<
1365  FunctorType, typename FunctorAnalysisReducerType::functor_type>>> {
1366  public:
1367  using functor_type = FunctorType;
1368  using reducer_type = FunctorAnalysisReducerType;
1369  CombinedFunctorReducer(const FunctorType& functor,
1370  const FunctorAnalysisReducerType&)
1371  : m_reducer(functor) {}
1372  KOKKOS_FUNCTION const FunctorType& get_functor() const {
1373  return m_reducer.get_functor();
1374  }
1375  KOKKOS_FUNCTION const FunctorAnalysisReducerType& get_reducer() const {
1376  return m_reducer;
1377  }
1378 
1379  private:
1380  FunctorAnalysisReducerType m_reducer;
1381 };
1382 
1383 template <class T, class ReturnType, class ValueTraits>
1384 struct ParallelReduceReturnValue;
1385 
1386 template <class ReturnType, class FunctorType>
1387 struct ParallelReduceReturnValue<
1388  std::enable_if_t<Kokkos::is_view<ReturnType>::value>, ReturnType,
1389  FunctorType> {
1390  using return_type = ReturnType;
1391  using reducer_type = InvalidType;
1392 
1393  using value_type_scalar = typename return_type::value_type;
1394  using value_type_array = typename return_type::value_type* const;
1395 
1396  using value_type = std::conditional_t<return_type::rank == 0,
1397  value_type_scalar, value_type_array>;
1398 
1399  static return_type& return_value(ReturnType& return_val, const FunctorType&) {
1400  return return_val;
1401  }
1402 };
1403 
1404 template <class ReturnType, class FunctorType>
1405 struct ParallelReduceReturnValue<
1406  std::enable_if_t<!Kokkos::is_view<ReturnType>::value &&
1407  (!std::is_array<ReturnType>::value &&
1408  !std::is_pointer<ReturnType>::value) &&
1409  !Kokkos::is_reducer<ReturnType>::value>,
1410  ReturnType, FunctorType> {
1411  using return_type =
1413 
1414  using reducer_type = InvalidType;
1415 
1416  using value_type = typename return_type::value_type;
1417 
1418  static return_type return_value(ReturnType& return_val, const FunctorType&) {
1419  return return_type(&return_val);
1420  }
1421 };
1422 
1423 template <class ReturnType, class FunctorType>
1424 struct ParallelReduceReturnValue<
1425  std::enable_if_t<(std::is_array<ReturnType>::value ||
1426  std::is_pointer<ReturnType>::value)>,
1427  ReturnType, FunctorType> {
1429  Kokkos::HostSpace, Kokkos::MemoryUnmanaged>;
1430 
1431  using reducer_type = InvalidType;
1432 
1433  using value_type = typename return_type::value_type[];
1434 
1435  static return_type return_value(ReturnType& return_val,
1436  const FunctorType& functor) {
1437  if (std::is_array<ReturnType>::value)
1438  return return_type(return_val);
1439  else
1440  return return_type(return_val, functor.value_count);
1441  }
1442 };
1443 
1444 template <class ReturnType, class FunctorType>
1445 struct ParallelReduceReturnValue<
1446  std::enable_if_t<Kokkos::is_reducer<ReturnType>::value>, ReturnType,
1447  FunctorType> {
1448  using return_type = typename ReturnType::result_view_type;
1449  using reducer_type = ReturnType;
1450  using value_type = typename return_type::value_type;
1451 
1452  static auto return_value(ReturnType& return_val, const FunctorType&) {
1453  return return_val.view();
1454  }
1455 };
1456 
1457 template <class T, class ReturnType, class FunctorType>
1458 struct ParallelReducePolicyType;
1459 
1460 template <class PolicyType, class FunctorType>
1461 struct ParallelReducePolicyType<
1462  std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value>,
1463  PolicyType, FunctorType> {
1464  using policy_type = PolicyType;
1465  static PolicyType policy(const PolicyType& policy_) { return policy_; }
1466 };
1467 
1468 template <class PolicyType, class FunctorType>
1469 struct ParallelReducePolicyType<
1470  std::enable_if_t<std::is_integral<PolicyType>::value>, PolicyType,
1471  FunctorType> {
1472  using execution_space =
1473  typename Impl::FunctorPolicyExecutionSpace<FunctorType,
1474  void>::execution_space;
1475 
1476  using policy_type = Kokkos::RangePolicy<execution_space>;
1477 
1478  static policy_type policy(const PolicyType& policy_) {
1479  return policy_type(0, policy_);
1480  }
1481 };
1482 
1483 template <class FunctorType, class ExecPolicy, class ValueType,
1484  class ExecutionSpace>
1485 struct ParallelReduceFunctorType {
1486  using functor_type = FunctorType;
1487  static const functor_type& functor(const functor_type& functor) {
1488  return functor;
1489  }
1490 };
1491 
1492 template <class PolicyType, class FunctorType, class ReturnType>
1493 struct ParallelReduceAdaptor {
1494  using return_value_adapter =
1495  Impl::ParallelReduceReturnValue<void, ReturnType, FunctorType>;
1496 
1497  static inline void execute_impl(const std::string& label,
1498  const PolicyType& policy,
1499  const FunctorType& functor,
1500  ReturnType& return_value) {
1501  using PassedReducerType = typename return_value_adapter::reducer_type;
1502  uint64_t kpID = 0;
1503 
1504  PolicyType inner_policy = policy;
1505  Kokkos::Tools::Impl::begin_parallel_reduce<PassedReducerType>(
1506  inner_policy, functor, label, kpID);
1507 
1508  using ReducerSelector =
1509  Kokkos::Impl::if_c<std::is_same<InvalidType, PassedReducerType>::value,
1510  FunctorType, PassedReducerType>;
1511  using Analysis = FunctorAnalysis<FunctorPatternInterface::REDUCE,
1512  PolicyType, typename ReducerSelector::type,
1513  typename return_value_adapter::value_type>;
1514 
1515  using CombinedFunctorReducerType =
1516  CombinedFunctorReducer<FunctorType, typename Analysis::Reducer>;
1517  auto closure = construct_with_shared_allocation_tracking_disabled<
1518  Impl::ParallelReduce<CombinedFunctorReducerType, PolicyType,
1519  typename Impl::FunctorPolicyExecutionSpace<
1520  FunctorType, PolicyType>::execution_space>>(
1521  CombinedFunctorReducerType(
1522  functor, typename Analysis::Reducer(
1523  ReducerSelector::select(functor, return_value))),
1524  inner_policy,
1525  return_value_adapter::return_value(return_value, functor));
1526  closure.execute();
1527 
1528  Kokkos::Tools::Impl::end_parallel_reduce<PassedReducerType>(
1529  inner_policy, functor, label, kpID);
1530  }
1531 
1532  static constexpr bool is_array_reduction =
1533  Impl::FunctorAnalysis<
1534  Impl::FunctorPatternInterface::REDUCE, PolicyType, FunctorType,
1535  typename return_value_adapter::value_type>::StaticValueSize == 0;
1536 
1537  template <typename Dummy = ReturnType>
1538  static inline std::enable_if_t<!(is_array_reduction &&
1539  std::is_pointer<Dummy>::value)>
1540  execute(const std::string& label, const PolicyType& policy,
1541  const FunctorType& functor, ReturnType& return_value) {
1542  execute_impl(label, policy, functor, return_value);
1543  }
1544 };
1545 } // namespace Impl
1546 
1547 //----------------------------------------------------------------------------
1548 
1560 // Parallel Reduce Blocking behavior
1561 
1562 namespace Impl {
1563 template <typename T>
1564 struct ReducerHasTestReferenceFunction {
1565  template <typename E>
1566  static std::true_type test_func(decltype(&E::references_scalar));
1567  template <typename E>
1568  static std::false_type test_func(...);
1569 
1570  enum {
1571  value = std::is_same<std::true_type, decltype(test_func<T>(nullptr))>::value
1572  };
1573 };
1574 
1575 template <class ExecutionSpace, class Arg>
1576 constexpr std::enable_if_t<
1577  // constraints only necessary because SFINAE lacks subsumption
1578  !ReducerHasTestReferenceFunction<Arg>::value &&
1579  !Kokkos::is_view<Arg>::value,
1580  // return type:
1581  bool>
1582 parallel_reduce_needs_fence(ExecutionSpace const&, Arg const&) {
1583  return true;
1584 }
1585 
1586 template <class ExecutionSpace, class Reducer>
1587 constexpr std::enable_if_t<
1588  // equivalent to:
1589  // (requires (Reducer const& r) {
1590  // { reducer.references_scalar() } -> std::convertible_to<bool>;
1591  // })
1592  ReducerHasTestReferenceFunction<Reducer>::value,
1593  // return type:
1594  bool>
1595 parallel_reduce_needs_fence(ExecutionSpace const&, Reducer const& reducer) {
1596  return reducer.references_scalar();
1597 }
1598 
1599 template <class ExecutionSpace, class ViewLike>
1600 constexpr std::enable_if_t<
1601  // requires Kokkos::ViewLike<ViewLike>
1602  Kokkos::is_view<ViewLike>::value,
1603  // return type:
1604  bool>
1605 parallel_reduce_needs_fence(ExecutionSpace const&, ViewLike const&) {
1606  return false;
1607 }
1608 
1609 template <class ExecutionSpace, class... Args>
1610 struct ParallelReduceFence {
1611  template <class... ArgsDeduced>
1612  static void fence(const ExecutionSpace& ex, const std::string& name,
1613  ArgsDeduced&&... args) {
1614  if (Impl::parallel_reduce_needs_fence(ex, (ArgsDeduced &&) args...)) {
1615  ex.fence(name);
1616  }
1617  }
1618 };
1619 
1620 } // namespace Impl
1621 
1660 // ReturnValue is scalar or array: take by reference
1661 
1662 template <class PolicyType, class FunctorType, class ReturnType>
1663 inline std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value &&
1664  !(Kokkos::is_view<ReturnType>::value ||
1665  Kokkos::is_reducer<ReturnType>::value ||
1666  std::is_pointer<ReturnType>::value)>
1667 parallel_reduce(const std::string& label, const PolicyType& policy,
1668  const FunctorType& functor, ReturnType& return_value) {
1669  static_assert(
1670  !std::is_const<ReturnType>::value,
1671  "A const reduction result type is only allowed for a View, pointer or "
1672  "reducer return type!");
1673 
1674  Impl::ParallelReduceAdaptor<PolicyType, FunctorType, ReturnType>::execute(
1675  label, policy, functor, return_value);
1676  Impl::ParallelReduceFence<typename PolicyType::execution_space, ReturnType>::
1677  fence(
1678  policy.space(),
1679  "Kokkos::parallel_reduce: fence due to result being value, not view",
1680  return_value);
1681 }
1682 
1683 template <class PolicyType, class FunctorType, class ReturnType>
1684 inline std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value &&
1685  !(Kokkos::is_view<ReturnType>::value ||
1686  Kokkos::is_reducer<ReturnType>::value ||
1687  std::is_pointer<ReturnType>::value)>
1688 parallel_reduce(const PolicyType& policy, const FunctorType& functor,
1689  ReturnType& return_value) {
1690  static_assert(
1691  !std::is_const<ReturnType>::value,
1692  "A const reduction result type is only allowed for a View, pointer or "
1693  "reducer return type!");
1694 
1695  Impl::ParallelReduceAdaptor<PolicyType, FunctorType, ReturnType>::execute(
1696  "", policy, functor, return_value);
1697  Impl::ParallelReduceFence<typename PolicyType::execution_space, ReturnType>::
1698  fence(
1699  policy.space(),
1700  "Kokkos::parallel_reduce: fence due to result being value, not view",
1701  return_value);
1702 }
1703 
1704 template <class FunctorType, class ReturnType>
1705 inline std::enable_if_t<!(Kokkos::is_view<ReturnType>::value ||
1706  Kokkos::is_reducer<ReturnType>::value ||
1707  std::is_pointer<ReturnType>::value)>
1708 parallel_reduce(const size_t& policy, const FunctorType& functor,
1709  ReturnType& return_value) {
1710  static_assert(
1711  !std::is_const<ReturnType>::value,
1712  "A const reduction result type is only allowed for a View, pointer or "
1713  "reducer return type!");
1714 
1715  using policy_type =
1716  typename Impl::ParallelReducePolicyType<void, size_t,
1717  FunctorType>::policy_type;
1718 
1719  Impl::ParallelReduceAdaptor<policy_type, FunctorType, ReturnType>::execute(
1720  "", policy_type(0, policy), functor, return_value);
1721  Impl::ParallelReduceFence<typename policy_type::execution_space, ReturnType>::
1722  fence(
1723  typename policy_type::execution_space(),
1724  "Kokkos::parallel_reduce: fence due to result being value, not view",
1725  return_value);
1726 }
1727 
1728 template <class FunctorType, class ReturnType>
1729 inline std::enable_if_t<!(Kokkos::is_view<ReturnType>::value ||
1730  Kokkos::is_reducer<ReturnType>::value ||
1731  std::is_pointer<ReturnType>::value)>
1732 parallel_reduce(const std::string& label, const size_t& policy,
1733  const FunctorType& functor, ReturnType& return_value) {
1734  static_assert(
1735  !std::is_const<ReturnType>::value,
1736  "A const reduction result type is only allowed for a View, pointer or "
1737  "reducer return type!");
1738 
1739  using policy_type =
1740  typename Impl::ParallelReducePolicyType<void, size_t,
1741  FunctorType>::policy_type;
1742  Impl::ParallelReduceAdaptor<policy_type, FunctorType, ReturnType>::execute(
1743  label, policy_type(0, policy), functor, return_value);
1744  Impl::ParallelReduceFence<typename policy_type::execution_space, ReturnType>::
1745  fence(
1746  typename policy_type::execution_space(),
1747  "Kokkos::parallel_reduce: fence due to result being value, not view",
1748  return_value);
1749 }
1750 
1751 // ReturnValue as View or Reducer: take by copy to allow for inline construction
1752 
1753 template <class PolicyType, class FunctorType, class ReturnType>
1754 inline std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value &&
1755  (Kokkos::is_view<ReturnType>::value ||
1756  Kokkos::is_reducer<ReturnType>::value ||
1757  std::is_pointer<ReturnType>::value)>
1758 parallel_reduce(const std::string& label, const PolicyType& policy,
1759  const FunctorType& functor, const ReturnType& return_value) {
1760  ReturnType return_value_impl = return_value;
1761  Impl::ParallelReduceAdaptor<PolicyType, FunctorType, ReturnType>::execute(
1762  label, policy, functor, return_value_impl);
1763  Impl::ParallelReduceFence<typename PolicyType::execution_space, ReturnType>::
1764  fence(
1765  policy.space(),
1766  "Kokkos::parallel_reduce: fence due to result being value, not view",
1767  return_value);
1768 }
1769 
1770 template <class PolicyType, class FunctorType, class ReturnType>
1771 inline std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value &&
1772  (Kokkos::is_view<ReturnType>::value ||
1773  Kokkos::is_reducer<ReturnType>::value ||
1774  std::is_pointer<ReturnType>::value)>
1775 parallel_reduce(const PolicyType& policy, const FunctorType& functor,
1776  const ReturnType& return_value) {
1777  ReturnType return_value_impl = return_value;
1778  Impl::ParallelReduceAdaptor<PolicyType, FunctorType, ReturnType>::execute(
1779  "", policy, functor, return_value_impl);
1780  Impl::ParallelReduceFence<typename PolicyType::execution_space, ReturnType>::
1781  fence(
1782  policy.space(),
1783  "Kokkos::parallel_reduce: fence due to result being value, not view",
1784  return_value);
1785 }
1786 
1787 template <class FunctorType, class ReturnType>
1788 inline std::enable_if_t<Kokkos::is_view<ReturnType>::value ||
1789  Kokkos::is_reducer<ReturnType>::value ||
1790  std::is_pointer<ReturnType>::value>
1791 parallel_reduce(const size_t& policy, const FunctorType& functor,
1792  const ReturnType& return_value) {
1793  using policy_type =
1794  typename Impl::ParallelReducePolicyType<void, size_t,
1795  FunctorType>::policy_type;
1796  ReturnType return_value_impl = return_value;
1797  Impl::ParallelReduceAdaptor<policy_type, FunctorType, ReturnType>::execute(
1798  "", policy_type(0, policy), functor, return_value_impl);
1799  Impl::ParallelReduceFence<typename policy_type::execution_space, ReturnType>::
1800  fence(
1801  typename policy_type::execution_space(),
1802  "Kokkos::parallel_reduce: fence due to result being value, not view",
1803  return_value);
1804 }
1805 
1806 template <class FunctorType, class ReturnType>
1807 inline std::enable_if_t<Kokkos::is_view<ReturnType>::value ||
1808  Kokkos::is_reducer<ReturnType>::value ||
1809  std::is_pointer<ReturnType>::value>
1810 parallel_reduce(const std::string& label, const size_t& policy,
1811  const FunctorType& functor, const ReturnType& return_value) {
1812  using policy_type =
1813  typename Impl::ParallelReducePolicyType<void, size_t,
1814  FunctorType>::policy_type;
1815  ReturnType return_value_impl = return_value;
1816  Impl::ParallelReduceAdaptor<policy_type, FunctorType, ReturnType>::execute(
1817  label, policy_type(0, policy), functor, return_value_impl);
1818  Impl::ParallelReduceFence<typename policy_type::execution_space, ReturnType>::
1819  fence(
1820  typename policy_type::execution_space(),
1821  "Kokkos::parallel_reduce: fence due to result being value, not view",
1822  return_value);
1823 }
1824 
1825 // No Return Argument
1826 
1827 template <class PolicyType, class FunctorType>
1828 inline void parallel_reduce(
1829  const std::string& label, const PolicyType& policy,
1830  const FunctorType& functor,
1831  std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value>* =
1832  nullptr) {
1833  using FunctorAnalysis =
1834  Impl::FunctorAnalysis<Impl::FunctorPatternInterface::REDUCE, PolicyType,
1835  FunctorType, void>;
1836  using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0),
1837  typename FunctorAnalysis::value_type,
1838  typename FunctorAnalysis::pointer_type>;
1839 
1840  static_assert(
1841  FunctorAnalysis::has_final_member_function,
1842  "Calling parallel_reduce without either return value or final function.");
1843 
1844  using result_view_type =
1846  result_view_type result_view;
1847 
1848  Impl::ParallelReduceAdaptor<PolicyType, FunctorType,
1849  result_view_type>::execute(label, policy, functor,
1850  result_view);
1851 }
1852 
1853 template <class PolicyType, class FunctorType>
1854 inline void parallel_reduce(
1855  const PolicyType& policy, const FunctorType& functor,
1856  std::enable_if_t<Kokkos::is_execution_policy<PolicyType>::value>* =
1857  nullptr) {
1858  using FunctorAnalysis =
1859  Impl::FunctorAnalysis<Impl::FunctorPatternInterface::REDUCE, PolicyType,
1860  FunctorType, void>;
1861  using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0),
1862  typename FunctorAnalysis::value_type,
1863  typename FunctorAnalysis::pointer_type>;
1864 
1865  static_assert(
1866  FunctorAnalysis::has_final_member_function,
1867  "Calling parallel_reduce without either return value or final function.");
1868 
1869  using result_view_type =
1871  result_view_type result_view;
1872 
1873  Impl::ParallelReduceAdaptor<PolicyType, FunctorType,
1874  result_view_type>::execute("", policy, functor,
1875  result_view);
1876 }
1877 
1878 template <class FunctorType>
1879 inline void parallel_reduce(const size_t& policy, const FunctorType& functor) {
1880  using policy_type =
1881  typename Impl::ParallelReducePolicyType<void, size_t,
1882  FunctorType>::policy_type;
1883  using FunctorAnalysis =
1884  Impl::FunctorAnalysis<Impl::FunctorPatternInterface::REDUCE, policy_type,
1885  FunctorType, void>;
1886  using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0),
1887  typename FunctorAnalysis::value_type,
1888  typename FunctorAnalysis::pointer_type>;
1889 
1890  static_assert(
1891  FunctorAnalysis::has_final_member_function,
1892  "Calling parallel_reduce without either return value or final function.");
1893 
1894  using result_view_type =
1896  result_view_type result_view;
1897 
1898  Impl::ParallelReduceAdaptor<policy_type, FunctorType,
1899  result_view_type>::execute("",
1900  policy_type(0, policy),
1901  functor, result_view);
1902 }
1903 
1904 template <class FunctorType>
1905 inline void parallel_reduce(const std::string& label, const size_t& policy,
1906  const FunctorType& functor) {
1907  using policy_type =
1908  typename Impl::ParallelReducePolicyType<void, size_t,
1909  FunctorType>::policy_type;
1910  using FunctorAnalysis =
1911  Impl::FunctorAnalysis<Impl::FunctorPatternInterface::REDUCE, policy_type,
1912  FunctorType, void>;
1913  using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0),
1914  typename FunctorAnalysis::value_type,
1915  typename FunctorAnalysis::pointer_type>;
1916 
1917  static_assert(
1918  FunctorAnalysis::has_final_member_function,
1919  "Calling parallel_reduce without either return value or final function.");
1920 
1921  using result_view_type =
1923  result_view_type result_view;
1924 
1925  Impl::ParallelReduceAdaptor<policy_type, FunctorType,
1926  result_view_type>::execute(label,
1927  policy_type(0, policy),
1928  functor, result_view);
1929 }
1930 
1931 } // namespace Kokkos
1932 
1933 #endif // KOKKOS_PARALLEL_REDUCE_HPP
Memory management for host memory.
ReturnType
Execution policy for work over a range of an integral type.