IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
ViewUtils.h
Go to the documentation of this file.
1//
2// View Utilities
3// Utility functions relating to Kokkos views
4//
5
6#ifndef IPPL_VIEW_UTILS_H
7#define IPPL_VIEW_UTILS_H
8
9#include <Kokkos_Core.hpp>
10
11#include "Types/ViewTypes.h"
12
13namespace ippl {
14 namespace detail {
31 template <unsigned Dim, unsigned Current = 0, class BeginFunctor, class EndFunctor,
32 class Functor, typename Check = std::nullptr_t>
33 constexpr void nestedLoop(BeginFunctor&& begin, EndFunctor&& end, Functor&& body,
34 Check&& check = nullptr) {
35 for (size_t i = begin(Current); i < end(Current); ++i) {
36 if constexpr (Dim - 1 == Current) {
37 body(i);
38 } else {
39 auto inner = [i, &body](auto... args) {
40 body(i, args...);
41 };
42 nestedLoop<Dim, Current + 1>(begin, end, inner, check);
43 }
44 }
45 if constexpr (!std::is_null_pointer_v<std::decay_t<Check>>) {
46 check(Current);
47 }
48 }
49
61 template <typename View, class Functor, typename Check = std::nullptr_t>
62 constexpr void nestedViewLoop(View& view, int shift, Functor&& body,
63 Check&& check = nullptr) {
65 [&](unsigned) {
66 return shift;
67 },
68 [&](unsigned d) {
69 return view.extent(d) - shift;
70 },
71 body, check);
72 }
73
83 template <typename T, unsigned Dim, class... Properties>
85 std::ostream& out = std::cout) {
86 using view_type = typename ViewType<T, Dim, Properties...>::view_type;
87 typename view_type::HostMirror hview = Kokkos::create_mirror_view(view);
88 Kokkos::deep_copy(hview, view);
89
91 view, 0,
92 [&]<typename... Args>(Args&&... args) {
93 out << hview(args...) << " ";
94 },
95 [&](unsigned axis) {
96 if (axis + 1 >= 2 || axis == 0) {
97 out << std::endl;
98 }
99 });
100 }
101
105 template <typename View, size_t... Idx>
106 decltype(auto) shrinkView_impl(std::string label, const View& view, int nghost,
107 const std::index_sequence<Idx...>&) {
108 using view_type = typename Kokkos::View<typename View::data_type, Kokkos::LayoutLeft,
109 typename View::memory_space>::uniform_type;
110 return view_type(label, (view.extent(Idx) - 2 * nghost)...);
111 }
112
121 template <typename View>
122 decltype(auto) shrinkView(std::string label, const View& view, int nghost) {
123 return shrinkView_impl(label, view, nghost, std::make_index_sequence<View::rank>{});
124 }
125 } // namespace detail
126} // namespace ippl
127
128#endif
typename ippl::detail::ViewType< ippl::Vector< double, Dim >, 1 >::view_type view_type
constexpr unsigned Dim
Definition Archive.h:20
decltype(auto) shrinkView(std::string label, const View &view, int nghost)
Definition ViewUtils.h:122
void write(const typename ViewType< T, Dim, Properties... >::view_type &view, std::ostream &out=std::cout)
Definition ViewUtils.h:84
constexpr void nestedViewLoop(View &view, int shift, Functor &&body, Check &&check=nullptr)
Definition ViewUtils.h:62
decltype(auto) shrinkView_impl(std::string label, const View &view, int nghost, const std::index_sequence< Idx... > &)
Definition ViewUtils.h:106
constexpr void nestedLoop(BeginFunctor &&begin, EndFunctor &&end, Functor &&body, Check &&check=nullptr)
Definition ViewUtils.h:33
Kokkos::View< typename NPtr< T, Dim >::type, Properties... > view_type
Definition ViewTypes.h:45