XGCa
vertex_list.hpp
Go to the documentation of this file.
1 #ifndef VERTEX_LIST_HPP
2 #define VERTEX_LIST_HPP
3 
4 #include "space_settings.hpp"
5 
6 template<typename Device, typename F>
7 int get_starts_and_ends(int n, F condition, View<bool*, CLayout, HostType>& starts_h, View<bool*, CLayout, HostType>& ends_h){
8  View<bool*, CLayout, Device> starts(NoInit("starts"), n);
9  View<bool*, CLayout, Device> ends(NoInit("ends"), n);
10 
11  // First, determine which vertices go on the list
12  View<bool*, CLayout, Device> is_in(NoInit("is_in"), n);
13  Kokkos::parallel_for("check_condition", Kokkos::RangePolicy<typename Device::execution_space>(0, n), KOKKOS_LAMBDA(const int i){
14  is_in(i) = condition(i);
15  });
16 
17  // Next, subtract to get list starts and stops
18  int nranges = 0;
19  Kokkos::parallel_reduce("check_condition", Kokkos::RangePolicy<typename Device::execution_space>(0, n), KOKKOS_LAMBDA(const int i, int& nranges_l){
20  starts(i) = false;
21  ends(i) = false;
22  if(is_in(i)){
23  if(i==n-1){
24  ends(i) = true;
25  nranges_l++;
26  }else{
27  if(!is_in(i+1)){ // List just ended
28  ends(i) = true;
29  nranges_l++;
30  }
31  }
32 
33  if(i==0){
34  starts(i) = true;
35  }else{
36  if(!is_in(i-1)){ // List just started
37  starts(i) = true;
38  }
39  }
40  }
41  }, nranges);
42 
43  ends_h = View<bool*, CLayout, HostType>(NoInit("ends"), n);
44  starts_h = View<bool*, CLayout, HostType>(NoInit("starts"), n);
45  Kokkos::deep_copy(ends_h, ends);
46  Kokkos::deep_copy(starts_h, starts);
47 
48  return nranges;
49 }
50 
51 // A compressed way to represent a list of vertices.
52 // Note: using C/python convention so range {0, 2} means vertices 0 and 1
53 class VertexList{
54  struct IntegerRange{
55  int start;
56  int end; // Note: using C/python convention so range {0, 2} means vertices 0 and 1
57 
58  KOKKOS_INLINE_FUNCTION bool is_in_range(int i) const{
59  return (i>=start && i<end);
60  }
61 
62  KOKKOS_INLINE_FUNCTION int size() const{
63  return (end - start);
64  }
65  };
66  public:
67  View<IntegerRange*,CLayout,HostType> ranges;
68 
69  View<int*,CLayout,DeviceType> mapping;
71 
72  public:
73 
74  void plan(){
75  // Start with assumption that we should just use mapping method
76  View<int*,CLayout,DeviceType> mapping_l(NoInit("mapping"), size());
77 
78  // Loop through ranges
79  int offset = 0;
80  for (int is=0; is<ranges.size(); is++){
81  int r_offset = ranges(is).start;
82  Kokkos::parallel_for("vertex_list_parallel_for", Kokkos::RangePolicy<ExSpace>( 0, ranges(is).size() ), KOKKOS_LAMBDA( const int i){
83  mapping_l(offset + i) = r_offset + i;
84  });
85  offset += ranges(is).size();
86  }
87  Kokkos::fence();
88  mapping = mapping_l;
89  use_mapping = true;
90  }
91 
93  : ranges("ranges", 0), use_mapping(false)
94  {}
95 
96  // Initialize with single range
97  VertexList(int start, int end);
98 
99  // Initialize with a condition, e.g.:
100  // VertexList(grid.nnode, KOKKOS_LAMBDA( const int i){return grid_h.is_in_psi_range(i)});
101  template<typename F>
102  VertexList(int n, F condition){
103  use_mapping=false;
104  View<bool*, CLayout, HostType> starts_h;
105  View<bool*, CLayout, HostType> ends_h;
106  int nranges = get_starts_and_ends<DeviceType>(n, condition, starts_h, ends_h);
107  ranges = View<IntegerRange*,CLayout,HostType>(NoInit("ranges"), nranges);
108 
109  // Now populate views
110  int i_range = 0;
111  for(int i=0; i<n; i++){
112  if(starts_h(i)){
113  ranges(i_range).start = i;
114  }
115  if(ends_h(i)){
116  ranges(i_range).end = i+1; // i is included, therefore the list ends at i+1
117  i_range += 1;
118  }
119  }
120  }
121 
122  // Initialize the list by passing in the full uncompressed list of vertices.
123  VertexList(const View<int*,CLayout,HostType>& unordered_full_list, bool one_indexed);
124 
132  KOKKOS_INLINE_FUNCTION bool is_in_list(int i) const {
133  // First check ranges
134  for (int is=0; is<ranges.size(); is++){
135  if(ranges(is).is_in_range(i)) return true;
136  }
137 
138  return false;
139  }
140 
141  // Number of vertices
142  KOKKOS_INLINE_FUNCTION int size() const{
143  int n = 0;
144 
145  // Add ranges
146  for (int is=0; is<ranges.size(); is++){
147  n += ranges(is).size();
148  }
149 
150  return n;
151  }
152 
158  template<typename DeviceExSpace, typename F>
159  inline void for_all(const std::string label, F lambda_func) const {
160  if(use_mapping && std::is_same<DeviceExSpace, ExSpace>::value){
161  // Sparse list - use a mapping instead of a separate parallel_for for each range
162  auto mapping_l = mapping;
163  Kokkos::parallel_for(label, Kokkos::RangePolicy<DeviceExSpace>(0, mapping.extent(0)), KOKKOS_LAMBDA(const int i){
164  lambda_func(mapping_l(i));
165  });
166  }else{
167  // Loop through ranges
168  for (int is=0; is<ranges.size(); is++){
169  Kokkos::parallel_for(label, Kokkos::RangePolicy<DeviceExSpace>(ranges(is).start, ranges(is).end), lambda_func);
170  }
171  }
172  }
173 
174  template<class Device>
175  void pack_contiguous(const View<double*, CLayout, Device>& input, const View<double*, CLayout, Device>& contiguous) const;
176 
177  template<class Device>
178  void unpack_contiguous(const View<double*, CLayout, Device>& contiguous, const View<double*, CLayout, Device>& output) const;
179 
180  VertexList operator&(const VertexList& list2) const;
181 
182  VertexList& operator|=(const VertexList& rhs); // Custom +=
183 
184  void shift(int shift_val);
185 
186  VertexList batch(int first_val, int nvals) const;
187 
188  void batch(int first_val, int nvals, const View<int*, CLayout, HostType>& view) const;
189 
190  int min() const;
191 
192  int max() const;
193 
194  bool is_contiguous() const;
195 
196  void print(const std::string& msg) const;
197 
198  bool is_subset(const VertexList& list2) const;
199 
200  bool intersects(const VertexList& list2) const;
201 
202  bool intersects(int start, int end) const;
203 
204  template<typename T>
205  void set_to_value(const T& view, double value) const{
206  for_all<typename T::execution_space>("vertex_list_set_to_value", KOKKOS_LAMBDA(const int i){
207  view(i) = value;
208  });
209  Kokkos::fence();
210  }
211 
212  template<typename T, typename F>
213  void set_to_value(const T& view, F func) const{
214  for_all<typename T::execution_space>("vertex_list_set_to_value", KOKKOS_LAMBDA(const int i){
215  view(i) = func(i);
216  });
217  Kokkos::fence();
218  }
219 
220  // Returns a int view of size nnode_in. Indices in the last are set to 0, while all others are set to 1.
221  // This function exists for historical reasons and usage should be replaced with get_view_in_list.
222  template<class Device>
223  View<int*, CLayout, Device> get_view_int_not_in_list(int nnode_in) const{
224  View<int*, CLayout, Device> view(NoInit("is_inside"), nnode_in);
225  // Set all values to 1
226  Kokkos::deep_copy(view, 1);
227 
228  // Set values inside list to 0
229  for_all<typename Device::execution_space>("set_is_inside", KOKKOS_LAMBDA(const int i){
230  view(i) = 0;
231  });
232  Kokkos::fence();
233 
234  return view;
235  }
236 
237  // Returns a bool view of size nnode_in that is true if the index appears in the list
238  template<class Device>
239  View<bool*, CLayout, Device> get_view_in_list(const std::string label, int nnode_in) const{
240  View<bool*, CLayout, Device> view(NoInit(label), nnode_in);
241  Kokkos::deep_copy(view, false);
242 
243  // Set values inside list to true
244  for_all<typename Device::execution_space>("set_view_in_list", KOKKOS_LAMBDA(const int i){
245  view(i) = true;
246  });
247  Kokkos::fence();
248 
249  return view;
250  }
251 };
252 
253 #endif
Definition: vertex_list.hpp:53
VertexList batch(int first_val, int nvals) const
Definition: vertex_list.cpp:255
VertexList operator&(const VertexList &list2) const
Definition: vertex_list.cpp:125
bool use_mapping
Definition: vertex_list.hpp:70
View< IntegerRange *, CLayout, HostType > ranges
Definition: vertex_list.hpp:67
void set_to_value(const T &view, double value) const
Definition: vertex_list.hpp:205
bool is_contiguous() const
Definition: vertex_list.cpp:327
VertexList()
Definition: vertex_list.hpp:92
KOKKOS_INLINE_FUNCTION int size() const
Definition: vertex_list.hpp:142
KOKKOS_INLINE_FUNCTION bool is_in_list(int i) const
Definition: vertex_list.hpp:132
bool intersects(const VertexList &list2) const
Definition: vertex_list.cpp:118
void plan()
Definition: vertex_list.hpp:74
void shift(int shift_val)
Definition: vertex_list.cpp:246
View< int *, CLayout, DeviceType > mapping
Definition: vertex_list.hpp:69
VertexList(int n, F condition)
Definition: vertex_list.hpp:102
bool is_subset(const VertexList &list2) const
Definition: vertex_list.cpp:342
void pack_contiguous(const View< double *, CLayout, Device > &input, const View< double *, CLayout, Device > &contiguous) const
Definition: vertex_list.cpp:58
void set_to_value(const T &view, F func) const
Definition: vertex_list.hpp:213
void unpack_contiguous(const View< double *, CLayout, Device > &contiguous, const View< double *, CLayout, Device > &output) const
Definition: vertex_list.cpp:74
View< bool *, CLayout, Device > get_view_in_list(const std::string label, int nnode_in) const
Definition: vertex_list.hpp:239
View< int *, CLayout, Device > get_view_int_not_in_list(int nnode_in) const
Definition: vertex_list.hpp:223
VertexList & operator|=(const VertexList &rhs)
Definition: vertex_list.cpp:173
int max() const
Definition: vertex_list.cpp:322
void print(const std::string &msg) const
Definition: vertex_list.cpp:331
void for_all(const std::string label, F lambda_func) const
Definition: vertex_list.hpp:159
int min() const
Definition: vertex_list.cpp:316
void parallel_for(const std::string name, int n_ptl, Function func, Option option, HostAoSoA aosoa_h, DeviceAoSoA aosoa_d)
Definition: streamed_parallel_for.hpp:252
logical false
Definition: module.F90:102
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:69
Definition: vertex_list.hpp:54
int end
Definition: vertex_list.hpp:56
KOKKOS_INLINE_FUNCTION bool is_in_range(int i) const
Definition: vertex_list.hpp:58
int start
Definition: vertex_list.hpp:55
KOKKOS_INLINE_FUNCTION int size() const
Definition: vertex_list.hpp:62
int get_starts_and_ends(int n, F condition, View< bool *, CLayout, HostType > &starts_h, View< bool *, CLayout, HostType > &ends_h)
Definition: vertex_list.hpp:7