Kokkos Core Kernels Package Version of the Day
Loading...
Searching...
No Matches
Kokkos_Tuners.hpp
1//@HEADER
2// ************************************************************************
3//
4// Kokkos v. 4.0
5// Copyright (2022) National Technology & Engineering
6// Solutions of Sandia, LLC (NTESS).
7//
8// Under the terms of Contract DE-NA0003525 with NTESS,
9// the U.S. Government retains certain rights in this software.
10//
11// Part of Kokkos, under the Apache License v2.0 with LLVM Exceptions.
12// See https://kokkos.org/LICENSE for license information.
13// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
14//
15//@HEADER
16
17#ifndef KOKKOS_IMPL_PUBLIC_INCLUDE
18#include <Kokkos_Macros.hpp>
19static_assert(false,
20 "Including non-public Kokkos header files is not allowed.");
21#endif
22#ifndef KOKKOS_KOKKOS_TUNERS_HPP
23#define KOKKOS_KOKKOS_TUNERS_HPP
24
25#include <Kokkos_Macros.hpp>
26#include <Kokkos_Core_fwd.hpp>
27#include <Kokkos_ExecPolicy.hpp>
28#include <KokkosExp_MDRangePolicy.hpp>
29#include <impl/Kokkos_Profiling_Interface.hpp>
30
31#include <array>
32#include <utility>
33#include <tuple>
34#include <string>
35#include <vector>
36#include <map>
37#include <cassert>
38
39namespace Kokkos {
40namespace Tools {
41
42namespace Experimental {
43
44// forward declarations
45SetOrRange make_candidate_set(size_t size, int64_t* data);
46bool have_tuning_tool();
47size_t declare_output_type(const std::string&,
48 Kokkos::Tools::Experimental::VariableInfo);
49void request_output_values(size_t, size_t,
50 Kokkos::Tools::Experimental::VariableValue*);
51VariableValue make_variable_value(size_t, int64_t);
52VariableValue make_variable_value(size_t, double);
53SetOrRange make_candidate_range(double lower, double upper, double step,
54 bool openLower, bool openUpper);
55SetOrRange make_candidate_range(int64_t lower, int64_t upper, int64_t step,
56 bool openLower, bool openUpper);
57size_t get_new_context_id();
58void begin_context(size_t context_id);
59void end_context(size_t context_id);
60namespace Impl {
61
67template <typename ValueType, typename ContainedType>
68struct ValueHierarchyNode;
69
70template <typename ValueType, typename ContainedType>
72 std::vector<ValueType> root_values;
73 std::vector<ContainedType> sub_values;
74 void add_root_value(const ValueType& in) { root_values.push_back(in); }
75 void add_sub_container(const ContainedType& in) { sub_values.push_back(in); }
76 const ValueType& get_root_value(const size_t index) const {
77 return root_values[index];
78 }
79 const ContainedType& get_sub_value(const size_t index) const {
80 return sub_values[index];
81 }
82};
83
84template <typename ValueType>
86 std::vector<ValueType> root_values;
87 explicit ValueHierarchyNode(std::vector<ValueType> rv)
88 : root_values(std::move(rv)) {}
89 void add_root_value(const ValueType& in) { root_values.push_back(in); }
90 const ValueType& get_root_value(const size_t index) const {
91 return root_values[index];
92 }
93};
94
100template <class NestedMap>
102
103// Vectors are our lowest-level, no nested values
104template <class T>
105struct MapTypeConverter<std::vector<T>> {
106 using type = ValueHierarchyNode<T, void>;
107};
108
109// Maps contain both the "root" types and sub-vectors
110template <class K, class V>
111struct MapTypeConverter<std::map<K, V>> {
113};
114
120template <class NestedMap>
122
123// Vectors are our lowest-level, no nested values. Just fill in the fundamental
124// values
125template <class T>
126struct ValueHierarchyConstructor<std::vector<T>> {
127 using return_type = typename MapTypeConverter<std::vector<T>>::type;
128 static return_type build(const std::vector<T>& in) { return return_type{in}; }
129};
130
131// For maps, we need to fill in the fundamental values, and construct child
132// nodes
133template <class K, class V>
134struct ValueHierarchyConstructor<std::map<K, V>> {
135 using return_type = typename MapTypeConverter<std::map<K, V>>::type;
136 static return_type build(const std::map<K, V>& in) {
137 return_type node_to_build;
138 for (auto& entry : in) {
139 node_to_build.add_root_value(entry.first);
140 node_to_build.add_sub_container(
141 ValueHierarchyConstructor<V>::build(entry.second));
142 }
143 return node_to_build;
144 }
145};
146
155template <class InspectForDepth>
157
158// The dimensionality of a vector is 1
159template <class T>
160struct get_space_dimensionality<std::vector<T>> {
161 static constexpr int value = 1;
162};
163
164// The dimensionality of a map is 1 (the map) plus the dimensionality
165// of the map's value type
166template <class K, class V>
167struct get_space_dimensionality<std::map<K, V>> {
168 static constexpr int value = 1 + get_space_dimensionality<V>::value;
169};
170
171template <class T, int N>
172struct n_dimensional_sparse_structure;
173
174template <class T>
175struct n_dimensional_sparse_structure<T, 1> {
176 using type = std::vector<T>;
177};
178
179template <class T, int N>
180struct n_dimensional_sparse_structure {
181 using type =
182 std::map<T, typename n_dimensional_sparse_structure<T, N - 1>::type>;
183};
184
191// First, a helper to get the value in one dimension
192template <class Container>
194
195// At any given level, just return your value at that level
196template <class RootType, class Subtype>
198 static RootType get(const ValueHierarchyNode<RootType, Subtype>& dimension,
199 double fraction_to_traverse) {
200 size_t index = dimension.root_values.size() * fraction_to_traverse;
201 return dimension.get_root_value(index);
202 }
203};
204
210// At the bottom level, we have one double and a base-level ValueHierarchyNode
211
212template <class HierarchyNode, class... InterpolationIndices>
214
215template <class ValueType>
217 using node_type = ValueHierarchyNode<ValueType, void>;
218 using return_type = std::tuple<ValueType>;
219 static return_type build(const node_type& in, double index) {
220 return std::make_tuple(DimensionValueExtractor<node_type>::get(in, index));
221 }
222};
223
224// At levels above the bottom, we tuple_cat the result of our child on the end
225// of our own tuple
226template <class ValueType, class Subtype, class... Indices>
227struct GetMultidimensionalPoint<ValueHierarchyNode<ValueType, Subtype>, double,
228 Indices...> {
230 using sub_tuple =
231 typename GetMultidimensionalPoint<Subtype, Indices...>::return_type;
232 using return_type = decltype(std::tuple_cat(
233 std::declval<std::tuple<ValueType>>(), std::declval<sub_tuple>()));
234 static return_type build(const node_type& in, double fraction_to_traverse,
235 Indices... indices) {
236 size_t index = in.sub_values.size() * fraction_to_traverse;
237 auto dimension_value = std::make_tuple(
238 DimensionValueExtractor<node_type>::get(in, fraction_to_traverse));
239 return std::tuple_cat(dimension_value,
240 GetMultidimensionalPoint<Subtype, Indices...>::build(
241 in.get_sub_value(index), indices...));
242 }
243};
244
245template <typename PointType, class ArrayType, size_t... Is>
246auto get_point_helper(const PointType& in, const ArrayType& indices,
247 std::index_sequence<Is...>) {
248 using helper = GetMultidimensionalPoint<
249 PointType,
250 decltype(std::get<Is>(std::declval<ArrayType>()).value.double_value)...>;
251 return helper::build(in, std::get<Is>(indices).value.double_value...);
252}
253
254template <typename PointType, typename ArrayType>
255struct GetPoint;
256
257template <typename PointType, size_t ArraySize>
258struct GetPoint<
259 PointType,
260 std::array<Kokkos::Tools::Experimental::VariableValue, ArraySize>> {
261 using index_set_type =
262 std::array<Kokkos::Tools::Experimental::VariableValue, ArraySize>;
263 static auto build(const PointType& in, const index_set_type& indices) {
264 return get_point_helper(in, indices, std::make_index_sequence<ArraySize>{});
265 }
266};
267
268template <typename PointType, typename ArrayType>
269auto get_point(const PointType& point, const ArrayType& indices) {
270 return GetPoint<PointType, ArrayType>::build(point, indices);
271}
272
273} // namespace Impl
274
275template <template <class...> class Container, size_t MaxDimensionSize = 100,
276 class... TemplateArguments>
277class MultidimensionalSparseTuningProblem {
278 public:
279 using ProblemSpaceInput = Container<TemplateArguments...>;
280 static constexpr int space_dimensionality =
281 Impl::get_space_dimensionality<ProblemSpaceInput>::value;
282 static constexpr size_t max_space_dimension_size = MaxDimensionSize;
283 static constexpr double tuning_min = 0.0;
284 static constexpr double tuning_max = 0.999;
285
286 // Not declared as static constexpr to work around the following compiler bug
287 // https://gcc.gnu.org/bugzilla/show_bug.cgi?id=96862
288 // where a floating-point expression cannot be constexpr under -frounding-math
289 double tuning_step = tuning_max / max_space_dimension_size;
290
291 using StoredProblemSpace =
292 typename Impl::MapTypeConverter<ProblemSpaceInput>::type;
293 using HierarchyConstructor =
294 typename Impl::ValueHierarchyConstructor<Container<TemplateArguments...>>;
295
296 using ValueArray = std::array<Kokkos::Tools::Experimental::VariableValue,
297 space_dimensionality>;
298 template <class Key, class Value>
299 using extended_map = std::map<Key, Value>;
300 template <typename Key>
301 using extended_problem =
302 MultidimensionalSparseTuningProblem<extended_map, MaxDimensionSize, Key,
303 ProblemSpaceInput>;
304 template <typename Key, typename Value>
305 using ExtendedProblemSpace =
306 typename Impl::MapTypeConverter<extended_map<Key, Value>>::type;
307
308 template <typename Key>
309 auto extend(const std::string& axis_name,
310 const std::vector<Key>& new_tuning_axis) const
311 -> extended_problem<Key> {
312 ExtendedProblemSpace<Key, ProblemSpaceInput> extended_space;
313 for (auto& key : new_tuning_axis) {
314 extended_space.add_root_value(key);
315 extended_space.add_sub_container(m_space);
316 }
317 std::vector<std::string> extended_names;
318 extended_names.reserve(m_variable_names.size() + 1);
319 extended_names.push_back(axis_name);
320 extended_names.insert(extended_names.end(), m_variable_names.begin(),
321 m_variable_names.end());
322 return extended_problem<Key>(extended_space, extended_names);
323 }
324
325 private:
326 StoredProblemSpace m_space;
327 std::array<size_t, space_dimensionality> variable_ids;
328 std::vector<std::string> m_variable_names;
329 size_t context;
330
331 public:
332 MultidimensionalSparseTuningProblem() = default;
333
334 MultidimensionalSparseTuningProblem(StoredProblemSpace space,
335 const std::vector<std::string>& names)
336 : m_space(std::move(space)), m_variable_names(names) {
337 KOKKOS_ASSERT(names.size() == space_dimensionality);
338 for (unsigned long x = 0; x < names.size(); ++x) {
339 VariableInfo info;
340 info.type = Kokkos::Tools::Experimental::ValueType::kokkos_value_double;
341 info.category = Kokkos::Tools::Experimental::StatisticalCategory::
342 kokkos_value_interval;
343 info.valueQuantity =
344 Kokkos::Tools::Experimental::CandidateValueType::kokkos_value_range;
345 info.candidates = Kokkos::Tools::Experimental::make_candidate_range(
346 tuning_min, tuning_max, tuning_step, true, true);
347 variable_ids[x] = declare_output_type(names[x], info);
348 }
349 }
350
351 MultidimensionalSparseTuningProblem(ProblemSpaceInput space,
352 const std::vector<std::string>& names)
353 : MultidimensionalSparseTuningProblem(HierarchyConstructor::build(space),
354 names) {}
355
356 template <typename... Coordinates>
357 auto get_point(Coordinates... coordinates) {
358 using ArrayType = std::array<Kokkos::Tools::Experimental::VariableValue,
359 sizeof...(coordinates)>;
360 return Impl::get_point(
361 m_space, ArrayType({Kokkos::Tools::Experimental::make_variable_value(
362 0, static_cast<double>(coordinates))...}));
363 }
364
365 auto begin() {
366 context = Kokkos::Tools::Experimental::get_new_context_id();
367 ValueArray values;
368 for (int x = 0; x < space_dimensionality; ++x) {
369 values[x] = Kokkos::Tools::Experimental::make_variable_value(
370 variable_ids[x], 0.0);
371 }
372 begin_context(context);
373 request_output_values(context, space_dimensionality, values.data());
374 return Impl::get_point(m_space, values);
375 }
376
377 auto end() { end_context(context); }
378};
379
380template <typename Tuner>
381struct ExtendableTunerMixin {
382 template <typename Key>
383 auto combine(const std::string& axis_name,
384 const std::vector<Key>& new_axis) const {
385 const auto& sub_tuner = static_cast<const Tuner*>(this)->get_tuner();
386 return sub_tuner.extend(axis_name, new_axis);
387 }
388
389 template <typename... Coordinates>
390 auto get_point(Coordinates... coordinates) {
391 const auto& sub_tuner = static_cast<const Tuner*>(this)->get_tuner();
392 return sub_tuner.get_point(coordinates...);
393 }
394
395 private:
396 ExtendableTunerMixin() = default;
397 friend Tuner;
398};
399
400template <size_t MaxDimensionSize = 100, template <class...> class Container,
401 class... TemplateArguments>
402auto make_multidimensional_sparse_tuning_problem(
403 const Container<TemplateArguments...>& in, std::vector<std::string> names) {
404 return MultidimensionalSparseTuningProblem<Container, MaxDimensionSize,
405 TemplateArguments...>(in, names);
406}
407
408class TeamSizeTuner : public ExtendableTunerMixin<TeamSizeTuner> {
409 private:
410 using SpaceDescription = std::map<int64_t, std::vector<int64_t>>;
411 using TunerType = decltype(make_multidimensional_sparse_tuning_problem<20>(
412 std::declval<SpaceDescription>(),
413 std::declval<std::vector<std::string>>()));
414 TunerType tuner;
415
416 public:
417 TeamSizeTuner() = default;
418 TeamSizeTuner& operator=(const TeamSizeTuner& other) = default;
419 TeamSizeTuner(const TeamSizeTuner& other) = default;
420 TeamSizeTuner& operator=(TeamSizeTuner&& other) = default;
421 TeamSizeTuner(TeamSizeTuner&& other) = default;
422 template <typename ViableConfigurationCalculator, typename Functor,
423 typename TagType, typename... Properties>
424 TeamSizeTuner(const std::string& name,
425 const Kokkos::TeamPolicy<Properties...>& policy_in,
426 const Functor& functor, const TagType& tag,
427 ViableConfigurationCalculator calc) {
428 using PolicyType = Kokkos::TeamPolicy<Properties...>;
429 PolicyType policy(policy_in);
430 auto initial_vector_length = policy.impl_vector_length();
431 if (initial_vector_length < 1) {
432 policy.impl_set_vector_length(1);
433 }
459 SpaceDescription space_description;
460
461 auto max_vector_length = PolicyType::vector_length_max();
462 std::vector<int64_t> allowed_vector_lengths;
463
464 if (policy.impl_auto_vector_length()) { // case 1 or 2
465 for (int vector_length = max_vector_length; vector_length >= 1;
466 vector_length /= 2) {
467 policy.impl_set_vector_length(vector_length);
480 auto max_team_size = calc.get_max_team_size(policy, functor, tag);
481 if ((policy.impl_auto_team_size()) ||
482 (policy.team_size() <= max_team_size)) {
483 allowed_vector_lengths.push_back(vector_length);
484 }
485 }
486 } else { // case 3, there's only one vector length to care about
487 allowed_vector_lengths.push_back(policy.impl_vector_length());
488 }
489
490 for (const auto vector_length : allowed_vector_lengths) {
491 std::vector<int64_t> allowed_team_sizes;
492 policy.impl_set_vector_length(vector_length);
493 auto max_team_size = calc.get_max_team_size(policy, functor, tag);
494 if (policy.impl_auto_team_size()) { // case 1 or 3, try all legal team
495 // sizes
496 for (int team_size = max_team_size; team_size >= 1; team_size /= 2) {
497 allowed_team_sizes.push_back(team_size);
498 }
499 } else { // case 2, just try the provided team size
500 allowed_team_sizes.push_back(policy.team_size());
501 }
502 space_description[vector_length] = allowed_team_sizes;
503 }
504 tuner = make_multidimensional_sparse_tuning_problem<20>(
505 space_description, {std::string(name + "_vector_length"),
506 std::string(name + "_team_size")});
507 policy.impl_set_vector_length(initial_vector_length);
508 }
509
510 template <typename... Properties>
511 auto tune(const Kokkos::TeamPolicy<Properties...>& policy_in) {
512 Kokkos::TeamPolicy<Properties...> policy(policy_in);
513 if (Kokkos::Tools::Experimental::have_tuning_tool()) {
514 auto configuration = tuner.begin();
515 auto team_size = std::get<1>(configuration);
516 auto vector_length = std::get<0>(configuration);
517 if (vector_length > 0) {
518 policy.impl_set_team_size(team_size);
519 policy.impl_set_vector_length(vector_length);
520 }
521 }
522 return policy;
523 }
524 void end() {
525 if (Kokkos::Tools::Experimental::have_tuning_tool()) {
526 tuner.end();
527 }
528 }
529
530 TunerType get_tuner() const { return tuner; }
531};
532namespace Impl {
533template <class T>
534struct tuning_type_for;
535
536template <>
537struct tuning_type_for<double> {
538 static constexpr Kokkos::Tools::Experimental::ValueType value =
539 Kokkos::Tools::Experimental::ValueType::kokkos_value_double;
540 static double get(
541 const Kokkos::Tools::Experimental::VariableValue& value_struct) {
542 return value_struct.value.double_value;
543 }
544};
545template <>
546struct tuning_type_for<int64_t> {
547 static constexpr Kokkos::Tools::Experimental::ValueType value =
548 Kokkos::Tools::Experimental::ValueType::kokkos_value_int64;
549 static int64_t get(
550 const Kokkos::Tools::Experimental::VariableValue& value_struct) {
551 return value_struct.value.int_value;
552 }
553};
554} // namespace Impl
555template <class Bound>
556class SingleDimensionalRangeTuner {
557 size_t id;
558 size_t context;
559 using tuning_util = Impl::tuning_type_for<Bound>;
560
561 Bound default_value;
562
563 public:
564 SingleDimensionalRangeTuner() = default;
565 SingleDimensionalRangeTuner(
566 const std::string& name,
567 Kokkos::Tools::Experimental::StatisticalCategory category,
568 Bound default_val, Bound lower, Bound upper, Bound step = (Bound)0) {
569 default_value = default_val;
570 Kokkos::Tools::Experimental::VariableInfo info;
571 info.category = category;
572 info.candidates = make_candidate_range(
573 static_cast<Bound>(lower), static_cast<Bound>(upper),
574 static_cast<Bound>(step), false, false);
575 info.valueQuantity =
576 Kokkos::Tools::Experimental::CandidateValueType::kokkos_value_range;
577 info.type = tuning_util::value;
578 id = Kokkos::Tools::Experimental::declare_output_type(name, info);
579 }
580
581 Bound begin() {
582 context = Kokkos::Tools::Experimental::get_new_context_id();
583 Kokkos::Tools::Experimental::begin_context(context);
584 auto tuned_value =
585 Kokkos::Tools::Experimental::make_variable_value(id, default_value);
586 Kokkos::Tools::Experimental::request_output_values(context, 1,
587 &tuned_value);
588 return tuning_util::get(tuned_value);
589 }
590
591 void end() { Kokkos::Tools::Experimental::end_context(context); }
592
593 template <typename Functor>
594 void with_tuned_value(Functor& func) {
595 func(begin());
596 end();
597 }
598};
599
600class RangePolicyOccupancyTuner {
601 private:
602 using TunerType = SingleDimensionalRangeTuner<int64_t>;
603 TunerType tuner;
604
605 public:
606 RangePolicyOccupancyTuner() = default;
607 template <typename ViableConfigurationCalculator, typename Functor,
608 typename TagType, typename... Properties>
609 RangePolicyOccupancyTuner(const std::string& name,
611 const Functor&, const TagType&,
612 ViableConfigurationCalculator)
613 : tuner(TunerType(name,
614 Kokkos::Tools::Experimental::StatisticalCategory::
615 kokkos_value_ratio,
616 100, 5, 100, 5)) {}
617
618 template <typename... Properties>
619 auto tune(const Kokkos::RangePolicy<Properties...>& policy_in) {
620 Kokkos::RangePolicy<Properties...> policy(policy_in);
621 if (Kokkos::Tools::Experimental::have_tuning_tool()) {
622 auto occupancy = tuner.begin();
623 policy.impl_set_desired_occupancy(
624 Kokkos::Experimental::DesiredOccupancy{static_cast<int>(occupancy)});
625 }
626 return policy;
627 }
628 void end() {
629 if (Kokkos::Tools::Experimental::have_tuning_tool()) {
630 tuner.end();
631 }
632 }
633
634 TunerType get_tuner() const { return tuner; }
635};
636
637namespace Impl {
638
639template <typename T>
640void fill_tile(std::vector<T>& cont, int tile_size) {
641 for (int x = 1; x < tile_size; x *= 2) {
642 cont.push_back(x);
643 }
644}
645template <typename T, typename Mapped>
646void fill_tile(std::map<T, Mapped>& cont, int tile_size) {
647 for (int x = 1; x < tile_size; x *= 2) {
648 fill_tile(cont[x], tile_size / x);
649 }
650}
651} // namespace Impl
652
653template <int MDRangeRank>
654struct MDRangeTuner : public ExtendableTunerMixin<MDRangeTuner<MDRangeRank>> {
655 private:
656 static constexpr int rank = MDRangeRank;
657 static constexpr int max_slices = 15;
658 using SpaceDescription =
659 typename Impl::n_dimensional_sparse_structure<int, rank>::type;
660 using TunerType =
661 decltype(make_multidimensional_sparse_tuning_problem<max_slices>(
662 std::declval<SpaceDescription>(),
663 std::declval<std::vector<std::string>>()));
664 TunerType tuner;
665
666 public:
667 MDRangeTuner() = default;
668 template <typename Functor, typename TagType, typename Calculator,
669 typename... Properties>
670 MDRangeTuner(const std::string& name,
672 const Functor& functor, const TagType& tag, Calculator calc) {
673 SpaceDescription desc;
674 int max_tile_size =
675 calc.get_mdrange_max_tile_size_product(policy, functor, tag);
676 Impl::fill_tile(desc, max_tile_size);
677 std::vector<std::string> feature_names;
678 for (int x = 0; x < rank; ++x) {
679 feature_names.push_back(name + "_tile_size_" + std::to_string(x));
680 }
681 tuner = make_multidimensional_sparse_tuning_problem<max_slices>(
682 desc, feature_names);
683 }
684 template <typename Policy, typename Tuple, size_t... Indices>
685 void set_policy_tile(Policy& policy, const Tuple& tuple,
686 const std::index_sequence<Indices...>&) {
687 policy.impl_change_tile_size({std::get<Indices>(tuple)...});
688 }
689 template <typename... Properties>
690 auto tune(const Kokkos::MDRangePolicy<Properties...>& policy_in) {
691 Kokkos::MDRangePolicy<Properties...> policy(policy_in);
692 if (Kokkos::Tools::Experimental::have_tuning_tool()) {
693 auto configuration = tuner.begin();
694 set_policy_tile(policy, configuration, std::make_index_sequence<rank>{});
695 }
696 return policy;
697 }
698 void end() {
699 if (Kokkos::Tools::Experimental::have_tuning_tool()) {
700 tuner.end();
701 }
702 }
703
704 TunerType get_tuner() const { return tuner; }
705};
706
707template <class Choice>
708struct CategoricalTuner {
709 using choice_list = std::vector<Choice>;
710 choice_list choices;
711 size_t context;
712 size_t tuning_variable_id;
713 CategoricalTuner(std::string name, choice_list m_choices)
714 : choices(m_choices) {
715 std::vector<int64_t> indices;
716 for (typename decltype(choices)::size_type x = 0; x < choices.size(); ++x) {
717 indices.push_back(x);
718 }
719 VariableInfo info;
720 info.category = StatisticalCategory::kokkos_value_categorical;
721 info.valueQuantity = CandidateValueType::kokkos_value_set;
722 info.type = ValueType::kokkos_value_int64;
723 info.candidates = make_candidate_set(indices.size(), indices.data());
724 tuning_variable_id = declare_output_type(name, info);
725 }
726 const Choice& begin() {
727 context = get_new_context_id();
728 begin_context(context);
729 VariableValue value = make_variable_value(tuning_variable_id, int64_t(0));
730 request_output_values(context, 1, &value);
731 return choices[value.value.int_value];
732 }
733 void end() { end_context(context); }
734};
735
736template <typename Choice>
737auto make_categorical_tuner(std::string name, std::vector<Choice> choices)
738 -> CategoricalTuner<Choice> {
739 return CategoricalTuner<Choice>(name, choices);
740}
741
742} // namespace Experimental
743} // namespace Tools
744} // namespace Kokkos
745
746#endif
A thread safe view to a bitset.
KOKKOS_FORCEINLINE_FUNCTION unsigned size() const
Execution policy for work over a range of an integral type.