1 #ifndef VERTEX_LIST_HPP
2 #define VERTEX_LIST_HPP
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);
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);
19 Kokkos::parallel_reduce(
"check_condition", Kokkos::RangePolicy<typename Device::execution_space>(0, n), KOKKOS_LAMBDA(
const int i,
int& nranges_l){
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);
62 KOKKOS_INLINE_FUNCTION
int size()
const{
67 View<IntegerRange*,CLayout,HostType>
ranges;
76 View<int*,CLayout,DeviceType> mapping_l(
NoInit(
"mapping"),
size());
80 for (
int is=0; is<
ranges.size(); is++){
81 int r_offset =
ranges(is).start;
83 mapping_l(offset + i) = r_offset + i;
85 offset +=
ranges(is).size();
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);
111 for(
int i=0; i<n; i++){
113 ranges(i_range).start = i;
116 ranges(i_range).end = i+1;
123 VertexList(
const View<int*,CLayout,HostType>& unordered_full_list,
bool one_indexed);
134 for (
int is=0; is<
ranges.size(); is++){
135 if(
ranges(is).is_in_range(i))
return true;
142 KOKKOS_INLINE_FUNCTION
int size()
const{
146 for (
int is=0; is<
ranges.size(); is++){
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){
164 lambda_func(mapping_l(i));
168 for (
int is=0; is<
ranges.size(); is++){
174 template<
class Device>
175 void pack_contiguous(
const View<double*, CLayout, Device>& input,
const View<double*, CLayout, Device>& contiguous)
const;
177 template<
class Device>
178 void unpack_contiguous(
const View<double*, CLayout, Device>& contiguous,
const View<double*, CLayout, Device>& output)
const;
184 void shift(
int shift_val);
188 void batch(
int first_val,
int nvals,
const View<int*, CLayout, HostType>& view)
const;
196 void print(
const std::string& msg)
const;
206 for_all<typename T::execution_space>(
"vertex_list_set_to_value", KOKKOS_LAMBDA(
const int i){
212 template<
typename T,
typename F>
214 for_all<typename T::execution_space>(
"vertex_list_set_to_value", KOKKOS_LAMBDA(
const int i){
222 template<
class Device>
224 View<int*, CLayout, Device> view(
NoInit(
"is_inside"), nnode_in);
226 Kokkos::deep_copy(view, 1);
229 for_all<typename Device::execution_space>(
"set_is_inside", KOKKOS_LAMBDA(
const int i){
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);
244 for_all<typename Device::execution_space>(
"set_view_in_list", KOKKOS_LAMBDA(
const int i){
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