Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana::LinkedCellList< MemorySpace, Scalar > Class Template Reference

Data describing the bin sizes and offsets resulting from a binning operation on a 3d regular Cartesian grid. More...

#include <Cabana_LinkedCellList.hpp>

Public Types

using memory_space = typename MemorySpace::memory_space
 Memory space.
 
using execution_space = typename memory_space::execution_space
 Default execution space.
 
using size_type = typename memory_space::size_type
 Memory space size type.
 
using CountView = Kokkos::View<int*, memory_space>
 Binning view type.
 
using OffsetView = Kokkos::View<size_type*, memory_space>
 Offset view type.
 
using stencil_type = Cabana::LinkedCellStencil<Scalar>
 Stencil type.
 

Public Member Functions

 LinkedCellList ()
 Default constructor.
 
template<class PositionType>
 LinkedCellList (PositionType positions, const Scalar grid_delta[3], const Scalar grid_min[3], const Scalar grid_max[3], typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type *=0)
 Simple constructor.
 
template<class PositionType>
 LinkedCellList (PositionType positions, const std::size_t begin, const std::size_t end, const Scalar grid_delta[3], const Scalar grid_min[3], const Scalar grid_max[3], typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type *=0)
 Partial range constructor.
 
template<class PositionType>
 LinkedCellList (PositionType positions, const Scalar grid_delta[3], const Scalar grid_min[3], const Scalar grid_max[3], const Scalar neighborhood_radius, const Scalar cell_size_ratio=1, typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type *=0)
 Explicit stencil constructor.
 
template<class PositionType>
 LinkedCellList (PositionType positions, const std::size_t begin, const std::size_t end, const Scalar grid_delta[3], const Scalar grid_min[3], const Scalar grid_max[3], const Scalar neighborhood_radius, const Scalar cell_size_ratio=1, typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type *=0)
 Explicit stencil and partial range constructor.
 
KOKKOS_INLINE_FUNCTION int numParticles () const
 Number of binned particles.
 
KOKKOS_INLINE_FUNCTION std::size_t getParticleBegin () const
 Beginning of binned range.
 
KOKKOS_INLINE_FUNCTION std::size_t getParticleEnd () const
 End of binned range.
 
KOKKOS_INLINE_FUNCTION int totalBins () const
 Get the total number of bins.
 
KOKKOS_INLINE_FUNCTION int numBin (const int dim) const
 Get the number of bins in a given dimension.
 
KOKKOS_INLINE_FUNCTION size_type cardinalBinIndex (const int i, const int j, const int k) const
 Given the ijk index of a bin get its cardinal index.
 
KOKKOS_INLINE_FUNCTION void ijkBinIndex (const int cardinal, int &i, int &j, int &k) const
 Given the cardinal index of a bin get its ijk indices.
 
KOKKOS_INLINE_FUNCTION int binSize (const int i, const int j, const int k) const
 Given a bin get the number of particles it contains.
 
KOKKOS_INLINE_FUNCTION size_type binOffset (const int i, const int j, const int k) const
 Given a bin get the particle index at which it sorts.
 
KOKKOS_INLINE_FUNCTION size_type permutation (const int particle_id) const
 Given a local particle id in the binned layout, get the id of the particle in the old (unbinned) layout.
 
KOKKOS_INLINE_FUNCTION std::size_t rangeBegin () const
 The beginning particle index binned by the linked cell list.
 
KOKKOS_INLINE_FUNCTION std::size_t rangeEnd () const
 The ending particle index binned by the linked cell list.
 
KOKKOS_INLINE_FUNCTION stencil_type cellStencil () const
 Get the linked cell stencil.
 
BinningData< MemorySpace > binningData () const
 Get the 1d bin data.
 
template<class ExecutionSpace, class PositionType>
void build (ExecutionSpace, PositionType positions, const std::size_t begin, const std::size_t end)
 Build the linked cell list with a subset of particles.
 
template<class PositionType>
void build (PositionType positions, const std::size_t begin, const std::size_t end)
 Build the linked cell list with a subset of particles.
 
template<class PositionType>
void build (PositionType positions)
 Build the linked cell list with all particles.
 
void storeParticleBins ()
 Store the bin cell index for each binned particle.
 
auto getParticleBins () const
 Get the bin cell index for each binned particle.
 
KOKKOS_INLINE_FUNCTION auto getParticleBin (const int particle_index) const
 Get the bin cell index of the input particle.
 
KOKKOS_FUNCTION void operator() (const int i) const
 Determines which particles belong to bin i.
 
void update (const bool sorted)
 Updates the status of the list sorting.
 
KOKKOS_INLINE_FUNCTION auto sorted () const
 Returns whether the list has been sorted or not.
 
KOKKOS_INLINE_FUNCTION void getStencilCells (const int cell, int &imin, int &imax, int &jmin, int &jmax, int &kmin, int &kmax) const
 Get the cell indices for the stencil about cell.
 
