Cabana 0.8.0-dev
Loading...
Searching...
No Matches
Cabana_Grid_GlobalParticleComm.hpp
Go to the documentation of this file.
1/****************************************************************************
2 * Copyright (c) 2018-2023 by the Cabana authors *
3 * All rights reserved. *
4 * *
5 * This file is part of the Cabana library. Cabana is distributed under a *
6 * BSD 3-clause license. For the licensing terms see the LICENSE file in *
7 * the top-level directory. *
8 * *
9 * SPDX-License-Identifier: BSD-3-Clause *
10 ****************************************************************************/
11
16#ifndef CABANA_GRID_GLOBALPARTICLECOMM_HPP
17#define CABANA_GRID_GLOBALPARTICLECOMM_HPP
18
20#include <Cabana_Grid_Types.hpp>
21#include <Cabana_Tags.hpp>
22
24#include <Cabana_Slice.hpp>
25
26#include <memory>
27#include <stdexcept>
28#include <type_traits>
29
30namespace Cabana
31{
32namespace Grid
33{
34//---------------------------------------------------------------------------//
38template <class MemorySpace, class LocalGridType, class CommSpaceType = Mpi>
40{
41 public:
43 static constexpr std::size_t num_space_dim = LocalGridType::num_space_dim;
45 using mesh_type = typename LocalGridType::mesh_type;
49 using memory_space = MemorySpace;
51 using commspace_type = CommSpaceType;
52
55 Kokkos::View<double* [num_space_dim][2], memory_space>;
57 using destination_view_type = Kokkos::View<int*, memory_space>;
60 Kokkos::View<int***, Kokkos::LayoutRight, memory_space>;
63 Kokkos::View<int***, Kokkos::LayoutRight, Kokkos::HostSpace>;
64
66 GlobalParticleComm( const LocalGridType local_grid )
67 {
68 auto global_grid = local_grid.globalGrid();
70 Kokkos::ViewAllocateWithoutInitializing( "global_destination" ),
71 0 );
72
73 int max_ranks_per_dim = -1;
74 for ( std::size_t d = 0; d < num_space_dim; ++d )
75 {
76 _ranks_per_dim[d] = global_grid.dimNumBlock( d );
77 if ( _ranks_per_dim[d] > max_ranks_per_dim )
78 max_ranks_per_dim = _ranks_per_dim[d];
79 }
80 copyRanks( global_grid );
81
82 // Purposely using zero-init. Some entries unused in non-cubic
83 // decompositions.
85 corner_view_type( "local_mpi_boundaries", max_ranks_per_dim );
86
87 _rank_1d = global_grid.blockId();
88 for ( std::size_t d = 0; d < num_space_dim; ++d )
89 _rank[d] = global_grid.dimBlockId( d );
90
91 auto local_mesh = createLocalMesh<memory_space>( local_grid );
92 storeRanks( local_mesh );
93
94 // Update local boundaries from all ranks.
95 auto comm = global_grid.comm();
96 // TODO: Could use subcommunicators instead.
97 MPI_Allreduce( MPI_IN_PLACE, _local_corners.data(),
98 _local_corners.size(), MPI_DOUBLE, MPI_SUM, comm );
99
100 scaleRanks();
101 }
102
104 template <class LocalMeshType>
105 void storeRanks( LocalMeshType local_mesh )
106 {
107 auto local_corners = _local_corners;
108 auto rank = _rank;
109 auto store_corners = KOKKOS_LAMBDA( const std::size_t d )
110 {
111 local_corners( rank[d], d, 0 ) =
112 local_mesh.lowCorner( Cabana::Grid::Own(), d );
113 local_corners( rank[d], d, 1 ) =
114 local_mesh.highCorner( Cabana::Grid::Own(), d );
115 };
116 using exec_space = typename memory_space::execution_space;
117 Kokkos::RangePolicy<exec_space> policy( 0, num_space_dim );
118 Kokkos::parallel_for( "Cabana::Grid::GlobalParticleComm::storeCorners",
119 policy, store_corners );
120 Kokkos::fence();
121 }
122
125 {
126 auto scale = getScaling();
127
128 auto local_corners = _local_corners;
129 auto ranks_per_dim = _ranks_per_dim;
130 auto scale_corners = KOKKOS_LAMBDA( const std::size_t d )
131 {
132 for ( int r = 0; r < ranks_per_dim[d]; ++r )
133 {
134 local_corners( r, d, 0 ) /= scale[d];
135 local_corners( r, d, 1 ) /= scale[d];
136 }
137 };
138 using exec_space = typename memory_space::execution_space;
139 Kokkos::RangePolicy<exec_space> policy( 0, num_space_dim );
140 Kokkos::parallel_for( "Cabana::Grid::GlobalParticleComm::scaleCorners",
141 policy, scale_corners );
142 Kokkos::fence();
143 }
144
146 template <std::size_t NSD = num_space_dim>
147 std::enable_if_t<3 == NSD, Kokkos::Array<int, num_space_dim>> getScaling()
148 {
149 Kokkos::Array<int, num_space_dim> scale;
150 scale[0] = _ranks_per_dim[1] * _ranks_per_dim[2];
151 scale[1] = _ranks_per_dim[0] * _ranks_per_dim[2];
152 scale[2] = _ranks_per_dim[0] * _ranks_per_dim[1];
153 return scale;
154 }
155
157 template <std::size_t NSD = num_space_dim>
158 std::enable_if_t<2 == NSD, Kokkos::Array<int, num_space_dim>> getScaling()
159 {
160 Kokkos::Array<int, num_space_dim> scale;
161 scale[0] = _ranks_per_dim[1];
162 scale[1] = _ranks_per_dim[0];
163 return scale;
164 }
165
167 template <class GlobalGridType, std::size_t NSD = num_space_dim>
168 std::enable_if_t<3 == NSD, void> copyRanks( GlobalGridType global_grid )
169 {
170 host_rank_view_type global_ranks_host(
171 Kokkos::ViewAllocateWithoutInitializing( "ranks_host" ),
173 for ( int i = 0; i < _ranks_per_dim[0]; ++i )
174 for ( int j = 0; j < _ranks_per_dim[1]; ++j )
175 for ( int k = 0; k < _ranks_per_dim[2]; ++k )
176 // Not device accessible (uses MPI), so must be copied.
177 global_ranks_host( i, j, k ) =
178 global_grid.blockRank( i, j, k );
179
180 _global_ranks = Kokkos::create_mirror_view_and_copy(
181 memory_space(), global_ranks_host );
182 }
183
185 template <class GlobalGridType, std::size_t NSD = num_space_dim>
186 std::enable_if_t<2 == NSD, void> copyRanks( GlobalGridType global_grid )
187 {
188 // Storing as 3d for convenience.
189 host_rank_view_type global_ranks_host(
190 Kokkos::ViewAllocateWithoutInitializing( "ranks_host" ),
191 _ranks_per_dim[0], _ranks_per_dim[1], 1 );
192 for ( int i = 0; i < _ranks_per_dim[0]; ++i )
193 for ( int j = 0; j < _ranks_per_dim[1]; ++j )
194 // Not device accessible (uses MPI), so must be copied.
195 global_ranks_host( i, j, 0 ) = global_grid.blockRank( i, j );
196
197 _global_ranks = Kokkos::create_mirror_view_and_copy(
198 memory_space(), global_ranks_host );
199 }
200
208 template <class ExecutionSpace, class PositionType>
209 void build( ExecutionSpace exec_space, PositionType positions,
210 const std::size_t begin, const std::size_t end )
211 {
212 Kokkos::Profiling::ScopedRegion region(
213 "Cabana::Grid::GlobalParticleComm::build" );
214
216 assert( end >= begin );
217 assert( end <= positions.size() );
218
219 // Must match the size of all particles, even if some can be ignored in
220 // this search.
221 Kokkos::resize( _destinations, positions.size() );
222 // Start with everything staying on this rank.
223 Kokkos::deep_copy( _destinations, _rank_1d );
224
225 // Local copies for lambda capture.
226 auto local_corners = _local_corners;
227 auto ranks_per_dim = _ranks_per_dim;
228 auto destinations = _destinations;
229 auto global_ranks = _global_ranks;
230 auto build_migrate = KOKKOS_LAMBDA( const std::size_t p )
231 {
232 // This is not num_space_dim because global_ranks is always rank-3
233 // out of convenience.
234 int ijk[3] = { 0, 0, 0 };
235
236 // Find the rank this particle should be moved to.
237 for ( std::size_t d = 0; d < num_space_dim; ++d )
238 {
239 int min = 0;
240 int max = ranks_per_dim[d];
241
242 // Check if outside the box in this dimension.
243 if ( ( positions( p, d ) < local_corners( min, d, 0 ) ) ||
244 ( positions( p, d ) > local_corners( max - 1, d, 1 ) ) )
245 destinations( p ) = -1;
246
247 // Do a binary search for this particle in this dimension.
248 while ( max - min > 1 )
249 {
250 int center = Kokkos::floor( ( max + min ) / 2.0 );
251 if ( positions( p, d ) < local_corners( center, d, 0 ) )
252 max = center;
253 else
254 min = center;
255 }
256 ijk[d] = min;
257 }
258 // Keep the destination rank for eventual migration.
259 destinations( p ) = global_ranks( ijk[0], ijk[1], ijk[2] );
260 };
261
262 Kokkos::RangePolicy<ExecutionSpace> policy( exec_space, begin, end );
263 Kokkos::parallel_for( "Cabana::Grid::GlobalParticleComm::build", policy,
264 build_migrate );
265 Kokkos::fence();
266 }
267
269 template <class ExecutionSpace, class PositionType>
270 void build( ExecutionSpace exec_space, PositionType positions )
271 {
272 build( exec_space, positions, 0, positions.size() );
273 }
274
276 template <class PositionType>
277 void build( PositionType positions )
278 {
279 using execution_space = typename memory_space::execution_space;
280 // TODO: enable views.
281 build( execution_space{}, positions, 0, positions.size() );
282 }
283
285 template <class AoSoAType>
286 void migrate( MPI_Comm comm, AoSoAType& aosoa )
287 {
289 comm, _destinations );
290 Cabana::migrate( distributor, aosoa );
291 }
292
293 protected:
297 Kokkos::Array<int, num_space_dim> _rank;
299 Kokkos::Array<int, num_space_dim> _ranks_per_dim;
306};
307
312template <class MemorySpace, class LocalGridType, class CommSpaceType = Mpi>
313auto createGlobalParticleComm( const LocalGridType& local_grid )
314{
315 return std::make_shared<
317 local_grid );
318}
319
320//---------------------------------------------------------------------------//
321
322} // namespace Grid
323} // namespace Cabana
324
325#endif // end CABANA_GRID_GLOBALPARTICLECOMM_HPP
Multi-node particle redistribution.
auto createGlobalParticleComm(const LocalGridType &local_grid)
Create global linked cell binning.
Definition Cabana_Grid_GlobalParticleComm.hpp:313
LocalMesh< MemorySpace, MeshType > createLocalMesh(const LocalGrid< MeshType > &local_grid)
Creation function for local mesh.
Definition Cabana_Grid_LocalMesh.hpp:779
Grid type tags.
Slice a single particle property from an AoSoA.
Type tags used in Cabana.
A communication plan for migrating data from one uniquely-owned decomposition to another uniquely own...
Definition Cabana_Distributor.hpp:64
Global logical grid.
Definition Cabana_Grid_GlobalGrid.hpp:39
Global particle communication based on the background grid.
Definition Cabana_Grid_GlobalParticleComm.hpp:40
Kokkos::View< double *[num_space_dim][2], memory_space > corner_view_type
Local boundary View type.
Definition Cabana_Grid_GlobalParticleComm.hpp:54
std::enable_if_t< 2==NSD, Kokkos::Array< int, num_space_dim > > getScaling()
Scaling factors due to double counting from MPI reduction.
Definition Cabana_Grid_GlobalParticleComm.hpp:158
void build(ExecutionSpace exec_space, PositionType positions, const std::size_t begin, const std::size_t end)
Bin particles across the global grid.
Definition Cabana_Grid_GlobalParticleComm.hpp:209
int _rank_1d
Current rank.
Definition Cabana_Grid_GlobalParticleComm.hpp:295
corner_view_type _local_corners
Local boundaries.
Definition Cabana_Grid_GlobalParticleComm.hpp:301
Cabana::Grid::GlobalGrid< mesh_type > global_grid_type
Global grid.
Definition Cabana_Grid_GlobalParticleComm.hpp:47
void migrate(MPI_Comm comm, AoSoAType &aosoa)
Migrate particles to the correct rank.
Definition Cabana_Grid_GlobalParticleComm.hpp:286
std::enable_if_t< 3==NSD, Kokkos::Array< int, num_space_dim > > getScaling()
Scaling factors due to double counting from MPI reduction.
Definition Cabana_Grid_GlobalParticleComm.hpp:147
std::enable_if_t< 2==NSD, void > copyRanks(GlobalGridType global_grid)
Store all cartesian MPI ranks.
Definition Cabana_Grid_GlobalParticleComm.hpp:186
std::enable_if_t< 3==NSD, void > copyRanks(GlobalGridType global_grid)
Store all cartesian MPI ranks.
Definition Cabana_Grid_GlobalParticleComm.hpp:168
destination_view_type _destinations
Particle destination ranks.
Definition Cabana_Grid_GlobalParticleComm.hpp:305
void storeRanks(LocalMeshType local_mesh)
Store local rank boundaries from the local mesh.
Definition Cabana_Grid_GlobalParticleComm.hpp:105
static constexpr std::size_t num_space_dim
Spatial dimension.
Definition Cabana_Grid_GlobalParticleComm.hpp:43
Kokkos::Array< int, num_space_dim > _ranks_per_dim
Total ranks per dimension.
Definition Cabana_Grid_GlobalParticleComm.hpp:299
CommSpaceType commspace_type
Communication backend.
Definition Cabana_Grid_GlobalParticleComm.hpp:51
void scaleRanks()
Scale local rank boundaries based on double counting from MPI reduction.
Definition Cabana_Grid_GlobalParticleComm.hpp:124
Kokkos::Array< int, num_space_dim > _rank
Current cartesian rank.
Definition Cabana_Grid_GlobalParticleComm.hpp:297
Kokkos::View< int ***, Kokkos::LayoutRight, Kokkos::HostSpace > host_rank_view_type
Cartesian rank View type (host).
Definition Cabana_Grid_GlobalParticleComm.hpp:62
void build(ExecutionSpace exec_space, PositionType positions)
Bin particles across the global grid.
Definition Cabana_Grid_GlobalParticleComm.hpp:270
Kokkos::View< int *, memory_space > destination_view_type
Particle destination ranks View type.
Definition Cabana_Grid_GlobalParticleComm.hpp:57
Kokkos::View< int ***, Kokkos::LayoutRight, memory_space > rank_view_type
Cartesian rank View type.
Definition Cabana_Grid_GlobalParticleComm.hpp:59
rank_view_type _global_ranks
All cartesian ranks.
Definition Cabana_Grid_GlobalParticleComm.hpp:303
void build(PositionType positions)
Bin particles across the global grid.
Definition Cabana_Grid_GlobalParticleComm.hpp:277
typename LocalGridType::mesh_type mesh_type
Mesh type.
Definition Cabana_Grid_GlobalParticleComm.hpp:45
MemorySpace memory_space
Kokkos memory space.
Definition Cabana_Grid_GlobalParticleComm.hpp:49
GlobalParticleComm(const LocalGridType local_grid)
Constructor.
Definition Cabana_Grid_GlobalParticleComm.hpp:66
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
void migrate(ExecutionSpace exec_space, const Distributor_t &distributor, const AoSoA_t &src, AoSoA_t &dst, typename std::enable_if<(is_distributor< Distributor_t >::value &&is_aosoa< AoSoA_t >::value), int >::type *=0)
Synchronously migrate data between two different decompositions using the distributor forward communi...
Definition Cabana_Distributor.hpp:203
Owned decomposition tag.
Definition Cabana_Grid_Types.hpp:190
Definition Cabana_Types.hpp:88