16#ifndef CABANA_VERLETLIST_HPP
17#define CABANA_VERLETLIST_HPP
23#include <Kokkos_Core.hpp>
24#include <Kokkos_Profiling_ScopedRegion.hpp>
46template <
class MemorySpace,
class LayoutTag>
50template <
class MemorySpace>
57 Kokkos::View<int*, memory_space>
counts;
66 KOKKOS_INLINE_FUNCTION
70 Kokkos::atomic_fetch_add( &
counts( pid ), 1 ) ) = nid;
74 KOKKOS_INLINE_FUNCTION
75 void setNeighbor(
const int pid,
const int nid,
const int new_id )
const
82template <
class MemorySpace>
89 Kokkos::View<int*, memory_space>
counts;
99 KOKKOS_INLINE_FUNCTION
102 std::size_t count = Kokkos::atomic_fetch_add( &
counts( pid ), 1 );
108 KOKKOS_INLINE_FUNCTION
109 void setNeighbor(
const int pid,
const int nid,
const int new_id )
const
126template <
class DeviceType,
class RandomAccessPositionType,
class RadiusType,
127 class AlgorithmTag,
class LayoutTag,
class BuildOpTag>
128struct VerletListBuilder
131 using device = DeviceType;
132 using PositionValueType =
typename RandomAccessPositionType::value_type;
133 using memory_space =
typename device::memory_space;
134 using execution_space =
typename device::execution_space;
137 VerletListData<memory_space, LayoutTag> _data;
139 PositionValueType rsqr;
144 RandomAccessPositionType _position;
145 std::size_t pid_begin, pid_end;
148 BinningData<memory_space> bin_data_1d;
149 LinkedCellList<memory_space, PositionValueType> linked_cell_list;
159 template <
class PositionType>
160 VerletListBuilder( PositionType positions,
const std::size_t begin,
161 const std::size_t end,
162 const RadiusType neighborhood_radius,
163 const PositionValueType cell_size_ratio,
164 const PositionValueType grid_min[3],
165 const PositionValueType grid_max[3],
166 const std::size_t max_neigh )
169 , alloc_n( max_neigh )
171 init( positions, neighborhood_radius, cell_size_ratio, grid_min,
175 radius = neighborhood_radius;
180 template <
class PositionType>
181 VerletListBuilder( PositionType positions,
const std::size_t begin,
182 const std::size_t end,
183 const PositionValueType background_radius,
184 const RadiusType neighborhood_radius,
185 const PositionValueType cell_size_ratio,
186 const PositionValueType grid_min[3],
187 const PositionValueType grid_max[3],
188 const std::size_t max_neigh )
191 , alloc_n( max_neigh )
193 assert( positions.size() == neighborhood_radius.size() );
194 init( positions, background_radius, cell_size_ratio, grid_min,
200 radius = neighborhood_radius;
203 template <
class PositionType>
204 void init( PositionType positions,
205 const PositionValueType neighborhood_radius,
206 const PositionValueType cell_size_ratio,
207 const PositionValueType grid_min[3],
208 const PositionValueType grid_max[3] )
214 _data.counts = Kokkos::View<int*, memory_space>(
"num_neighbors",
218 initCounts( LayoutTag() );
221 _position = positions;
227 double grid_size = cell_size_ratio * neighborhood_radius;
228 PositionValueType grid_delta[3] = { grid_size, grid_size, grid_size };
229 linked_cell_list = createLinkedCellList<memory_space>(
230 _position, grid_delta, grid_min, grid_max, neighborhood_radius,
232 bin_data_1d = linked_cell_list.binningData();
235 rsqr = neighborhood_radius * neighborhood_radius;
238 KOKKOS_INLINE_FUNCTION
auto cutoff( [[maybe_unused]]
const int p )
const
242 if constexpr ( is_slice<RadiusType>::value ||
243 Kokkos::is_view<RadiusType>::value )
244 return radius( p ) * radius( p );
251 struct CountNeighborsTag
254 using CountNeighborsPolicy =
255 Kokkos::TeamPolicy<execution_space, CountNeighborsTag,
256 Kokkos::IndexType<int>,
257 Kokkos::Schedule<Kokkos::Dynamic>>;
259 KOKKOS_INLINE_FUNCTION
261 operator()(
const CountNeighborsTag&,
262 const typename CountNeighborsPolicy::member_type& team )
const
266 int cell = team.league_rank();
269 int imin, imax, jmin, jmax, kmin, kmax;
270 linked_cell_list.getStencilCells( cell, imin, imax, jmin, jmax, kmin,
274 std::size_t b_offset = bin_data_1d.binOffset( cell );
275 Kokkos::parallel_for(
276 Kokkos::TeamThreadRange( team, 0, bin_data_1d.binSize( cell ) ),
281 std::size_t pid = linked_cell_list.permutation( bi + b_offset );
283 if ( ( pid >= pid_begin ) && ( pid < pid_end ) )
286 double x_p = _position( pid, 0 );
287 double y_p = _position( pid, 1 );
288 double z_p = _position( pid, 2 );
291 int stencil_count = 0;
292 for ( int i = imin; i < imax; ++i )
293 for ( int j = jmin; j < jmax; ++j )
294 for ( int k = kmin; k < kmax; ++k )
298 if ( linked_cell_list.cellStencil()
299 .grid.minDistanceToPoint(
300 x_p, y_p, z_p, i, j, k ) <= rsqr )
302 std::size_t n_offset =
303 linked_cell_list.binOffset( i, j, k );
305 linked_cell_list.binSize( i, j, k );
311 neighbor_reduce( team, pid, x_p, y_p, z_p,
313 cell_count, BuildOpTag() );
314 stencil_count += cell_count;
317 Kokkos::single( Kokkos::PerThread( team ), [&]()
318 { _data.counts( pid ) = stencil_count; } );
324 KOKKOS_INLINE_FUNCTION
void
325 neighbor_reduce(
const typename CountNeighborsPolicy::member_type& team,
326 const std::size_t pid,
const double x_p,
const double y_p,
327 const double z_p,
const int n_offset,
const int num_n,
328 int& cell_count, TeamVectorOpTag )
const
330 Kokkos::parallel_reduce(
331 Kokkos::ThreadVectorRange( team, num_n ),
332 [&](
const int n,
int& local_count ) {
333 neighbor_kernel( pid, x_p, y_p, z_p, n_offset, n, local_count );
339 KOKKOS_INLINE_FUNCTION
340 void neighbor_reduce(
const typename CountNeighborsPolicy::member_type,
341 const std::size_t pid,
const double x_p,
342 const double y_p,
const double z_p,
343 const int n_offset,
const int num_n,
int& cell_count,
346 for (
int n = 0; n < num_n; n++ )
347 neighbor_kernel( pid, x_p, y_p, z_p, n_offset, n, cell_count );
351 KOKKOS_INLINE_FUNCTION
352 void neighbor_kernel(
const int pid,
const double x_p,
const double y_p,
353 const double z_p,
const int n_offset,
const int n,
354 int& local_count )
const
357 std::size_t nid = linked_cell_list.permutation( n_offset + n );
360 double x_n = _position( nid, 0 );
361 double y_n = _position( nid, 1 );
362 double z_n = _position( nid, 2 );
365 if ( NeighborDiscriminator<AlgorithmTag>::isValid(
366 pid, x_p, y_p, z_p, nid, x_n, y_n, z_n ) )
370 PositionValueType dx = x_p - x_n;
371 PositionValueType dy = y_p - y_n;
372 PositionValueType dz = z_p - z_n;
373 PositionValueType dist_sqr = dx * dx + dy * dy + dz * dz;
376 if ( dist_sqr <= cutoff( pid ) )
383 template <
class KokkosMemorySpace>
386 using kokkos_mem_space = KokkosMemorySpace;
387 Kokkos::View<int*, kokkos_mem_space> counts;
388 Kokkos::View<int*, kokkos_mem_space> offsets;
389 KOKKOS_INLINE_FUNCTION
390 void operator()(
const int i,
int& update,
const bool final_pass )
const
398 void initCounts( VerletLayoutCSR ) {}
400 void initCounts( VerletLayout2D )
406 _data.neighbors = Kokkos::View<int**, memory_space>(
407 Kokkos::ViewAllocateWithoutInitializing(
"neighbors" ),
408 _data.counts.size(), alloc_n );
412 void processCounts( VerletLayoutCSR )
415 _data.offsets = Kokkos::View<int*, memory_space>(
416 Kokkos::ViewAllocateWithoutInitializing(
"neighbor_offsets" ),
417 _data.counts.size() );
420 OffsetScanOp<memory_space> offset_op;
421 offset_op.counts = _data.counts;
422 offset_op.offsets = _data.offsets;
423 int total_num_neighbor;
424 Kokkos::RangePolicy<execution_space> range_policy(
425 0, _data.counts.extent( 0 ) );
426 Kokkos::parallel_scan(
"Cabana::VerletListBuilder::offset_scan",
427 range_policy, offset_op, total_num_neighbor );
431 _data.neighbors = Kokkos::View<int*, memory_space>(
432 Kokkos::ViewAllocateWithoutInitializing(
"neighbors" ),
433 total_num_neighbor );
436 Kokkos::deep_copy( _data.counts, 0 );
441 void processCounts( VerletLayout2D )
444 auto counts = _data.counts;
446 Kokkos::Max<int> max_reduce( max );
447 Kokkos::parallel_reduce(
448 "Cabana::VerletListBuilder::reduce_max",
449 Kokkos::RangePolicy<execution_space>( 0, _data.counts.size() ),
450 KOKKOS_LAMBDA(
const int i,
int& value ) {
451 if ( counts( i ) > value )
456 _data.max_n =
static_cast<std::size_t
>( max );
459 if ( count or _data.max_n > _data.neighbors.extent( 1 ) )
462 Kokkos::deep_copy( _data.counts, 0 );
463 _data.neighbors = Kokkos::View<int**, memory_space>(
464 Kokkos::ViewAllocateWithoutInitializing(
"neighbors" ),
465 _data.counts.size(), _data.max_n );
470 struct FillNeighborsTag
473 using FillNeighborsPolicy =
474 Kokkos::TeamPolicy<execution_space, FillNeighborsTag,
475 Kokkos::IndexType<int>,
476 Kokkos::Schedule<Kokkos::Dynamic>>;
477 KOKKOS_INLINE_FUNCTION
479 operator()(
const FillNeighborsTag&,
480 const typename FillNeighborsPolicy::member_type& team )
const
484 int cell = team.league_rank();
487 int imin, imax, jmin, jmax, kmin, kmax;
488 linked_cell_list.getStencilCells( cell, imin, imax, jmin, jmax, kmin,
492 std::size_t b_offset = bin_data_1d.binOffset( cell );
493 Kokkos::parallel_for(
494 Kokkos::TeamThreadRange( team, 0, bin_data_1d.binSize( cell ) ),
499 std::size_t pid = linked_cell_list.permutation( bi + b_offset );
501 if ( ( pid >= pid_begin ) && ( pid < pid_end ) )
504 double x_p = _position( pid, 0 );
505 double y_p = _position( pid, 1 );
506 double z_p = _position( pid, 2 );
509 for ( int i = imin; i < imax; ++i )
510 for ( int j = jmin; j < jmax; ++j )
511 for ( int k = kmin; k < kmax; ++k )
515 if ( linked_cell_list.cellStencil()
516 .grid.minDistanceToPoint(
517 x_p, y_p, z_p, i, j, k ) <= rsqr )
521 std::size_t n_offset =
522 linked_cell_list.binOffset( i, j, k );
524 linked_cell_list.binSize( i, j, k );
525 neighbor_for( team, pid, x_p, y_p, z_p,
535 KOKKOS_INLINE_FUNCTION
void
536 neighbor_for(
const typename FillNeighborsPolicy::member_type& team,
537 const std::size_t pid,
const double x_p,
const double y_p,
538 const double z_p,
const int n_offset,
const int num_n,
541 Kokkos::parallel_for(
542 Kokkos::ThreadVectorRange( team, num_n ), [&](
const int n )
543 { neighbor_kernel( pid, x_p, y_p, z_p, n_offset, n ); } );
547 KOKKOS_INLINE_FUNCTION
548 void neighbor_for(
const typename FillNeighborsPolicy::member_type team,
549 const std::size_t pid,
const double x_p,
550 const double y_p,
const double z_p,
const int n_offset,
553 for (
int n = 0; n < num_n; n++ )
555 Kokkos::PerThread( team ),
556 [&]() { neighbor_kernel( pid, x_p, y_p, z_p, n_offset, n ); } );
560 KOKKOS_INLINE_FUNCTION
561 void neighbor_kernel(
const int pid,
const double x_p,
const double y_p,
562 const double z_p,
const int n_offset,
566 std::size_t nid = linked_cell_list.permutation( n_offset + n );
569 double x_n = _position( nid, 0 );
570 double y_n = _position( nid, 1 );
571 double z_n = _position( nid, 2 );
575 pid, x_p, y_p, z_p, nid, x_n, y_n, z_n ) )
579 PositionValueType dx = x_p - x_n;
580 PositionValueType dy = y_p - y_n;
581 PositionValueType dz = z_p - z_n;
582 PositionValueType dist_sqr = dx * dx + dy * dy + dz * dz;
586 if ( dist_sqr <= cutoff( pid ) )
588 _data.addNeighbor( pid, nid );
596template <
class DeviceType,
class AlgorithmTag,
class LayoutTag,
597 class BuildOpTag,
class PositionType>
598auto createVerletListBuilder(
599 PositionType x,
const std::size_t begin,
const std::size_t end,
600 const typename PositionType::value_type radius,
601 const typename PositionType::value_type cell_size_ratio,
602 const typename PositionType::value_type grid_min[3],
603 const typename PositionType::value_type grid_max[3],
604 const std::size_t max_neigh,
605 typename std::enable_if<( is_slice<PositionType>::value ),
int>::type* = 0 )
607 using RandomAccessPositionType =
typename PositionType::random_access_slice;
608 return VerletListBuilder<DeviceType, RandomAccessPositionType,
609 typename PositionType::value_type, AlgorithmTag,
610 LayoutTag, BuildOpTag>(
611 x, begin, end, radius, cell_size_ratio, grid_min, grid_max, max_neigh );
614template <
class DeviceType,
class AlgorithmTag,
class LayoutTag,
615 class BuildOpTag,
class PositionType>
616auto createVerletListBuilder(
617 PositionType x,
const std::size_t begin,
const std::size_t end,
618 const typename PositionType::value_type radius,
619 const typename PositionType::value_type cell_size_ratio,
620 const typename PositionType::value_type grid_min[3],
621 const typename PositionType::value_type grid_max[3],
622 const std::size_t max_neigh,
623 typename std::enable_if<( Kokkos::is_view<PositionType>::value ),
626 using RandomAccessPositionType =
627 Kokkos::View<
typename PositionType::value_type**, DeviceType,
628 Kokkos::MemoryTraits<Kokkos::RandomAccess>>;
629 return VerletListBuilder<DeviceType, RandomAccessPositionType,
630 typename PositionType::value_type, AlgorithmTag,
631 LayoutTag, BuildOpTag>(
632 x, begin, end, radius, cell_size_ratio, grid_min, grid_max, max_neigh );
635template <
class DeviceType,
class AlgorithmTag,
class LayoutTag,
636 class BuildOpTag,
class PositionType,
class RadiusType>
637auto createVerletListBuilder(
638 PositionType x,
const std::size_t begin,
const std::size_t end,
639 const typename PositionType::value_type background_radius,
640 const RadiusType radius,
641 const typename PositionType::value_type cell_size_ratio,
642 const typename PositionType::value_type grid_min[3],
643 const typename PositionType::value_type grid_max[3],
644 const std::size_t max_neigh,
645 typename std::enable_if<( is_slice<PositionType>::value ),
int>::type* = 0 )
647 using RandomAccessPositionType =
typename PositionType::random_access_slice;
648 return VerletListBuilder<DeviceType, RandomAccessPositionType, RadiusType,
649 AlgorithmTag, LayoutTag, BuildOpTag>(
650 x, begin, end, background_radius, radius, cell_size_ratio, grid_min,
651 grid_max, max_neigh );
654template <
class DeviceType,
class AlgorithmTag,
class LayoutTag,
655 class BuildOpTag,
class PositionType,
class RadiusType>
656auto createVerletListBuilder(
657 PositionType x,
const std::size_t begin,
const std::size_t end,
658 const typename PositionType::value_type background_radius,
659 const RadiusType radius,
660 const typename PositionType::value_type cell_size_ratio,
661 const typename PositionType::value_type grid_min[3],
662 const typename PositionType::value_type grid_max[3],
663 const std::size_t max_neigh,
664 typename std::enable_if<( Kokkos::is_view<PositionType>::value ),
667 using RandomAccessPositionType =
668 Kokkos::View<
typename PositionType::value_type**, DeviceType,
669 Kokkos::MemoryTraits<Kokkos::RandomAccess>>;
670 return VerletListBuilder<DeviceType, RandomAccessPositionType, RadiusType,
671 AlgorithmTag, LayoutTag, BuildOpTag>(
672 x, begin, end, background_radius, radius, cell_size_ratio, grid_min,
673 grid_max, max_neigh );
699template <
class MemorySpace,
class AlgorithmTag,
class LayoutTag,
700 class BuildTag = TeamVectorOpTag>
704 static_assert( Kokkos::is_memory_space<MemorySpace>::value,
"" );
754 template <
class PositionType>
756 PositionType x,
const std::size_t begin,
const std::size_t end,
757 const typename PositionType::value_type neighborhood_radius,
758 const typename PositionType::value_type cell_size_ratio,
759 const typename PositionType::value_type grid_min[3],
760 const typename PositionType::value_type grid_max[3],
761 const std::size_t max_neigh = 0,
763 Kokkos::is_view<PositionType>::value ),
766 build( x, begin, end, neighborhood_radius, cell_size_ratio, grid_min,
767 grid_max, max_neigh );
806 template <
class PositionSlice,
class RadiusSlice>
807 VerletList( PositionSlice x,
const std::size_t begin,
const std::size_t end,
808 const typename PositionSlice::value_type background_radius,
809 RadiusSlice neighborhood_radius,
810 const typename PositionSlice::value_type cell_size_ratio,
811 const typename PositionSlice::value_type grid_min[3],
812 const typename PositionSlice::value_type grid_max[3],
813 const std::size_t max_neigh = 0,
817 build( x, begin, end, background_radius, neighborhood_radius,
818 cell_size_ratio, grid_min, grid_max, max_neigh );
825 template <
class PositionType>
827 build( PositionType x,
const std::size_t begin,
const std::size_t end,
828 const typename PositionType::value_type neighborhood_radius,
829 const typename PositionType::value_type cell_size_ratio,
830 const typename PositionType::value_type grid_min[3],
831 const typename PositionType::value_type grid_max[3],
832 const std::size_t max_neigh = 0,
834 Kokkos::is_view<PositionType>::value ),
839 cell_size_ratio, grid_min, grid_max, max_neigh );
845 template <
class PositionType,
class ExecutionSpace>
847 build( ExecutionSpace, PositionType x,
const std::size_t begin,
848 const std::size_t end,
849 const typename PositionType::value_type neighborhood_radius,
850 const typename PositionType::value_type cell_size_ratio,
851 const typename PositionType::value_type grid_min[3],
852 const typename PositionType::value_type grid_max[3],
853 const std::size_t max_neigh = 0,
855 Kokkos::is_view<PositionType>::value ),
858 Kokkos::Profiling::ScopedRegion region(
"Cabana::VerletList::build" );
862 assert( end >= begin );
863 assert( end <=
size( x ) );
865 using device_type = Kokkos::Device<ExecutionSpace, memory_space>;
867 auto builder = Impl::createVerletListBuilder<device_type, AlgorithmTag,
868 LayoutTag, BuildTag>(
869 x, begin, end, neighborhood_radius, cell_size_ratio, grid_min,
870 grid_max, max_neigh );
871 buildImpl( builder );
878 template <
class PositionSlice,
class RadiusSlice>
879 void build( PositionSlice x,
const std::size_t begin,
const std::size_t end,
880 const typename PositionSlice::value_type background_radius,
881 RadiusSlice neighborhood_radius,
882 const typename PositionSlice::value_type cell_size_ratio,
883 const typename PositionSlice::value_type grid_min[3],
884 const typename PositionSlice::value_type grid_max[3],
885 const std::size_t max_neigh = 0 )
889 neighborhood_radius, cell_size_ratio, grid_min, grid_max,
896 template <
class PositionSlice,
class RadiusSlice,
class ExecutionSpace>
897 void build( ExecutionSpace, PositionSlice x,
const std::size_t begin,
898 const std::size_t end,
899 const typename PositionSlice::value_type background_radius,
900 RadiusSlice neighborhood_radius,
901 const typename PositionSlice::value_type cell_size_ratio,
902 const typename PositionSlice::value_type grid_min[3],
903 const typename PositionSlice::value_type grid_max[3],
904 const std::size_t max_neigh = 0 )
906 Kokkos::Profiling::ScopedRegion region(
"Cabana::VerletList::build" );
910 assert( end >= begin );
911 assert( end <= x.size() );
914 using device_type = Kokkos::Device<ExecutionSpace, memory_space>;
915 auto builder = Impl::createVerletListBuilder<device_type, AlgorithmTag,
916 LayoutTag, BuildTag>(
917 x, begin, end, background_radius, neighborhood_radius,
918 cell_size_ratio, grid_min, grid_max, max_neigh );
919 buildImpl( builder );
923 template <
class BuilderType>
924 void buildImpl( BuilderType builder )
934 typename BuilderType::FillNeighborsPolicy fill_policy(
935 builder.bin_data_1d.numBin(), Kokkos::AUTO, 4 );
938 typename BuilderType::CountNeighborsPolicy count_policy(
939 builder.bin_data_1d.numBin(), Kokkos::AUTO, 4 );
940 Kokkos::parallel_for(
"Cabana::VerletList::count_neighbors",
941 count_policy, builder );
945 builder.processCounts( LayoutTag() );
946 Kokkos::parallel_for(
"Cabana::VerletList::fill_neighbors",
947 fill_policy, builder );
953 builder.processCounts( LayoutTag() );
957 if ( builder.count or builder.refill )
959 Kokkos::parallel_for(
"Cabana::VerletList::fill_neighbors",
960 fill_policy, builder );
965 _data = builder._data;
970 KOKKOS_INLINE_FUNCTION
972 const std::size_t neighbor_index,
973 const int new_index )
const
975 _data.setNeighbor( particle_index, neighbor_index, new_index );
983template <
class MemorySpace,
class AlgorithmTag,
class BuildTag>
995 KOKKOS_INLINE_FUNCTION
999 return list.
_data.neighbors.extent( 0 );
1003 KOKKOS_INLINE_FUNCTION
1006 std::size_t num_p = list.
_data.counts.size();
1011 KOKKOS_INLINE_FUNCTION
1013 const std::size_t particle_index )
1015 return list.
_data.counts( particle_index );
1020 KOKKOS_INLINE_FUNCTION
1022 const std::size_t particle_index,
1023 const std::size_t neighbor_index )
1025 return list.
_data.neighbors( list.
_data.offsets( particle_index ) +
1032template <
class MemorySpace,
class AlgorithmTag,
class BuildTag>
1044 KOKKOS_INLINE_FUNCTION
1047 std::size_t num_p = list.
_data.counts.size();
1052 KOKKOS_INLINE_FUNCTION
1056 return list.
_data.max_n;
1060 KOKKOS_INLINE_FUNCTION
1062 const std::size_t particle_index )
1064 return list.
_data.counts( particle_index );
1069 KOKKOS_INLINE_FUNCTION
1071 const std::size_t particle_index,
1072 const std::size_t neighbor_index )
1074 return list.
_data.neighbors( particle_index, neighbor_index );
std::enable_if_t< 3==Array_t::num_space_dim, void > update(Array_t &a, const typename Array_t::value_type alpha, const Array_t &b, const typename Array_t::value_type beta, DecompositionTag tag)
Update two vectors such that a = alpha * a + beta * b. 3D specialization.
Definition Cabana_Grid_Array.hpp:595
KOKKOS_INLINE_FUNCTION std::enable_if_t< 3==SplineDataType::num_space_dim, void > value(const ViewType &view, const SplineDataType &sd, PointDataType &result, typename std::enable_if<(std::rank< PointDataType >::value==0), void * >::type=0)
Interpolate a scalar value to a point. 3D specialization.
Definition Cabana_Grid_Interpolation.hpp:56
Linked cell list binning (spatial sorting) and neighbor iteration.
KOKKOS_INLINE_FUNCTION std::size_t maxNeighbor(const ListType &list, const std::size_t num_particles)
Iterate to find the maximum number of neighbors.
Definition Cabana_NeighborList.hpp:164
KOKKOS_INLINE_FUNCTION std::size_t totalNeighbor(const ListType &list, const std::size_t num_particles)
Iterate to get the total number of neighbors.
Definition Cabana_NeighborList.hpp:152
SIMD and neighbor extension of Kokkos parallel iteration.
Neighborhood discriminator.
Definition Cabana_NeighborList.hpp:55
static KOKKOS_INLINE_FUNCTION std::size_t getNeighbor(const list_type &list, const std::size_t particle_index, const std::size_t neighbor_index)
Definition Cabana_VerletList.hpp:1070
VerletList< MemorySpace, AlgorithmTag, VerletLayout2D, BuildTag > list_type
Neighbor list type.
Definition Cabana_VerletList.hpp:1040
MemorySpace memory_space
Kokkos memory space.
Definition Cabana_VerletList.hpp:1038
static KOKKOS_INLINE_FUNCTION std::size_t numNeighbor(const list_type &list, const std::size_t particle_index)
Get the number of neighbors for a given particle index.
Definition Cabana_VerletList.hpp:1061
static KOKKOS_INLINE_FUNCTION std::size_t totalNeighbor(const list_type &list)
Get the total number of neighbors across all particles.
Definition Cabana_VerletList.hpp:1045
static KOKKOS_INLINE_FUNCTION std::size_t maxNeighbor(const list_type &list)
Get the maximum number of neighbors per particle.
Definition Cabana_VerletList.hpp:1053
static KOKKOS_INLINE_FUNCTION std::size_t numNeighbor(const list_type &list, const std::size_t particle_index)
Get the number of neighbors for a given particle index.
Definition Cabana_VerletList.hpp:1012
VerletList< MemorySpace, AlgorithmTag, VerletLayoutCSR, BuildTag > list_type
Neighbor list type.
Definition Cabana_VerletList.hpp:991
static KOKKOS_INLINE_FUNCTION std::size_t maxNeighbor(const list_type &list)
Get the maximum number of neighbors across all particles.
Definition Cabana_VerletList.hpp:1004
static KOKKOS_INLINE_FUNCTION std::size_t totalNeighbor(const list_type &list)
Get the total number of neighbors across all particles.
Definition Cabana_VerletList.hpp:996
static KOKKOS_INLINE_FUNCTION std::size_t getNeighbor(const list_type &list, const std::size_t particle_index, const std::size_t neighbor_index)
Definition Cabana_VerletList.hpp:1021
MemorySpace memory_space
Kokkos memory space.
Definition Cabana_VerletList.hpp:989
Neighbor list interface. Provides an interface callable at the functor level that gives access to nei...
Definition Cabana_NeighborList.hpp:114
Neighbor operations are executed with team parallelism.
Definition Cabana_Parallel.hpp:206
Neighbor operations are executed with team vector parallelism.
Definition Cabana_Parallel.hpp:211
Neighbor list implementation based on binning particles on a 3d Cartesian grid with cells of the same...
Definition Cabana_VerletList.hpp:702
void build(PositionSlice x, const std::size_t begin, const std::size_t end, const typename PositionSlice::value_type background_radius, RadiusSlice neighborhood_radius, const typename PositionSlice::value_type cell_size_ratio, const typename PositionSlice::value_type grid_min[3], const typename PositionSlice::value_type grid_max[3], const std::size_t max_neigh=0)
Given a list of particle positions and a neighborhood radius calculate the neighbor list.
Definition Cabana_VerletList.hpp:879
VerletListData< memory_space, VerletLayoutCSR > _data
Definition Cabana_VerletList.hpp:713
VerletList(PositionType x, const std::size_t begin, const std::size_t end, const typename PositionType::value_type neighborhood_radius, const typename PositionType::value_type cell_size_ratio, const typename PositionType::value_type grid_min[3], const typename PositionType::value_type grid_max[3], const std::size_t max_neigh=0, typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type *=0)
VerletList constructor. Given a list of particle positions and a neighborhood radius calculate the ne...
Definition Cabana_VerletList.hpp:755
KOKKOS_INLINE_FUNCTION void setNeighbor(const std::size_t particle_index, const std::size_t neighbor_index, const int new_index) const
Modify a neighbor in the list; for example, mark it as a broken bond.
Definition Cabana_VerletList.hpp:971
typename memory_space::execution_space execution_space
Kokkos default execution space for this memory space.
Definition Cabana_VerletList.hpp:710
void build(ExecutionSpace, PositionType x, const std::size_t begin, const std::size_t end, const typename PositionType::value_type neighborhood_radius, const typename PositionType::value_type cell_size_ratio, const typename PositionType::value_type grid_min[3], const typename PositionType::value_type grid_max[3], const std::size_t max_neigh=0, typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type *=0)
Given a list of particle positions and a neighborhood radius calculate the neighbor list.
Definition Cabana_VerletList.hpp:847
MemorySpace memory_space
Kokkos memory space in which the neighbor list data resides.
Definition Cabana_VerletList.hpp:707
VerletList()
Default constructor.
Definition Cabana_VerletList.hpp:718
void build(ExecutionSpace, PositionSlice x, const std::size_t begin, const std::size_t end, const typename PositionSlice::value_type background_radius, RadiusSlice neighborhood_radius, const typename PositionSlice::value_type cell_size_ratio, const typename PositionSlice::value_type grid_min[3], const typename PositionSlice::value_type grid_max[3], const std::size_t max_neigh=0)
Given a list of particle positions and a neighborhood radius calculate the neighbor list.
Definition Cabana_VerletList.hpp:897
VerletList(PositionSlice x, const std::size_t begin, const std::size_t end, const typename PositionSlice::value_type background_radius, RadiusSlice neighborhood_radius, const typename PositionSlice::value_type cell_size_ratio, const typename PositionSlice::value_type grid_min[3], const typename PositionSlice::value_type grid_max[3], const std::size_t max_neigh=0, typename std::enable_if<(is_slice< PositionSlice >::value), int >::type *=0)
VerletList constructor. Given a list of particle positions and a neighborhood radius calculate the ne...
Definition Cabana_VerletList.hpp:807
void build(PositionType x, const std::size_t begin, const std::size_t end, const typename PositionType::value_type neighborhood_radius, const typename PositionType::value_type cell_size_ratio, const typename PositionType::value_type grid_min[3], const typename PositionType::value_type grid_max[3], const std::size_t max_neigh=0, typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type *=0)
Definition Cabana_VerletList.hpp:827
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
auto size(SliceType slice, typename std::enable_if< is_slice< SliceType >::value, int >::type *=0)
Check slice size (differs from Kokkos View).
Definition Cabana_Slice.hpp:1012
2D array neighbor list layout.
Definition Cabana_VerletList.hpp:40
CSR (compressed sparse row) neighbor list layout.
Definition Cabana_VerletList.hpp:35
KOKKOS_INLINE_FUNCTION void addNeighbor(const int pid, const int nid) const
Add a neighbor to the list.
Definition Cabana_VerletList.hpp:100
MemorySpace memory_space
Kokkos memory space.
Definition Cabana_VerletList.hpp:86
Kokkos::View< int *, memory_space > counts
Number of neighbors per particle.
Definition Cabana_VerletList.hpp:89
Kokkos::View< int **, memory_space > neighbors
Neighbor list.
Definition Cabana_VerletList.hpp:92
KOKKOS_INLINE_FUNCTION void setNeighbor(const int pid, const int nid, const int new_id) const
Modify a neighbor in the list.
Definition Cabana_VerletList.hpp:109
std::size_t max_n
Definition Cabana_VerletList.hpp:96
KOKKOS_INLINE_FUNCTION void setNeighbor(const int pid, const int nid, const int new_id) const
Modify a neighbor in the list.
Definition Cabana_VerletList.hpp:75
Kokkos::View< int *, memory_space > offsets
Offsets into the neighbor list.
Definition Cabana_VerletList.hpp:60
MemorySpace memory_space
Kokkos memory space.
Definition Cabana_VerletList.hpp:54
Kokkos::View< int *, memory_space > counts
Number of neighbors per particle.
Definition Cabana_VerletList.hpp:57
KOKKOS_INLINE_FUNCTION void addNeighbor(const int pid, const int nid) const
Add a neighbor to the list.
Definition Cabana_VerletList.hpp:67
Kokkos::View< int *, memory_space > neighbors
Neighbor list.
Definition Cabana_VerletList.hpp:63
Definition Cabana_VerletList.hpp:47
Definition Cabana_Types.hpp:88
Slice static type checker.
Definition Cabana_Slice.hpp:861