KOKKOS_INLINE_FUNCTION auto getParticle (const int offset) const
 Get a candidate neighbor particle at a given binned offset.
 

Detailed Description

template<class MemorySpace, class Scalar = double>
class Cabana::LinkedCellList< MemorySpace, Scalar >

Data describing the bin sizes and offsets resulting from a binning operation on a 3d regular Cartesian grid.

Note
Defaults to double precision for backwards compatibility.

Constructor & Destructor Documentation

◆ LinkedCellList() [1/4]

template<class MemorySpace, class Scalar = double>
template<class PositionType>
Cabana::LinkedCellList< MemorySpace, Scalar >::LinkedCellList ( PositionType positions,
const Scalar grid_delta[3],
const Scalar grid_min[3],
const Scalar grid_max[3],
typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type * = 0 )
inline

Simple constructor.

Template Parameters
PositionTypeType for positions.
Parameters
positionsParticle positions.
grid_deltaGrid sizes in each cardinal direction.
grid_minGrid minimum value in each direction.
grid_maxGrid maximum value in each direction.

◆ LinkedCellList() [2/4]

template<class MemorySpace, class Scalar = double>
template<class PositionType>
Cabana::LinkedCellList< MemorySpace, Scalar >::LinkedCellList ( PositionType positions,
const std::size_t begin,
const std::size_t end,
const Scalar grid_delta[3],
const Scalar grid_min[3],
const Scalar grid_max[3],
typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type * = 0 )
inline

Partial range constructor.

Template Parameters
PositionTypeType for positions.
Parameters
positionsParticle positions.
beginThe beginning index of particles to bin or find neighbors for. Particles outside this range will NOT be considered as candidate neighbors.
endThe end index of particles to bin or find neighbors for. Particles outside this range will NOT be considered as candidate neighbors.
grid_deltaGrid sizes in each cardinal direction.
grid_minGrid minimum value in each direction.
grid_maxGrid maximum value in each direction.

◆ LinkedCellList() [3/4]

template<class MemorySpace, class Scalar = double>
template<class PositionType>
Cabana::LinkedCellList< MemorySpace, Scalar >::LinkedCellList ( PositionType positions,
const Scalar grid_delta[3],
const Scalar grid_min[3],
const Scalar grid_max[3],
const Scalar neighborhood_radius,
const Scalar cell_size_ratio = 1,
typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type * = 0 )
inline

Explicit stencil constructor.

Template Parameters
PositionTypeType for positions.
Parameters
positionsParticle positions.
grid_deltaGrid sizes in each cardinal direction.
grid_minGrid minimum value in each direction.
grid_maxGrid maximum value in each direction.
neighborhood_radiusRadius for neighbors.
cell_size_ratioRatio of the cell size to the neighborhood size.

◆ LinkedCellList() [4/4]

template<class MemorySpace, class Scalar = double>
template<class PositionType>
Cabana::LinkedCellList< MemorySpace, Scalar >::LinkedCellList ( PositionType positions,
const std::size_t begin,
const std::size_t end,
const Scalar grid_delta[3],
const Scalar grid_min[3],
const Scalar grid_max[3],
const Scalar neighborhood_radius,
const Scalar cell_size_ratio = 1,
typename std::enable_if<(is_slice< PositionType >::value||Kokkos::is_view< PositionType >::value), int >::type * = 0 )
inline

Explicit stencil and partial range constructor.

Template Parameters
PositionTypeType for positions.
Parameters
positionsParticle positions.
beginThe beginning index of particles to bin or find neighbors for. Particles outside this range will NOT be considered as candidate neighbors.
endThe end index of particles to bin or find neighbors for. Particles outside this range will NOT be considered as candidate neighbors.
grid_deltaGrid sizes in each cardinal direction.
grid_minGrid minimum value in each direction.
grid_maxGrid maximum value in each direction.
neighborhood_radiusRadius for neighbors.
cell_size_ratioRatio of the cell size to the neighborhood size.

Member Function Documentation

◆ binningData()

template<class MemorySpace, class Scalar = double>
BinningData< MemorySpace > Cabana::LinkedCellList< MemorySpace, Scalar >::binningData ( ) const
inline

Get the 1d bin data.

Returns
The 1d bin data.

◆ binOffset()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION size_type Cabana::LinkedCellList< MemorySpace, Scalar >::binOffset ( const int i,
const int j,
const int k ) const
inline

Given a bin get the particle index at which it sorts.

Parameters
iThe i bin index (x).
jThe j bin index (y).
kThe k bin index (z).
Returns
The starting particle index of the bin.

◆ binSize()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION int Cabana::LinkedCellList< MemorySpace, Scalar >::binSize ( const int i,
const int j,
const int k ) const
inline

Given a bin get the number of particles it contains.

