48 #ifndef KOKKOS_SERIAL_HPP
49 #define KOKKOS_SERIAL_HPP
51 #include <Kokkos_Macros.hpp>
52 #if defined(KOKKOS_ENABLE_SERIAL)
57 #include <Kokkos_TaskScheduler.hpp>
59 #include <Kokkos_HostSpace.hpp>
60 #include <Kokkos_ScratchSpace.hpp>
61 #include <Kokkos_MemoryTraits.hpp>
62 #include <impl/Kokkos_Tags.hpp>
63 #include <impl/Kokkos_HostThreadTeam.hpp>
64 #include <impl/Kokkos_FunctorAnalysis.hpp>
65 #include <impl/Kokkos_FunctorAdapter.hpp>
66 #include <impl/Kokkos_Profiling_Interface.hpp>
68 #include <KokkosExp_MDRangePolicy.hpp>
70 #include <Kokkos_UniqueToken.hpp>
92 typedef Serial execution_space;
94 typedef HostSpace::size_type size_type;
96 typedef HostSpace memory_space;
98 typedef Kokkos::Device<execution_space, memory_space> device_type;
101 typedef LayoutRight array_layout;
104 typedef ScratchMemorySpace<Kokkos::Serial> scratch_memory_space;
114 inline static int in_parallel() {
return false; }
122 static void impl_static_fence() {}
124 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
125 static void fence() {}
127 void fence()
const {}
131 static int concurrency() {
return 1; }
135 const bool =
false) {}
137 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
141 static void initialize(
unsigned threads_count = 1,
142 unsigned use_numa_count = 0,
143 unsigned use_cores_per_numa = 0,
144 bool allow_asynchronous_threadpool =
false);
146 static bool is_initialized();
153 inline static int thread_pool_size(
int = 0) {
return 1; }
154 KOKKOS_INLINE_FUNCTION
static int thread_pool_rank() {
return 0; }
158 KOKKOS_INLINE_FUNCTION
static unsigned hardware_thread_id() {
159 return thread_pool_rank();
161 inline static unsigned max_hardware_threads() {
return thread_pool_size(0); }
163 static void impl_initialize();
165 static bool impl_is_initialized();
168 static void impl_finalize();
172 inline static int impl_thread_pool_size(
int = 0) {
return 1; }
173 KOKKOS_INLINE_FUNCTION
static int impl_thread_pool_rank() {
return 0; }
177 KOKKOS_INLINE_FUNCTION
static unsigned impl_hardware_thread_id() {
178 return impl_thread_pool_rank();
180 inline static unsigned impl_max_hardware_threads() {
181 return impl_thread_pool_size(0);
184 uint32_t impl_instance_id() const noexcept {
return 0; }
186 static const char* name();
190 namespace Profiling {
191 namespace Experimental {
193 struct DeviceTypeTraits<Serial> {
194 static constexpr DeviceType
id = DeviceType::Serial;
207 struct MemorySpaceAccess<Kokkos::Serial::memory_space,
208 Kokkos::Serial::scratch_memory_space> {
209 enum { assignable =
false };
210 enum { accessible =
true };
211 enum { deepcopy =
false };
215 struct VerifyExecutionCanAccessMemorySpace<
216 Kokkos::Serial::memory_space, Kokkos::Serial::scratch_memory_space> {
217 enum { value =
true };
218 inline static void verify(
void) {}
219 inline static void verify(
const void*) {}
232 void serial_resize_thread_team_data(
size_t pool_reduce_bytes,
233 size_t team_reduce_bytes,
234 size_t team_shared_bytes,
235 size_t thread_local_bytes);
237 HostThreadTeamData* serial_get_thread_team_data();
251 template <
class... Properties>
252 class TeamPolicyInternal<Kokkos::Serial, Properties...>
253 :
public PolicyTraits<Properties...> {
255 size_t m_team_scratch_size[2];
256 size_t m_thread_scratch_size[2];
262 typedef TeamPolicyInternal execution_policy;
264 typedef PolicyTraits<Properties...> traits;
267 typedef Kokkos::Serial execution_space;
269 const typename traits::execution_space& space()
const {
270 static typename traits::execution_space m_space;
274 TeamPolicyInternal& operator=(
const TeamPolicyInternal& p) {
275 m_league_size = p.m_league_size;
276 m_team_scratch_size[0] = p.m_team_scratch_size[0];
277 m_thread_scratch_size[0] = p.m_thread_scratch_size[0];
278 m_team_scratch_size[1] = p.m_team_scratch_size[1];
279 m_thread_scratch_size[1] = p.m_thread_scratch_size[1];
280 m_chunk_size = p.m_chunk_size;
284 template <
class ExecSpace,
class... OtherProperties>
285 friend class TeamPolicyInternal;
287 template <
class... OtherProperties>
289 const TeamPolicyInternal<Kokkos::Serial, OtherProperties...>& p) {
290 m_league_size = p.m_league_size;
291 m_team_scratch_size[0] = p.m_team_scratch_size[0];
292 m_thread_scratch_size[0] = p.m_thread_scratch_size[0];
293 m_team_scratch_size[1] = p.m_team_scratch_size[1];
294 m_thread_scratch_size[1] = p.m_thread_scratch_size[1];
295 m_chunk_size = p.m_chunk_size;
299 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
300 template <
class FunctorType>
301 static int team_size_max(
const FunctorType&) {
305 template <
class FunctorType>
306 static int team_size_recommended(
const FunctorType&) {
310 template <
class FunctorType>
311 static int team_size_recommended(
const FunctorType&,
const int&) {
316 template <
class FunctorType>
317 int team_size_max(
const FunctorType&,
const ParallelForTag&)
const {
320 template <
class FunctorType>
321 int team_size_max(
const FunctorType&,
const ParallelReduceTag&)
const {
324 template <
class FunctorType,
class ReducerType>
325 int team_size_max(
const FunctorType&,
const ReducerType&,
326 const ParallelReduceTag&)
const {
329 template <
class FunctorType>
330 int team_size_recommended(
const FunctorType&,
const ParallelForTag&)
const {
333 template <
class FunctorType>
334 int team_size_recommended(
const FunctorType&,
335 const ParallelReduceTag&)
const {
338 template <
class FunctorType,
class ReducerType>
339 int team_size_recommended(
const FunctorType&,
const ReducerType&,
340 const ParallelReduceTag&)
const {
346 inline int team_size()
const {
return 1; }
347 inline int league_size()
const {
return m_league_size; }
348 inline size_t scratch_size(
const int& level,
int = 0)
const {
349 return m_team_scratch_size[level] + m_thread_scratch_size[level];
352 inline static int vector_length_max() {
356 inline static int scratch_size_max(
int level) {
357 return (level == 0 ? 1024 * 32 : 20 * 1024 * 1024);
360 TeamPolicyInternal(
const execution_space&,
int league_size_request
361 #ifndef KOKKOS_ENABLE_DEPRECATED_CODE
363 int team_size_request
370 : m_team_scratch_size{0, 0},
371 m_thread_scratch_size{0, 0},
372 m_league_size(league_size_request),
374 #ifndef KOKKOS_ENABLE_DEPRECATED_CODE
375 if (team_size_request > 1)
376 Kokkos::abort(
"Kokkos::abort: Requested Team Size is too large!");
380 TeamPolicyInternal(
const execution_space&,
int league_size_request,
381 const Kokkos::AUTO_t&
384 : m_team_scratch_size{0, 0},
385 m_thread_scratch_size{0, 0},
386 m_league_size(league_size_request),
389 TeamPolicyInternal(
int league_size_request
390 #ifndef KOKKOS_ENABLE_DEPRECATED_CODE
392 int team_size_request
399 : m_team_scratch_size{0, 0},
400 m_thread_scratch_size{0, 0},
401 m_league_size(league_size_request),
403 #ifndef KOKKOS_ENABLE_DEPRECATED_CODE
404 if (team_size_request > 1)
405 Kokkos::abort(
"Kokkos::abort: Requested Team Size is too large!");
409 TeamPolicyInternal(
int league_size_request,
410 const Kokkos::AUTO_t&
413 : m_team_scratch_size{0, 0},
414 m_thread_scratch_size{0, 0},
415 m_league_size(league_size_request),
418 inline int chunk_size()
const {
return m_chunk_size; }
420 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
422 inline TeamPolicyInternal set_chunk_size(
423 typename traits::index_type chunk_size_)
const {
424 TeamPolicyInternal p = *
this;
425 p.m_chunk_size = chunk_size_;
431 inline TeamPolicyInternal set_scratch_size(
432 const int& level,
const PerTeamValue& per_team)
const {
433 TeamPolicyInternal p = *
this;
434 p.m_team_scratch_size[level] = per_team.value;
440 inline TeamPolicyInternal set_scratch_size(
441 const int& level,
const PerThreadValue& per_thread)
const {
442 TeamPolicyInternal p = *
this;
443 p.m_thread_scratch_size[level] = per_thread.value;
449 inline TeamPolicyInternal set_scratch_size(
450 const int& level,
const PerTeamValue& per_team,
451 const PerThreadValue& per_thread)
const {
452 TeamPolicyInternal p = *
this;
453 p.m_team_scratch_size[level] = per_team.value;
454 p.m_thread_scratch_size[level] = per_thread.value;
459 inline TeamPolicyInternal& set_chunk_size(
460 typename traits::index_type chunk_size_) {
461 m_chunk_size = chunk_size_;
467 inline TeamPolicyInternal& set_scratch_size(
const int& level,
468 const PerTeamValue& per_team) {
469 m_team_scratch_size[level] = per_team.value;
475 inline TeamPolicyInternal& set_scratch_size(
476 const int& level,
const PerThreadValue& per_thread) {
477 m_thread_scratch_size[level] = per_thread.value;
483 inline TeamPolicyInternal& set_scratch_size(
484 const int& level,
const PerTeamValue& per_team,
485 const PerThreadValue& per_thread) {
486 m_team_scratch_size[level] = per_team.value;
487 m_thread_scratch_size[level] = per_thread.value;
492 typedef Impl::HostThreadTeamMember<Kokkos::Serial> member_type;
495 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
497 inline TeamPolicyInternal internal_set_chunk_size(
498 typename traits::index_type chunk_size_) {
499 m_chunk_size = chunk_size_;
505 inline TeamPolicyInternal internal_set_scratch_size(
506 const int& level,
const PerTeamValue& per_team) {
507 m_team_scratch_size[level] = per_team.value;
513 inline TeamPolicyInternal internal_set_scratch_size(
514 const int& level,
const PerThreadValue& per_thread) {
515 m_thread_scratch_size[level] = per_thread.value;
521 inline TeamPolicyInternal internal_set_scratch_size(
522 const int& level,
const PerTeamValue& per_team,
523 const PerThreadValue& per_thread) {
524 m_team_scratch_size[level] = per_team.value;
525 m_thread_scratch_size[level] = per_thread.value;
540 template <
class FunctorType,
class... Traits>
541 class ParallelFor<FunctorType, Kokkos::RangePolicy<Traits...>, Kokkos::Serial> {
545 const FunctorType m_functor;
546 const Policy m_policy;
548 template <
class TagType>
549 typename std::enable_if<std::is_same<TagType, void>::value>::type exec()
551 const typename Policy::member_type e = m_policy.end();
552 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
557 template <
class TagType>
558 typename std::enable_if<!std::is_same<TagType, void>::value>::type exec()
561 const typename Policy::member_type e = m_policy.end();
562 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
568 inline void execute()
const {
569 this->
template exec<typename Policy::work_tag>();
572 inline ParallelFor(
const FunctorType& arg_functor,
const Policy& arg_policy)
573 : m_functor(arg_functor), m_policy(arg_policy) {}
578 template <
class FunctorType,
class ReducerType,
class... Traits>
579 class ParallelReduce<FunctorType, Kokkos::RangePolicy<Traits...>, ReducerType,
583 typedef typename Policy::work_tag WorkTag;
585 typedef Kokkos::Impl::if_c<std::is_same<InvalidType, ReducerType>::value,
586 FunctorType, ReducerType>
589 typedef typename ReducerConditional::type ReducerTypeFwd;
591 typename Kokkos::Impl::if_c<std::is_same<InvalidType, ReducerType>::value,
592 WorkTag,
void>::type WorkTagFwd;
594 typedef FunctorAnalysis<FunctorPatternInterface::REDUCE, Policy, FunctorType>
597 typedef Kokkos::Impl::FunctorValueInit<ReducerTypeFwd, WorkTagFwd> ValueInit;
599 typedef typename Analysis::pointer_type pointer_type;
600 typedef typename Analysis::reference_type reference_type;
602 const FunctorType m_functor;
603 const Policy m_policy;
604 const ReducerType m_reducer;
605 const pointer_type m_result_ptr;
607 template <
class TagType>
608 inline typename std::enable_if<std::is_same<TagType, void>::value>::type exec(
609 reference_type update)
const {
610 const typename Policy::member_type e = m_policy.end();
611 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
612 m_functor(i, update);
616 template <
class TagType>
617 inline typename std::enable_if<!std::is_same<TagType, void>::value>::type
618 exec(reference_type update)
const {
621 const typename Policy::member_type e = m_policy.end();
622 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
623 m_functor(t, i, update);
628 inline void execute()
const {
629 const size_t pool_reduce_size =
630 Analysis::value_size(ReducerConditional::select(m_functor, m_reducer));
631 const size_t team_reduce_size = 0;
632 const size_t team_shared_size = 0;
633 const size_t thread_local_size = 0;
635 serial_resize_thread_team_data(pool_reduce_size, team_reduce_size,
636 team_shared_size, thread_local_size);
638 HostThreadTeamData& data = *serial_get_thread_team_data();
641 m_result_ptr ? m_result_ptr : pointer_type(data.pool_reduce_local());
643 reference_type update =
644 ValueInit::init(ReducerConditional::select(m_functor, m_reducer), ptr);
646 this->
template exec<WorkTag>(update);
648 Kokkos::Impl::FunctorFinal<ReducerTypeFwd, WorkTagFwd>::final(
649 ReducerConditional::select(m_functor, m_reducer), ptr);
652 template <
class HostViewType>
654 const FunctorType& arg_functor,
const Policy& arg_policy,
655 const HostViewType& arg_result_view,
656 typename std::enable_if<Kokkos::is_view<HostViewType>::value &&
657 !Kokkos::is_reducer_type<ReducerType>::value,
658 void*>::type =
nullptr)
659 : m_functor(arg_functor),
660 m_policy(arg_policy),
661 m_reducer(InvalidType()),
662 m_result_ptr(arg_result_view.data()) {
663 static_assert(Kokkos::is_view<HostViewType>::value,
664 "Kokkos::Serial reduce result must be a View");
667 std::is_same<typename HostViewType::memory_space, HostSpace>::value,
668 "Kokkos::Serial reduce result must be a View in HostSpace");
671 inline ParallelReduce(
const FunctorType& arg_functor, Policy arg_policy,
672 const ReducerType& reducer)
673 : m_functor(arg_functor),
674 m_policy(arg_policy),
676 m_result_ptr(reducer.view().data()) {
686 template <
class FunctorType,
class... Traits>
687 class ParallelScan<FunctorType, Kokkos::RangePolicy<Traits...>,
691 typedef typename Policy::work_tag WorkTag;
693 typedef FunctorAnalysis<FunctorPatternInterface::SCAN, Policy, FunctorType>
696 typedef Kokkos::Impl::FunctorValueInit<FunctorType, WorkTag> ValueInit;
698 typedef typename Analysis::pointer_type pointer_type;
699 typedef typename Analysis::reference_type reference_type;
701 const FunctorType m_functor;
702 const Policy m_policy;
704 template <
class TagType>
705 inline typename std::enable_if<std::is_same<TagType, void>::value>::type exec(
706 reference_type update)
const {
707 const typename Policy::member_type e = m_policy.end();
708 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
709 m_functor(i, update,
true);
713 template <
class TagType>
714 inline typename std::enable_if<!std::is_same<TagType, void>::value>::type
715 exec(reference_type update)
const {
717 const typename Policy::member_type e = m_policy.end();
718 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
719 m_functor(t, i, update,
true);
724 inline void execute()
const {
725 const size_t pool_reduce_size = Analysis::value_size(m_functor);
726 const size_t team_reduce_size = 0;
727 const size_t team_shared_size = 0;
728 const size_t thread_local_size = 0;
730 serial_resize_thread_team_data(pool_reduce_size, team_reduce_size,
731 team_shared_size, thread_local_size);
733 HostThreadTeamData& data = *serial_get_thread_team_data();
735 reference_type update =
736 ValueInit::init(m_functor, pointer_type(data.pool_reduce_local()));
738 this->
template exec<WorkTag>(update);
741 inline ParallelScan(
const FunctorType& arg_functor,
const Policy& arg_policy)
742 : m_functor(arg_functor), m_policy(arg_policy) {}
746 template <
class FunctorType,
class ReturnType,
class... Traits>
747 class ParallelScanWithTotal<FunctorType, Kokkos::RangePolicy<Traits...>,
751 typedef typename Policy::work_tag WorkTag;
753 typedef FunctorAnalysis<FunctorPatternInterface::SCAN, Policy, FunctorType>
756 typedef Kokkos::Impl::FunctorValueInit<FunctorType, WorkTag> ValueInit;
758 typedef typename Analysis::pointer_type pointer_type;
759 typedef typename Analysis::reference_type reference_type;
761 const FunctorType m_functor;
762 const Policy m_policy;
765 template <
class TagType>
766 inline typename std::enable_if<std::is_same<TagType, void>::value>::type exec(
767 reference_type update)
const {
768 const typename Policy::member_type e = m_policy.end();
769 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
770 m_functor(i, update,
true);
774 template <
class TagType>
775 inline typename std::enable_if<!std::is_same<TagType, void>::value>::type
776 exec(reference_type update)
const {
778 const typename Policy::member_type e = m_policy.end();
779 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
780 m_functor(t, i, update,
true);
785 inline void execute() {
786 const size_t pool_reduce_size = Analysis::value_size(m_functor);
787 const size_t team_reduce_size = 0;
788 const size_t team_shared_size = 0;
789 const size_t thread_local_size = 0;
791 serial_resize_thread_team_data(pool_reduce_size, team_reduce_size,
792 team_shared_size, thread_local_size);
794 HostThreadTeamData& data = *serial_get_thread_team_data();
796 reference_type update =
797 ValueInit::init(m_functor, pointer_type(data.pool_reduce_local()));
799 this->
template exec<WorkTag>(update);
801 m_returnvalue = update;
804 inline ParallelScanWithTotal(
const FunctorType& arg_functor,
805 const Policy& arg_policy,
807 : m_functor(arg_functor),
808 m_policy(arg_policy),
809 m_returnvalue(arg_returnvalue) {}
822 template <
class FunctorType,
class... Traits>
823 class ParallelFor<FunctorType, Kokkos::MDRangePolicy<Traits...>,
826 typedef Kokkos::MDRangePolicy<Traits...> MDRangePolicy;
827 typedef typename MDRangePolicy::impl_range_policy Policy;
829 typedef typename Kokkos::Impl::HostIterateTile<
830 MDRangePolicy, FunctorType,
typename MDRangePolicy::work_tag,
void>
833 const FunctorType m_functor;
834 const MDRangePolicy m_mdr_policy;
835 const Policy m_policy;
838 const typename Policy::member_type e = m_policy.end();
839 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
840 iterate_type(m_mdr_policy, m_functor)(i);
845 inline void execute()
const { this->exec(); }
847 inline ParallelFor(
const FunctorType& arg_functor,
848 const MDRangePolicy& arg_policy)
849 : m_functor(arg_functor),
850 m_mdr_policy(arg_policy),
851 m_policy(Policy(0, m_mdr_policy.m_num_tiles).set_chunk_size(1)) {}
854 template <
class FunctorType,
class ReducerType,
class... Traits>
855 class ParallelReduce<FunctorType, Kokkos::MDRangePolicy<Traits...>, ReducerType,
858 typedef Kokkos::MDRangePolicy<Traits...> MDRangePolicy;
859 typedef typename MDRangePolicy::impl_range_policy Policy;
861 typedef typename MDRangePolicy::work_tag WorkTag;
863 typedef Kokkos::Impl::if_c<std::is_same<InvalidType, ReducerType>::value,
864 FunctorType, ReducerType>
866 typedef typename ReducerConditional::type ReducerTypeFwd;
868 typename Kokkos::Impl::if_c<std::is_same<InvalidType, ReducerType>::value,
869 WorkTag,
void>::type WorkTagFwd;
871 typedef FunctorAnalysis<FunctorPatternInterface::REDUCE, MDRangePolicy,
875 typedef Kokkos::Impl::FunctorValueInit<ReducerTypeFwd, WorkTagFwd> ValueInit;
877 typedef typename Analysis::pointer_type pointer_type;
878 typedef typename Analysis::value_type value_type;
879 typedef typename Analysis::reference_type reference_type;
882 typename Kokkos::Impl::HostIterateTile<MDRangePolicy, FunctorType,
883 WorkTag, reference_type>;
885 const FunctorType m_functor;
886 const MDRangePolicy m_mdr_policy;
887 const Policy m_policy;
888 const ReducerType m_reducer;
889 const pointer_type m_result_ptr;
891 inline void exec(reference_type update)
const {
892 const typename Policy::member_type e = m_policy.end();
893 for (
typename Policy::member_type i = m_policy.begin(); i < e; ++i) {
894 iterate_type(m_mdr_policy, m_functor, update)(i);
899 inline void execute()
const {
900 const size_t pool_reduce_size =
901 Analysis::value_size(ReducerConditional::select(m_functor, m_reducer));
902 const size_t team_reduce_size = 0;
903 const size_t team_shared_size = 0;
904 const size_t thread_local_size = 0;
906 serial_resize_thread_team_data(pool_reduce_size, team_reduce_size,
907 team_shared_size, thread_local_size);
909 HostThreadTeamData& data = *serial_get_thread_team_data();
912 m_result_ptr ? m_result_ptr : pointer_type(data.pool_reduce_local());
914 reference_type update =
915 ValueInit::init(ReducerConditional::select(m_functor, m_reducer), ptr);
919 Kokkos::Impl::FunctorFinal<ReducerTypeFwd, WorkTagFwd>::final(
920 ReducerConditional::select(m_functor, m_reducer), ptr);
923 template <
class HostViewType>
925 const FunctorType& arg_functor,
const MDRangePolicy& arg_policy,
926 const HostViewType& arg_result_view,
927 typename std::enable_if<Kokkos::is_view<HostViewType>::value &&
928 !Kokkos::is_reducer_type<ReducerType>::value,
929 void*>::type =
nullptr)
930 : m_functor(arg_functor),
931 m_mdr_policy(arg_policy),
932 m_policy(Policy(0, m_mdr_policy.m_num_tiles).set_chunk_size(1)),
933 m_reducer(InvalidType()),
934 m_result_ptr(arg_result_view.data()) {
935 static_assert(Kokkos::is_view<HostViewType>::value,
936 "Kokkos::Serial reduce result must be a View");
939 std::is_same<typename HostViewType::memory_space, HostSpace>::value,
940 "Kokkos::Serial reduce result must be a View in HostSpace");
943 inline ParallelReduce(
const FunctorType& arg_functor,
944 MDRangePolicy arg_policy,
const ReducerType& reducer)
945 : m_functor(arg_functor),
946 m_mdr_policy(arg_policy),
947 m_policy(Policy(0, m_mdr_policy.m_num_tiles).set_chunk_size(1)),
949 m_result_ptr(reducer.view().data()) {
967 template <
class FunctorType,
class... Properties>
968 class ParallelFor<FunctorType, Kokkos::TeamPolicy<Properties...>,
971 enum { TEAM_REDUCE_SIZE = 512 };
973 typedef TeamPolicyInternal<Kokkos::Serial, Properties...> Policy;
974 typedef typename Policy::member_type Member;
976 const FunctorType m_functor;
980 template <
class TagType>
981 inline typename std::enable_if<std::is_same<TagType, void>::value>::type exec(
982 HostThreadTeamData& data)
const {
983 for (
int ileague = 0; ileague < m_league; ++ileague) {
984 m_functor(Member(data, ileague, m_league));
988 template <
class TagType>
989 inline typename std::enable_if<!std::is_same<TagType, void>::value>::type
990 exec(HostThreadTeamData& data)
const {
992 for (
int ileague = 0; ileague < m_league; ++ileague) {
993 m_functor(t, Member(data, ileague, m_league));
998 inline void execute()
const {
999 const size_t pool_reduce_size = 0;
1000 const size_t team_reduce_size = TEAM_REDUCE_SIZE;
1001 const size_t team_shared_size = m_shared;
1002 const size_t thread_local_size = 0;
1004 serial_resize_thread_team_data(pool_reduce_size, team_reduce_size,
1005 team_shared_size, thread_local_size);
1007 HostThreadTeamData& data = *serial_get_thread_team_data();
1009 this->
template exec<typename Policy::work_tag>(data);
1012 ParallelFor(
const FunctorType& arg_functor,
const Policy& arg_policy)
1013 : m_functor(arg_functor),
1014 m_league(arg_policy.league_size()),
1015 m_shared(arg_policy.scratch_size(0) + arg_policy.scratch_size(1) +
1016 FunctorTeamShmemSize<FunctorType>::value(arg_functor, 1)) {}
1021 template <
class FunctorType,
class ReducerType,
class... Properties>
1022 class ParallelReduce<FunctorType, Kokkos::TeamPolicy<Properties...>,
1023 ReducerType, Kokkos::Serial> {
1025 enum { TEAM_REDUCE_SIZE = 512 };
1027 typedef TeamPolicyInternal<Kokkos::Serial, Properties...> Policy;
1029 typedef FunctorAnalysis<FunctorPatternInterface::REDUCE, Policy, FunctorType>
1032 typedef typename Policy::member_type Member;
1033 typedef typename Policy::work_tag WorkTag;
1035 typedef Kokkos::Impl::if_c<std::is_same<InvalidType, ReducerType>::value,
1036 FunctorType, ReducerType>
1038 typedef typename ReducerConditional::type ReducerTypeFwd;
1040 typename Kokkos::Impl::if_c<std::is_same<InvalidType, ReducerType>::value,
1041 WorkTag,
void>::type WorkTagFwd;
1043 typedef Kokkos::Impl::FunctorValueInit<ReducerTypeFwd, WorkTagFwd> ValueInit;
1045 typedef typename Analysis::pointer_type pointer_type;
1046 typedef typename Analysis::reference_type reference_type;
1048 const FunctorType m_functor;
1050 const ReducerType m_reducer;
1051 pointer_type m_result_ptr;
1054 template <
class TagType>
1055 inline typename std::enable_if<std::is_same<TagType, void>::value>::type exec(
1056 HostThreadTeamData& data, reference_type update)
const {
1057 for (
int ileague = 0; ileague < m_league; ++ileague) {
1058 m_functor(Member(data, ileague, m_league), update);
1062 template <
class TagType>
1063 inline typename std::enable_if<!std::is_same<TagType, void>::value>::type
1064 exec(HostThreadTeamData& data, reference_type update)
const {
1067 for (
int ileague = 0; ileague < m_league; ++ileague) {
1068 m_functor(t, Member(data, ileague, m_league), update);
1073 inline void execute()
const {
1074 const size_t pool_reduce_size =
1075 Analysis::value_size(ReducerConditional::select(m_functor, m_reducer));
1077 const size_t team_reduce_size = TEAM_REDUCE_SIZE;
1078 const size_t team_shared_size = m_shared;
1079 const size_t thread_local_size = 0;
1081 serial_resize_thread_team_data(pool_reduce_size, team_reduce_size,
1082 team_shared_size, thread_local_size);
1084 HostThreadTeamData& data = *serial_get_thread_team_data();
1087 m_result_ptr ? m_result_ptr : pointer_type(data.pool_reduce_local());
1089 reference_type update =
1090 ValueInit::init(ReducerConditional::select(m_functor, m_reducer), ptr);
1092 this->
template exec<WorkTag>(data, update);
1094 Kokkos::Impl::FunctorFinal<ReducerTypeFwd, WorkTagFwd>::final(
1095 ReducerConditional::select(m_functor, m_reducer), ptr);
1098 template <
class ViewType>
1100 const FunctorType& arg_functor,
const Policy& arg_policy,
1101 const ViewType& arg_result,
1102 typename std::enable_if<Kokkos::is_view<ViewType>::value &&
1103 !Kokkos::is_reducer_type<ReducerType>::value,
1104 void*>::type =
nullptr)
1105 : m_functor(arg_functor),
1106 m_league(arg_policy.league_size()),
1107 m_reducer(InvalidType()),
1108 m_result_ptr(arg_result.data()),
1109 m_shared(arg_policy.scratch_size(0) + arg_policy.scratch_size(1) +
1110 FunctorTeamShmemSize<FunctorType>::value(m_functor, 1)) {
1111 static_assert(Kokkos::is_view<ViewType>::value,
1112 "Reduction result on Kokkos::Serial must be a Kokkos::View");
1115 std::is_same<typename ViewType::memory_space, Kokkos::HostSpace>::value,
1116 "Reduction result on Kokkos::Serial must be a Kokkos::View in "
1120 inline ParallelReduce(
const FunctorType& arg_functor, Policy arg_policy,
1121 const ReducerType& reducer)
1122 : m_functor(arg_functor),
1123 m_league(arg_policy.league_size()),
1125 m_result_ptr(reducer.view().data()),
1126 m_shared(arg_policy.scratch_size(0) + arg_policy.scratch_size(1) +
1127 FunctorTeamShmemSize<FunctorType>::value(arg_functor, 1)) {
1142 namespace Experimental {
1145 class UniqueToken<Serial, UniqueTokenScope::Instance> {
1147 using execution_space = Serial;
1148 using size_type = int;
1153 UniqueToken(execution_space
const& = execution_space()) noexcept {}
1156 KOKKOS_INLINE_FUNCTION
1157 int size() const noexcept {
return 1; }
1160 KOKKOS_INLINE_FUNCTION
1161 int acquire() const noexcept {
return 0; }
1164 KOKKOS_INLINE_FUNCTION
1165 void release(
int) const noexcept {}
1169 class UniqueToken<Serial, UniqueTokenScope::Global> {
1171 using execution_space = Serial;
1172 using size_type = int;
1177 UniqueToken(execution_space
const& = execution_space()) noexcept {}
1180 KOKKOS_INLINE_FUNCTION
1181 int size() const noexcept {
return 1; }
1184 KOKKOS_INLINE_FUNCTION
1185 int acquire() const noexcept {
return 0; }
1188 KOKKOS_INLINE_FUNCTION
1189 void release(
int) const noexcept {}
1195 #include <impl/Kokkos_Serial_Task.hpp>
1197 #endif // defined( KOKKOS_ENABLE_SERIAL )
void print_configuration(std::ostream &, const bool detail=false)
Print "Bill of Materials".
Declaration of various MemoryLayout options.
Declaration of parallel operators.
void finalize()
Finalize the spaces that were initialized via Kokkos::initialize.
Execution policy for work over a range of an integral type.