Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana::Grid::GlobalParticleComm< MemorySpace, LocalGridType > Class Template Reference

Global particle communication based on the background grid. More...

#include <Cabana_Grid_GlobalParticleComm.hpp>

Public Types

using mesh_type = typename LocalGridType::mesh_type
 Mesh type.
 
using global_grid_type = Cabana::Grid::GlobalGrid<mesh_type>
 Global grid.
 
using memory_space = MemorySpace
 Kokkos memory space.
 
using corner_view_type
 Local boundary View type.
 
using destination_view_type = Kokkos::View<int*, memory_space>
 Particle destination ranks View type.
 
using rank_view_type
 Cartesian rank View type.
 
using host_rank_view_type
 Cartesian rank View type (host).
 

Public Member Functions

 GlobalParticleComm (const LocalGridType local_grid)
 Constructor.
 
template<class LocalMeshType>
void storeRanks (LocalMeshType local_mesh)
 Store local rank boundaries from the local mesh.
 
void scaleRanks ()
 Scale local rank boundaries based on double counting from MPI reduction.
 
template<std::size_t NSD = num_space_dim>
std::enable_if_t< 3==NSD, Kokkos::Array< int, num_space_dim > > getScaling ()
 Scaling factors due to double counting from MPI reduction.
 
template<std::size_t NSD = num_space_dim>
std::enable_if_t< 2==NSD, Kokkos::Array< int, num_space_dim > > getScaling ()
 Scaling factors due to double counting from MPI reduction.
 
template<class GlobalGridType, std::size_t NSD = num_space_dim>
std::enable_if_t< 3==NSD, void > copyRanks (GlobalGridType global_grid)
 Store all cartesian MPI ranks.
 
template<class GlobalGridType, std::size_t NSD = num_space_dim>
std::enable_if_t< 2==NSD, void > copyRanks (GlobalGridType global_grid)
 Store all cartesian MPI ranks.
 
template<class ExecutionSpace, class PositionType>
void build (ExecutionSpace exec_space, PositionType positions, const std::size_t begin, const std::size_t end)
 Bin particles across the global grid.
 
template<class ExecutionSpace, class PositionType>
void build (ExecutionSpace exec_space, PositionType positions)
 Bin particles across the global grid.
 
template<class PositionType>
void build (PositionType positions)
 Bin particles across the global grid.
 
template<class AoSoAType>
void migrate (MPI_Comm comm, AoSoAType &aosoa)
 Migrate particles to the correct rank.
 

Static Public Attributes

static constexpr std::size_t num_space_dim = LocalGridType::num_space_dim
 Spatial dimension.
 

Protected Attributes

int _rank_1d
 Current rank.
 
Kokkos::Array< int, num_space_dim_rank
 Current cartesian rank.
 
Kokkos::Array< int, num_space_dim_ranks_per_dim
 Total ranks per dimension.
 
corner_view_type _local_corners
 Local boundaries.
 
rank_view_type _global_ranks
 All cartesian ranks.
 
destination_view_type _destinations
 Particle destination ranks.
 

Detailed Description

template<class MemorySpace, class LocalGridType>
class Cabana::Grid::GlobalParticleComm< MemorySpace, LocalGridType >

Global particle communication based on the background grid.

Member Typedef Documentation

◆ corner_view_type

template<class MemorySpace, class LocalGridType>
using Cabana::Grid::GlobalParticleComm< MemorySpace, LocalGridType >::corner_view_type
Initial value:
Kokkos::View<double* [num_space_dim][2], memory_space>

Local boundary View type.

◆ host_rank_view_type

template<class MemorySpace, class LocalGridType>
using Cabana::Grid::GlobalParticleComm< MemorySpace, LocalGridType >::host_rank_view_type
Initial value:
Kokkos::View<int***, Kokkos::LayoutRight, Kokkos::HostSpace>

Cartesian rank View type (host).

◆ rank_view_type

template<class MemorySpace, class LocalGridType>
using Cabana::Grid::GlobalParticleComm< MemorySpace, LocalGridType >::rank_view_type
Initial value:
Kokkos::View<int***, Kokkos::LayoutRight, memory_space>

Cartesian rank View type.

Member Function Documentation

◆ build()

template<class MemorySpace, class LocalGridType>
template<class ExecutionSpace, class PositionType>
void Cabana::Grid::GlobalParticleComm< MemorySpace, LocalGridType >::build ( ExecutionSpace exec_space,
PositionType positions,
const std::size_t begin,
const std::size_t end )
inline

Bin particles across the global grid.

Because of MPI partitioning, this is not a perfect grid (as the Core LinkedCellList is), so we use binary search instead of direct 3d->1d indexing.


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