Parameters
iThe i bin index (x).
jThe j bin index (y).
kThe k bin index (z).
Returns
The number of particles in the bin.

◆ build() [1/3]

template<class MemorySpace, class Scalar = double>
template<class ExecutionSpace, class PositionType>
void Cabana::LinkedCellList< MemorySpace, Scalar >::build ( ExecutionSpace ,
PositionType positions,
const std::size_t begin,
const std::size_t end )
inline

Build the linked cell list with a subset of particles.

Template Parameters
ExecutionSpaceKokkos execution space.
PositionTypeType for positions.
Parameters
positionsParticle positions.
beginThe beginning index of particles to bin or find neighbors for. Particles outside this range will NOT be considered as candidate neighbors.
endThe end index of particles to bin or find neighbors for. Particles outside this range will NOT be considered as candidate neighbors.

◆ build() [2/3]

template<class MemorySpace, class Scalar = double>
template<class PositionType>
void Cabana::LinkedCellList< MemorySpace, Scalar >::build ( PositionType positions)
inline

Build the linked cell list with all particles.

Template Parameters
PositionTypeType for positions.
Parameters
positionsParticle positions.

◆ build() [3/3]

template<class MemorySpace, class Scalar = double>
template<class PositionType>
void Cabana::LinkedCellList< MemorySpace, Scalar >::build ( PositionType positions,
const std::size_t begin,
const std::size_t end )
inline

Build the linked cell list with a subset of particles.

Template Parameters
PositionTypeType for positions.
Parameters
positionsParticle positions.
beginThe beginning index of particles to bin or find neighbors for. Particles outside this range will NOT be considered as candidate neighbors.
endThe end index of particles to bin or find neighbors for. Particles outside this range will NOT be considered as candidate neighbors.

◆ cardinalBinIndex()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION size_type Cabana::LinkedCellList< MemorySpace, Scalar >::cardinalBinIndex ( const int i,
const int j,
const int k ) const
inline

Given the ijk index of a bin get its cardinal index.

Parameters
iThe i bin index (x).
jThe j bin index (y).
kThe k bin index (z).
Returns
The cardinal bin index.

Note that the Kokkos sort orders the bins such that the i index moves the slowest and the k index mvoes the fastest.

◆ cellStencil()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION stencil_type Cabana::LinkedCellList< MemorySpace, Scalar >::cellStencil ( ) const
inline

Get the linked cell stencil.

Returns
Cell stencil.

◆ getParticle()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION auto Cabana::LinkedCellList< MemorySpace, Scalar >::getParticle ( const int offset) const
inline

Get a candidate neighbor particle at a given binned offset.

Parameters
offsetParticle offset in the binned layout.

◆ getParticleBins()

template<class MemorySpace, class Scalar = double>
auto Cabana::LinkedCellList< MemorySpace, Scalar >::getParticleBins ( ) const
inline

Get the bin cell index for each binned particle.

Note
This View only stores bins for particles which are binned.

◆ ijkBinIndex()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION void Cabana::LinkedCellList< MemorySpace, Scalar >::ijkBinIndex ( const int cardinal,
int & i,
int & j,
int & k ) const
inline

Given the cardinal index of a bin get its ijk indices.

Parameters
cardinalThe cardinal bin index.
iThe i bin index (x).
jThe j bin index (y).
kThe k bin index (z).

Note that the Kokkos sort orders the bins such that the i index moves the slowest and the k index mvoes the fastest.

◆ numBin()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION int Cabana::LinkedCellList< MemorySpace, Scalar >::numBin ( const int dim) const
inline

Get the number of bins in a given dimension.

Parameters
dimThe dimension to get the number of bins for.
Returns
The number of bins.

◆ operator()()

template<class MemorySpace, class Scalar = double>
KOKKOS_FUNCTION void Cabana::LinkedCellList< MemorySpace, Scalar >::operator() ( const int i) const
inline

Determines which particles belong to bin i.

Parameters
ithe cardinal bin index of the current bin

◆ permutation()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION size_type Cabana::LinkedCellList< MemorySpace, Scalar >::permutation ( const int particle_id) const
inline

Given a local particle id in the binned layout, get the id of the particle in the old (unbinned) layout.

Parameters
particle_idThe id of the particle in the binned layout.
Returns
The particle id in the old (unbinned) layout.

◆ totalBins()

template<class MemorySpace, class Scalar = double>
KOKKOS_INLINE_FUNCTION int Cabana::LinkedCellList< MemorySpace, Scalar >::totalBins ( ) const
inline

Get the total number of bins.

Returns
the total number of bins.

◆ update()

template<class MemorySpace, class Scalar = double>
void Cabana::LinkedCellList< MemorySpace, Scalar >::update ( const bool sorted)
inline

Updates the status of the list sorting.

Parameters
sortedBool that is true if the list has been sorted

The documentation for this class was generated from the following file: