//@HEADER // ************************************************************************ // // Kokkos v. 4.0 // Copyright (2022) National Technology & Engineering // Solutions of Sandia, LLC (NTESS). // // Under the terms of Contract DE-NA0003525 with NTESS, // the U.S. Government retains certain rights in this software. // // Part of Kokkos, under the Apache License v2.0 with LLVM Exceptions. // See https://kokkos.org/LICENSE for license information. // SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception // //@HEADER #ifndef KOKKOS_IMPL_PUBLIC_INCLUDE #include static_assert(false, "Including non-public Kokkos header files is not allowed."); #endif #ifndef KOKKOS_PARALLEL_REDUCE_HPP #define KOKKOS_PARALLEL_REDUCE_HPP #include #include #include #include #include namespace Kokkos { template struct Sum { public: // Required using reducer = Sum; using value_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION Sum(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION Sum(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest += src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val = reduction_identity::sum(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE Sum(View const&) ->Sum::memory_space>; template struct Prod { public: // Required using reducer = Prod; using value_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION Prod(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION Prod(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest *= src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val = reduction_identity::prod(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE Prod(View const&) ->Prod::memory_space>; template struct Min { public: // Required using reducer = Min; using value_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION Min(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION Min(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (src < dest) dest = src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE Min(View const&) ->Min::memory_space>; template struct Max { public: // Required using reducer = Max; using value_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION Max(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION Max(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (src > dest) dest = src; } // Required KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val = reduction_identity::max(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE Max(View const&) ->Max::memory_space>; template struct LAnd { public: // Required using reducer = LAnd; using value_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION LAnd(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION LAnd(const result_view_type& value_) : value(value_), references_scalar_v(false) {} KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest = dest && src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val = reduction_identity::land(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE LAnd(View const&) ->LAnd::memory_space>; template struct LOr { public: // Required using reducer = LOr; using value_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION LOr(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION LOr(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest = dest || src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val = reduction_identity::lor(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE LOr(View const&) ->LOr::memory_space>; template struct BAnd { public: // Required using reducer = BAnd; using value_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION BAnd(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION BAnd(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest = dest & src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val = reduction_identity::band(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE BAnd(View const&) ->BAnd::memory_space>; template struct BOr { public: // Required using reducer = BOr; using value_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION BOr(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION BOr(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest = dest | src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val = reduction_identity::bor(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE BOr(View const&) ->BOr::memory_space>; template struct ValLocScalar { Scalar val; Index loc; }; template struct MinLoc { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); public: // Required using reducer = MinLoc; using value_type = ValLocScalar; using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION MinLoc(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION MinLoc(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (src.val < dest.val) dest = src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.val = reduction_identity::min(); val.loc = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MinLoc( View, Properties...> const&) ->MinLoc, Properties...>::memory_space>; template struct MaxLoc { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); public: // Required using reducer = MaxLoc; using value_type = ValLocScalar; using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION MaxLoc(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION MaxLoc(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (src.val > dest.val) dest = src; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.val = reduction_identity::max(); val.loc = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MaxLoc( View, Properties...> const&) ->MaxLoc, Properties...>::memory_space>; template struct MinMaxScalar { Scalar min_val, max_val; }; template struct MinMax { private: using scalar_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); public: // Required using reducer = MinMax; using value_type = MinMaxScalar; using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION MinMax(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION MinMax(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (src.min_val < dest.min_val) { dest.min_val = src.min_val; } if (src.max_val > dest.max_val) { dest.max_val = src.max_val; } } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.max_val = reduction_identity::max(); val.min_val = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MinMax(View, Properties...> const&) ->MinMax, Properties...>::memory_space>; template struct MinMaxLocScalar { Scalar min_val, max_val; Index min_loc, max_loc; }; template struct MinMaxLoc { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); public: // Required using reducer = MinMaxLoc; using value_type = MinMaxLocScalar; using result_view_type = Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION MinMaxLoc(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION MinMaxLoc(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (src.min_val < dest.min_val) { dest.min_val = src.min_val; dest.min_loc = src.min_loc; } if (src.max_val > dest.max_val) { dest.max_val = src.max_val; dest.max_loc = src.max_loc; } } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.max_val = reduction_identity::max(); val.min_val = reduction_identity::min(); val.max_loc = reduction_identity::min(); val.min_loc = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MinMaxLoc( View, Properties...> const&) ->MinMaxLoc, Properties...>::memory_space>; // -------------------------------------------------- // reducers added to support std algorithms // -------------------------------------------------- // // MaxFirstLoc // template struct MaxFirstLoc { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); static_assert(std::is_integral_v); public: // Required using reducer = MaxFirstLoc; using value_type = ::Kokkos::ValLocScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION MaxFirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION MaxFirstLoc(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (dest.val < src.val) { dest = src; } else if (!(src.val < dest.val)) { dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc; } } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.val = reduction_identity::max(); val.loc = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MaxFirstLoc( View, Properties...> const&) ->MaxFirstLoc, Properties...>::memory_space>; // // MaxFirstLocCustomComparator // recall that comp(a,b) returns true is a < b // template struct MaxFirstLocCustomComparator { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); static_assert(std::is_integral_v); public: // Required using reducer = MaxFirstLocCustomComparator; using value_type = ::Kokkos::ValLocScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; ComparatorType m_comp; public: KOKKOS_INLINE_FUNCTION MaxFirstLocCustomComparator(value_type& value_, ComparatorType comp_) : value(&value_), references_scalar_v(true), m_comp(comp_) {} KOKKOS_INLINE_FUNCTION MaxFirstLocCustomComparator(const result_view_type& value_, ComparatorType comp_) : value(value_), references_scalar_v(false), m_comp(comp_) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (m_comp(dest.val, src.val)) { dest = src; } else if (!m_comp(src.val, dest.val)) { dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc; } } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.val = reduction_identity::max(); val.loc = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MaxFirstLocCustomComparator( View, Properties...> const&, ComparatorType) ->MaxFirstLocCustomComparator, Properties...>::memory_space>; // // MinFirstLoc // template struct MinFirstLoc { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); static_assert(std::is_integral_v); public: // Required using reducer = MinFirstLoc; using value_type = ::Kokkos::ValLocScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION MinFirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION MinFirstLoc(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (src.val < dest.val) { dest = src; } else if (!(dest.val < src.val)) { dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc; } } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.val = reduction_identity::min(); val.loc = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MinFirstLoc( View, Properties...> const&) ->MinFirstLoc, Properties...>::memory_space>; // // MinFirstLocCustomComparator // recall that comp(a,b) returns true is a < b // template struct MinFirstLocCustomComparator { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); static_assert(std::is_integral_v); public: // Required using reducer = MinFirstLocCustomComparator; using value_type = ::Kokkos::ValLocScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; ComparatorType m_comp; public: KOKKOS_INLINE_FUNCTION MinFirstLocCustomComparator(value_type& value_, ComparatorType comp_) : value(&value_), references_scalar_v(true), m_comp(comp_) {} KOKKOS_INLINE_FUNCTION MinFirstLocCustomComparator(const result_view_type& value_, ComparatorType comp_) : value(value_), references_scalar_v(false), m_comp(comp_) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (m_comp(src.val, dest.val)) { dest = src; } else if (!m_comp(dest.val, src.val)) { dest.loc = (src.loc < dest.loc) ? src.loc : dest.loc; } } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.val = reduction_identity::min(); val.loc = reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MinFirstLocCustomComparator( View, Properties...> const&, ComparatorType) ->MinFirstLocCustomComparator, Properties...>::memory_space>; // // MinMaxFirstLastLoc // template struct MinMaxFirstLastLoc { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); static_assert(std::is_integral_v); public: // Required using reducer = MinMaxFirstLastLoc; using value_type = ::Kokkos::MinMaxLocScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION MinMaxFirstLastLoc(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION MinMaxFirstLastLoc(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (src.min_val < dest.min_val) { dest.min_val = src.min_val; dest.min_loc = src.min_loc; } else if (!(dest.min_val < src.min_val)) { dest.min_loc = (src.min_loc < dest.min_loc) ? src.min_loc : dest.min_loc; } if (dest.max_val < src.max_val) { dest.max_val = src.max_val; dest.max_loc = src.max_loc; } else if (!(src.max_val < dest.max_val)) { dest.max_loc = (src.max_loc > dest.max_loc) ? src.max_loc : dest.max_loc; } } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.max_val = ::Kokkos::reduction_identity::max(); val.min_val = ::Kokkos::reduction_identity::min(); val.max_loc = ::Kokkos::reduction_identity::max(); val.min_loc = ::Kokkos::reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MinMaxFirstLastLoc( View, Properties...> const&) ->MinMaxFirstLastLoc, Properties...>::memory_space>; // // MinMaxFirstLastLocCustomComparator // recall that comp(a,b) returns true is a < b // template struct MinMaxFirstLastLocCustomComparator { private: using scalar_type = std::remove_cv_t; using index_type = std::remove_cv_t; static_assert(!std::is_pointer_v && !std::is_array_v); static_assert(std::is_integral_v); public: // Required using reducer = MinMaxFirstLastLocCustomComparator; using value_type = ::Kokkos::MinMaxLocScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; ComparatorType m_comp; public: KOKKOS_INLINE_FUNCTION MinMaxFirstLastLocCustomComparator(value_type& value_, ComparatorType comp_) : value(&value_), references_scalar_v(true), m_comp(comp_) {} KOKKOS_INLINE_FUNCTION MinMaxFirstLastLocCustomComparator(const result_view_type& value_, ComparatorType comp_) : value(value_), references_scalar_v(false), m_comp(comp_) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { if (m_comp(src.min_val, dest.min_val)) { dest.min_val = src.min_val; dest.min_loc = src.min_loc; } else if (!m_comp(dest.min_val, src.min_val)) { dest.min_loc = (src.min_loc < dest.min_loc) ? src.min_loc : dest.min_loc; } if (m_comp(dest.max_val, src.max_val)) { dest.max_val = src.max_val; dest.max_loc = src.max_loc; } else if (!m_comp(src.max_val, dest.max_val)) { dest.max_loc = (src.max_loc > dest.max_loc) ? src.max_loc : dest.max_loc; } } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.max_val = ::Kokkos::reduction_identity::max(); val.min_val = ::Kokkos::reduction_identity::min(); val.max_loc = ::Kokkos::reduction_identity::max(); val.min_loc = ::Kokkos::reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE MinMaxFirstLastLocCustomComparator( View, Properties...> const&, ComparatorType) ->MinMaxFirstLastLocCustomComparator< Scalar, Index, ComparatorType, typename View, Properties...>::memory_space>; // // FirstLoc // template struct FirstLocScalar { Index min_loc_true; }; template struct FirstLoc { private: using index_type = std::remove_cv_t; static_assert(std::is_integral_v); public: // Required using reducer = FirstLoc; using value_type = FirstLocScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION FirstLoc(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION FirstLoc(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest.min_loc_true = (src.min_loc_true < dest.min_loc_true) ? src.min_loc_true : dest.min_loc_true; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.min_loc_true = ::Kokkos::reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE FirstLoc( View, Properties...> const&) ->FirstLoc, Properties...>::memory_space>; // // LastLoc // template struct LastLocScalar { Index max_loc_true; }; template struct LastLoc { private: using index_type = std::remove_cv_t; static_assert(std::is_integral_v); public: // Required using reducer = LastLoc; using value_type = LastLocScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION LastLoc(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION LastLoc(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest.max_loc_true = (src.max_loc_true > dest.max_loc_true) ? src.max_loc_true : dest.max_loc_true; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.max_loc_true = ::Kokkos::reduction_identity::max(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE LastLoc(View, Properties...> const&) ->LastLoc, Properties...>::memory_space>; template struct StdIsPartScalar { Index max_loc_true, min_loc_false; }; // // StdIsPartitioned // template struct StdIsPartitioned { private: using index_type = std::remove_cv_t; static_assert(std::is_integral_v); public: // Required using reducer = StdIsPartitioned; using value_type = StdIsPartScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION StdIsPartitioned(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION StdIsPartitioned(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest.max_loc_true = (dest.max_loc_true < src.max_loc_true) ? src.max_loc_true : dest.max_loc_true; dest.min_loc_false = (dest.min_loc_false < src.min_loc_false) ? dest.min_loc_false : src.min_loc_false; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.max_loc_true = ::Kokkos::reduction_identity::max(); val.min_loc_false = ::Kokkos::reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE StdIsPartitioned( View, Properties...> const&) ->StdIsPartitioned, Properties...>::memory_space>; template struct StdPartPointScalar { Index min_loc_false; }; // // StdPartitionPoint // template struct StdPartitionPoint { private: using index_type = std::remove_cv_t; static_assert(std::is_integral_v); public: // Required using reducer = StdPartitionPoint; using value_type = StdPartPointScalar; using result_view_type = ::Kokkos::View; private: result_view_type value; bool references_scalar_v; public: KOKKOS_INLINE_FUNCTION StdPartitionPoint(value_type& value_) : value(&value_), references_scalar_v(true) {} KOKKOS_INLINE_FUNCTION StdPartitionPoint(const result_view_type& value_) : value(value_), references_scalar_v(false) {} // Required KOKKOS_INLINE_FUNCTION void join(value_type& dest, const value_type& src) const { dest.min_loc_false = (dest.min_loc_false < src.min_loc_false) ? dest.min_loc_false : src.min_loc_false; } KOKKOS_INLINE_FUNCTION void init(value_type& val) const { val.min_loc_false = ::Kokkos::reduction_identity::min(); } KOKKOS_INLINE_FUNCTION value_type& reference() const { return *value.data(); } KOKKOS_INLINE_FUNCTION result_view_type view() const { return value; } KOKKOS_INLINE_FUNCTION bool references_scalar() const { return references_scalar_v; } }; template KOKKOS_DEDUCTION_GUIDE StdPartitionPoint( View, Properties...> const&) ->StdPartitionPoint, Properties...>::memory_space>; } // namespace Kokkos namespace Kokkos { namespace Impl { template class CombinedFunctorReducer { public: using functor_type = FunctorType; using reducer_type = FunctorAnalysisReducerType; CombinedFunctorReducer(const FunctorType& functor, const FunctorAnalysisReducerType& reducer) : m_functor(functor), m_reducer(reducer) {} KOKKOS_FUNCTION const FunctorType& get_functor() const { return m_functor; } KOKKOS_FUNCTION const FunctorAnalysisReducerType& get_reducer() const { return m_reducer; } private: FunctorType m_functor; FunctorAnalysisReducerType m_reducer; }; template class CombinedFunctorReducer< FunctorType, FunctorAnalysisReducerType, std::enable_if_t>> { public: using functor_type = FunctorType; using reducer_type = FunctorAnalysisReducerType; CombinedFunctorReducer(const FunctorType& functor, const FunctorAnalysisReducerType&) : m_reducer(functor) {} KOKKOS_FUNCTION const FunctorType& get_functor() const { return m_reducer.get_functor(); } KOKKOS_FUNCTION const FunctorAnalysisReducerType& get_reducer() const { return m_reducer; } private: FunctorAnalysisReducerType m_reducer; }; template struct ParallelReduceReturnValue; template struct ParallelReduceReturnValue< std::enable_if_t::value>, ReturnType, FunctorType> { using return_type = ReturnType; using reducer_type = InvalidType; using value_type_scalar = typename return_type::value_type; using value_type_array = typename return_type::value_type* const; using value_type = std::conditional_t; static return_type& return_value(ReturnType& return_val, const FunctorType&) { return return_val; } }; template struct ParallelReduceReturnValue< std::enable_if_t::value && (!std::is_array::value && !std::is_pointer::value) && !Kokkos::is_reducer::value>, ReturnType, FunctorType> { using return_type = Kokkos::View; using reducer_type = InvalidType; using value_type = typename return_type::value_type; static return_type return_value(ReturnType& return_val, const FunctorType&) { return return_type(&return_val); } }; template struct ParallelReduceReturnValue< std::enable_if_t<(std::is_array::value || std::is_pointer::value)>, ReturnType, FunctorType> { using return_type = Kokkos::View, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>; using reducer_type = InvalidType; using value_type = typename return_type::value_type[]; static return_type return_value(ReturnType& return_val, const FunctorType& functor) { if (std::is_array::value) return return_type(return_val); else return return_type(return_val, functor.value_count); } }; template struct ParallelReduceReturnValue< std::enable_if_t::value>, ReturnType, FunctorType> { using return_type = typename ReturnType::result_view_type; using reducer_type = ReturnType; using value_type = typename return_type::value_type; static auto return_value(ReturnType& return_val, const FunctorType&) { return return_val.view(); } }; template struct ParallelReducePolicyType; template struct ParallelReducePolicyType< std::enable_if_t::value>, PolicyType, FunctorType> { using policy_type = PolicyType; static PolicyType policy(const PolicyType& policy_) { return policy_; } }; template struct ParallelReducePolicyType< std::enable_if_t::value>, PolicyType, FunctorType> { using execution_space = typename Impl::FunctorPolicyExecutionSpace::execution_space; using policy_type = Kokkos::RangePolicy; static policy_type policy(const PolicyType& policy_) { return policy_type(0, policy_); } }; template struct ParallelReduceFunctorType { using functor_type = FunctorType; static const functor_type& functor(const functor_type& functor) { return functor; } }; template struct ParallelReduceAdaptor { using return_value_adapter = Impl::ParallelReduceReturnValue; static inline void execute_impl(const std::string& label, const PolicyType& policy, const FunctorType& functor, ReturnType& return_value) { using PassedReducerType = typename return_value_adapter::reducer_type; uint64_t kpID = 0; PolicyType inner_policy = policy; Kokkos::Tools::Impl::begin_parallel_reduce( inner_policy, functor, label, kpID); using ReducerSelector = Kokkos::Impl::if_c::value, FunctorType, PassedReducerType>; using Analysis = FunctorAnalysis; using CombinedFunctorReducerType = CombinedFunctorReducer; auto closure = construct_with_shared_allocation_tracking_disabled< Impl::ParallelReduce::execution_space>>( CombinedFunctorReducerType( functor, typename Analysis::Reducer( ReducerSelector::select(functor, return_value))), inner_policy, return_value_adapter::return_value(return_value, functor)); closure.execute(); Kokkos::Tools::Impl::end_parallel_reduce( inner_policy, functor, label, kpID); } static constexpr bool is_array_reduction = Impl::FunctorAnalysis< Impl::FunctorPatternInterface::REDUCE, PolicyType, FunctorType, typename return_value_adapter::value_type>::StaticValueSize == 0; template static inline std::enable_if_t::value)> execute(const std::string& label, const PolicyType& policy, const FunctorType& functor, ReturnType& return_value) { execute_impl(label, policy, functor, return_value); } }; } // namespace Impl //---------------------------------------------------------------------------- /*! \fn void parallel_reduce(label,policy,functor,return_argument) \brief Perform a parallel reduction. \param label An optional Label giving the call name. Must be able to construct a std::string from the argument. \param policy A Kokkos Execution Policy, such as an integer, a RangePolicy or a TeamPolicy. \param functor A functor with a reduction operator, and optional init, join and final functions. \param return_argument A return argument which can be a scalar, a View, or a ReducerStruct. This argument can be left out if the functor has a final function. */ // Parallel Reduce Blocking behavior namespace Impl { template struct ReducerHasTestReferenceFunction { template static std::true_type test_func(decltype(&E::references_scalar)); template static std::false_type test_func(...); enum { value = std::is_same(nullptr))>::value }; }; template constexpr std::enable_if_t< // constraints only necessary because SFINAE lacks subsumption !ReducerHasTestReferenceFunction::value && !Kokkos::is_view::value, // return type: bool> parallel_reduce_needs_fence(ExecutionSpace const&, Arg const&) { return true; } template constexpr std::enable_if_t< // equivalent to: // (requires (Reducer const& r) { // { reducer.references_scalar() } -> std::convertible_to; // }) ReducerHasTestReferenceFunction::value, // return type: bool> parallel_reduce_needs_fence(ExecutionSpace const&, Reducer const& reducer) { return reducer.references_scalar(); } template constexpr std::enable_if_t< // requires Kokkos::ViewLike Kokkos::is_view::value, // return type: bool> parallel_reduce_needs_fence(ExecutionSpace const&, ViewLike const&) { return false; } template struct ParallelReduceFence { template static void fence(const ExecutionSpace& ex, const std::string& name, ArgsDeduced&&... args) { if (Impl::parallel_reduce_needs_fence(ex, (ArgsDeduced &&) args...)) { ex.fence(name); } } }; } // namespace Impl /** \brief Parallel reduction * * parallel_reduce performs parallel reductions with arbitrary functions - i.e. * it is not solely data based. The call expects up to 4 arguments: * * * Example of a parallel_reduce functor for a POD (plain old data) value type: * \code * class FunctorType { // For POD value type * public: * using execution_space = ...; * using value_type = ; * void operator()( iwork , & update ) const ; * void init( & update ) const ; * void join( & update , * const & input ) const ; * * void final( & update ) const ; * }; * \endcode * * Example of a parallel_reduce functor for an array of POD (plain old data) * values: * \code * class FunctorType { // For array of POD value * public: * using execution_space = ...; * using value_type = []; * void operator()( , update[] ) const ; * void init( update[] ) const ; * void join( update[] , * const input[] ) const ; * * void final( update[] ) const ; * }; * \endcode */ // ReturnValue is scalar or array: take by reference template inline std::enable_if_t::value && !(Kokkos::is_view::value || Kokkos::is_reducer::value || std::is_pointer::value)> parallel_reduce(const std::string& label, const PolicyType& policy, const FunctorType& functor, ReturnType& return_value) { static_assert( !std::is_const::value, "A const reduction result type is only allowed for a View, pointer or " "reducer return type!"); Impl::ParallelReduceAdaptor::execute( label, policy, functor, return_value); Impl::ParallelReduceFence:: fence( policy.space(), "Kokkos::parallel_reduce: fence due to result being value, not view", return_value); } template inline std::enable_if_t::value && !(Kokkos::is_view::value || Kokkos::is_reducer::value || std::is_pointer::value)> parallel_reduce(const PolicyType& policy, const FunctorType& functor, ReturnType& return_value) { static_assert( !std::is_const::value, "A const reduction result type is only allowed for a View, pointer or " "reducer return type!"); Impl::ParallelReduceAdaptor::execute( "", policy, functor, return_value); Impl::ParallelReduceFence:: fence( policy.space(), "Kokkos::parallel_reduce: fence due to result being value, not view", return_value); } template inline std::enable_if_t::value || Kokkos::is_reducer::value || std::is_pointer::value)> parallel_reduce(const size_t& policy, const FunctorType& functor, ReturnType& return_value) { static_assert( !std::is_const::value, "A const reduction result type is only allowed for a View, pointer or " "reducer return type!"); using policy_type = typename Impl::ParallelReducePolicyType::policy_type; Impl::ParallelReduceAdaptor::execute( "", policy_type(0, policy), functor, return_value); Impl::ParallelReduceFence:: fence( typename policy_type::execution_space(), "Kokkos::parallel_reduce: fence due to result being value, not view", return_value); } template inline std::enable_if_t::value || Kokkos::is_reducer::value || std::is_pointer::value)> parallel_reduce(const std::string& label, const size_t& policy, const FunctorType& functor, ReturnType& return_value) { static_assert( !std::is_const::value, "A const reduction result type is only allowed for a View, pointer or " "reducer return type!"); using policy_type = typename Impl::ParallelReducePolicyType::policy_type; Impl::ParallelReduceAdaptor::execute( label, policy_type(0, policy), functor, return_value); Impl::ParallelReduceFence:: fence( typename policy_type::execution_space(), "Kokkos::parallel_reduce: fence due to result being value, not view", return_value); } // ReturnValue as View or Reducer: take by copy to allow for inline construction template inline std::enable_if_t::value && (Kokkos::is_view::value || Kokkos::is_reducer::value || std::is_pointer::value)> parallel_reduce(const std::string& label, const PolicyType& policy, const FunctorType& functor, const ReturnType& return_value) { ReturnType return_value_impl = return_value; Impl::ParallelReduceAdaptor::execute( label, policy, functor, return_value_impl); Impl::ParallelReduceFence:: fence( policy.space(), "Kokkos::parallel_reduce: fence due to result being value, not view", return_value); } template inline std::enable_if_t::value && (Kokkos::is_view::value || Kokkos::is_reducer::value || std::is_pointer::value)> parallel_reduce(const PolicyType& policy, const FunctorType& functor, const ReturnType& return_value) { ReturnType return_value_impl = return_value; Impl::ParallelReduceAdaptor::execute( "", policy, functor, return_value_impl); Impl::ParallelReduceFence:: fence( policy.space(), "Kokkos::parallel_reduce: fence due to result being value, not view", return_value); } template inline std::enable_if_t::value || Kokkos::is_reducer::value || std::is_pointer::value> parallel_reduce(const size_t& policy, const FunctorType& functor, const ReturnType& return_value) { using policy_type = typename Impl::ParallelReducePolicyType::policy_type; ReturnType return_value_impl = return_value; Impl::ParallelReduceAdaptor::execute( "", policy_type(0, policy), functor, return_value_impl); Impl::ParallelReduceFence:: fence( typename policy_type::execution_space(), "Kokkos::parallel_reduce: fence due to result being value, not view", return_value); } template inline std::enable_if_t::value || Kokkos::is_reducer::value || std::is_pointer::value> parallel_reduce(const std::string& label, const size_t& policy, const FunctorType& functor, const ReturnType& return_value) { using policy_type = typename Impl::ParallelReducePolicyType::policy_type; ReturnType return_value_impl = return_value; Impl::ParallelReduceAdaptor::execute( label, policy_type(0, policy), functor, return_value_impl); Impl::ParallelReduceFence:: fence( typename policy_type::execution_space(), "Kokkos::parallel_reduce: fence due to result being value, not view", return_value); } // No Return Argument template inline void parallel_reduce( const std::string& label, const PolicyType& policy, const FunctorType& functor, std::enable_if_t::value>* = nullptr) { using FunctorAnalysis = Impl::FunctorAnalysis; using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0), typename FunctorAnalysis::value_type, typename FunctorAnalysis::pointer_type>; static_assert( FunctorAnalysis::has_final_member_function, "Calling parallel_reduce without either return value or final function."); using result_view_type = Kokkos::View; result_view_type result_view; Impl::ParallelReduceAdaptor::execute(label, policy, functor, result_view); } template inline void parallel_reduce( const PolicyType& policy, const FunctorType& functor, std::enable_if_t::value>* = nullptr) { using FunctorAnalysis = Impl::FunctorAnalysis; using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0), typename FunctorAnalysis::value_type, typename FunctorAnalysis::pointer_type>; static_assert( FunctorAnalysis::has_final_member_function, "Calling parallel_reduce without either return value or final function."); using result_view_type = Kokkos::View; result_view_type result_view; Impl::ParallelReduceAdaptor::execute("", policy, functor, result_view); } template inline void parallel_reduce(const size_t& policy, const FunctorType& functor) { using policy_type = typename Impl::ParallelReducePolicyType::policy_type; using FunctorAnalysis = Impl::FunctorAnalysis; using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0), typename FunctorAnalysis::value_type, typename FunctorAnalysis::pointer_type>; static_assert( FunctorAnalysis::has_final_member_function, "Calling parallel_reduce without either return value or final function."); using result_view_type = Kokkos::View; result_view_type result_view; Impl::ParallelReduceAdaptor::execute("", policy_type(0, policy), functor, result_view); } template inline void parallel_reduce(const std::string& label, const size_t& policy, const FunctorType& functor) { using policy_type = typename Impl::ParallelReducePolicyType::policy_type; using FunctorAnalysis = Impl::FunctorAnalysis; using value_type = std::conditional_t<(FunctorAnalysis::StaticValueSize != 0), typename FunctorAnalysis::value_type, typename FunctorAnalysis::pointer_type>; static_assert( FunctorAnalysis::has_final_member_function, "Calling parallel_reduce without either return value or final function."); using result_view_type = Kokkos::View; result_view_type result_view; Impl::ParallelReduceAdaptor::execute(label, policy_type(0, policy), functor, result_view); } } // namespace Kokkos #endif // KOKKOS_PARALLEL_REDUCE_HPP