Loading...
Searching...
No Matches
communication_plan.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <algorithm>
4#include <map>
5#include <tuple>
6#include <vector>
7
8#include "dense/vec.hpp"
9#include "grid/grid_types.hpp"
13#include "util/timer.hpp"
14
16
17
18// --------------------------------------------------------------------------------------
19// Reusable, precomputed plan
20// --------------------------------------------------------------------------------------
21//
22// Goal: avoid rebuilding send/recv pair lists, sorting, chunk layout, and per-rank buffer sizes
23// on every halo exchange. We do that once in the ctor, then `exchange_and_reduce(...)` just runs
24// the hot path: local copies, pack, fence, post isends/irecvs, wait, scatter, unpack.
25//
26// Notes:
27// - This keeps your existing per-boundary recv buffers intact and still used by unpack.
28// - It keeps your per-rank aggregation optimization.
29// - It does not depend on send_buffers_ (your argument is unused currently anyway).
30//
31template < class GridDataType >
33{
34 public:
35 using ScalarType = typename GridDataType::value_type;
36 static constexpr int VecDim = grid::grid_data_vec_dim< GridDataType >();
37 using memory_space = typename GridDataType::memory_space;
38 using rank_buffer_view = Kokkos::View< ScalarType*, memory_space >;
39
40 explicit ShellBoundaryCommPlan( const grid::shell::DistributedDomain& domain, bool enable_local_comm = true )
41 : domain_( &domain ), enable_local_comm_( enable_local_comm )
42 {
43 build_plan_();
44 allocate_rank_buffers_();
45 }
46
47 // Call this each timestep/iteration.
49 const GridDataType& data,
52 {
53 util::Timer timer_all( "shell_boundary_exchange_and_reduce" );
54
55 post_irecvs_();
56
57 local_comm_copy_into_recv_buffers_( data, boundary_recv_buffers );
58
59 pack_remote_sends_( data );
60
61 post_isends_();
62
63 wait_all_();
64
65 scatter_recvs_into_boundary_buffers_( boundary_recv_buffers );
66
67 unpack_and_reduce_( data, boundary_recv_buffers, reduction );
68 }
69
70 // Optional: if domain topology changes (rare), rebuild everything.
71 void rebuild()
72 {
73 build_plan_();
74 allocate_rank_buffers_();
75 }
76
77 private:
78 struct SendRecvPair
79 {
80 int boundary_type = -1; // 0 vertex, 1 edge, 2 face
81 mpi::MPIRank local_rank;
82 grid::shell::SubdomainInfo local_subdomain;
83 int local_subdomain_boundary;
84 int local_subdomain_id;
85
86 mpi::MPIRank neighbor_rank;
87 grid::shell::SubdomainInfo neighbor_subdomain;
88 int neighbor_subdomain_boundary;
89 };
90
91 struct ChunkInfo
92 {
93 SendRecvPair pair;
94 int offset = 0; // in scalars
95 int size = 0; // in scalars
96 };
97
98 // --------------------------
99 // Plan build / layout
100 // --------------------------
101 int piece_num_scalars_( const SendRecvPair& p ) const
102 {
103 const auto& domain = *domain_;
104
105 if ( p.boundary_type == 0 )
106 {
107 return VecDim;
108 }
109 else if ( p.boundary_type == 1 )
110 {
111 const auto local_edge_boundary = static_cast< grid::BoundaryEdge >( p.local_subdomain_boundary );
112 const int n_nodes = grid::is_edge_boundary_radial( local_edge_boundary ) ?
114 domain.domain_info().subdomain_num_nodes_per_side_laterally();
115 return n_nodes * VecDim;
116 }
117 else if ( p.boundary_type == 2 )
118 {
119 const auto local_face_boundary = static_cast< grid::BoundaryFace >( p.local_subdomain_boundary );
120 const int ni = domain.domain_info().subdomain_num_nodes_per_side_laterally();
121 const int nj = grid::is_face_boundary_normal_to_radial_direction( local_face_boundary ) ?
122 domain.domain_info().subdomain_num_nodes_per_side_laterally() :
123 domain.domain_info().subdomain_num_nodes_radially();
124 return ni * nj * VecDim;
125 }
126 Kokkos::abort( "Unknown boundary type" );
127 return 0;
128 }
129
130 void build_plan_()
131 {
132 util::Timer timer( "ShellBoundaryCommPlan::build_plan" );
133
134 const auto& domain = *domain_;
135
136 send_recv_pairs_.clear();
137 send_recv_pairs_.reserve( 1024 );
138
139 // Build the full (unsorted) pair list once.
140 for ( const auto& [local_subdomain_info, idx_and_neighborhood] : domain.subdomains() )
141 {
142 const auto& [local_subdomain_id, neighborhood] = idx_and_neighborhood;
143
144 for ( const auto& [local_vertex_boundary, neighbors] : neighborhood.neighborhood_vertex() )
145 {
146 for ( const auto& neighbor : neighbors )
147 {
148 const auto& [neighbor_subdomain_info, neighbor_local_boundary, neighbor_rank] = neighbor;
149 send_recv_pairs_.push_back( SendRecvPair{
150 .boundary_type = 0,
151 .local_rank = mpi::rank(),
152 .local_subdomain = local_subdomain_info,
153 .local_subdomain_boundary = static_cast< int >( local_vertex_boundary ),
154 .local_subdomain_id = local_subdomain_id,
155 .neighbor_rank = neighbor_rank,
156 .neighbor_subdomain = neighbor_subdomain_info,
157 .neighbor_subdomain_boundary = static_cast< int >( neighbor_local_boundary ) } );
158 }
159 }
160
161 for ( const auto& [local_edge_boundary, neighbors] : neighborhood.neighborhood_edge() )
162 {
163 for ( const auto& neighbor : neighbors )
164 {
165 const auto& [neighbor_subdomain_info, neighbor_local_boundary, _, neighbor_rank] = neighbor;
166 send_recv_pairs_.push_back( SendRecvPair{
167 .boundary_type = 1,
168 .local_rank = mpi::rank(),
169 .local_subdomain = local_subdomain_info,
170 .local_subdomain_boundary = static_cast< int >( local_edge_boundary ),
171 .local_subdomain_id = local_subdomain_id,
172 .neighbor_rank = neighbor_rank,
173 .neighbor_subdomain = neighbor_subdomain_info,
174 .neighbor_subdomain_boundary = static_cast< int >( neighbor_local_boundary ) } );
175 }
176 }
177
178 for ( const auto& [local_face_boundary, neighbor] : neighborhood.neighborhood_face() )
179 {
180 const auto& [neighbor_subdomain_info, neighbor_local_boundary, _, neighbor_rank] = neighbor;
181 send_recv_pairs_.push_back( SendRecvPair{
182 .boundary_type = 2,
183 .local_rank = mpi::rank(),
184 .local_subdomain = local_subdomain_info,
185 .local_subdomain_boundary = static_cast< int >( local_face_boundary ),
186 .local_subdomain_id = local_subdomain_id,
187 .neighbor_rank = neighbor_rank,
188 .neighbor_subdomain = neighbor_subdomain_info,
189 .neighbor_subdomain_boundary = static_cast< int >( neighbor_local_boundary ) } );
190 }
191 }
192
193 // Precompute local-comm subset (fixed list).
194 local_pairs_.clear();
195 local_pairs_.reserve( send_recv_pairs_.size() );
196 for ( const auto& p : send_recv_pairs_ )
197 {
198 if ( enable_local_comm_ && p.local_rank == p.neighbor_rank )
199 local_pairs_.push_back( p );
200 }
201
202 // SEND layout (sorted and chunked per rank, remote only)
203 {
204 auto send_pairs = send_recv_pairs_;
205 std::sort( send_pairs.begin(), send_pairs.end(), []( const SendRecvPair& a, const SendRecvPair& b ) {
206 if ( a.boundary_type != b.boundary_type ) return a.boundary_type < b.boundary_type;
207 if ( a.local_subdomain != b.local_subdomain ) return a.local_subdomain < b.local_subdomain;
208 if ( a.local_subdomain_boundary != b.local_subdomain_boundary )
209 return a.local_subdomain_boundary < b.local_subdomain_boundary;
210 if ( a.neighbor_subdomain != b.neighbor_subdomain ) return a.neighbor_subdomain < b.neighbor_subdomain;
211 return a.neighbor_subdomain_boundary < b.neighbor_subdomain_boundary;
212 } );
213
214 send_chunks_by_rank_.clear();
215 send_total_by_rank_.clear();
216
217 for ( const auto& p : send_pairs )
218 {
219 if ( enable_local_comm_ && p.local_rank == p.neighbor_rank )
220 continue;
221
222 const int sz = piece_num_scalars_( p );
223 auto& chunks = send_chunks_by_rank_[p.neighbor_rank];
224
225 const int off = send_total_by_rank_[p.neighbor_rank];
226 send_total_by_rank_[p.neighbor_rank] += sz;
227
228 chunks.push_back( ChunkInfo{ .pair = p, .offset = off, .size = sz } );
229 }
230 }
231
232 // RECV layout (sorted and chunked per rank, remote only)
233 {
234 auto recv_pairs = send_recv_pairs_;
235 std::sort( recv_pairs.begin(), recv_pairs.end(), []( const SendRecvPair& a, const SendRecvPair& b ) {
236 if ( a.boundary_type != b.boundary_type ) return a.boundary_type < b.boundary_type;
237 if ( a.neighbor_subdomain != b.neighbor_subdomain ) return a.neighbor_subdomain < b.neighbor_subdomain;
238 if ( a.neighbor_subdomain_boundary != b.neighbor_subdomain_boundary )
239 return a.neighbor_subdomain_boundary < b.neighbor_subdomain_boundary;
240 if ( a.local_subdomain != b.local_subdomain ) return a.local_subdomain < b.local_subdomain;
241 return a.local_subdomain_boundary < b.local_subdomain_boundary;
242 } );
243
244 recv_chunks_by_rank_.clear();
245 recv_total_by_rank_.clear();
246
247 for ( const auto& p : recv_pairs )
248 {
249 if ( enable_local_comm_ && p.local_rank == p.neighbor_rank )
250 continue;
251
252 const int sz = piece_num_scalars_( p );
253 auto& chunks = recv_chunks_by_rank_[p.neighbor_rank];
254
255 const int off = recv_total_by_rank_[p.neighbor_rank];
256 recv_total_by_rank_[p.neighbor_rank] += sz;
257
258 chunks.push_back( ChunkInfo{ .pair = p, .offset = off, .size = sz } );
259 }
260 }
261 }
262
263 void allocate_rank_buffers_()
264 {
265 util::Timer timer( "ShellBoundaryCommPlan::allocate_rank_buffers" );
266
267 send_rank_buffers_.clear();
268 recv_rank_buffers_.clear();
269
270 for ( const auto& [rank, total] : send_total_by_rank_ )
271 {
272 if ( total > 0 )
273 send_rank_buffers_[rank] = rank_buffer_view( "rank_send_buffer", total );
274 }
275 for ( const auto& [rank, total] : recv_total_by_rank_ )
276 {
277 if ( total > 0 )
278 recv_rank_buffers_[rank] = rank_buffer_view( "rank_recv_buffer", total );
279 }
280
281 data_send_requests_.resize( send_rank_buffers_.size() );
282 data_recv_requests_.resize( recv_rank_buffers_.size() );
283 }
284
285 // --------------------------
286 // Hot path
287 // --------------------------
288 void post_irecvs_() const
289 {
290 util::Timer timer( "ShellBoundaryCommPlan::post_irecvs" );
291
292 int i = 0;
293 for ( const auto& [rank, buf] : recv_rank_buffers_ )
294 {
295 const int total_sz = static_cast< int >( buf.extent( 0 ) );
296 MPI_Irecv(
297 buf.data(),
298 total_sz,
299 mpi::mpi_datatype< ScalarType >(),
300 rank,
302 MPI_COMM_WORLD,
303 &data_recv_requests_[i] );
304 ++i;
305 }
306 recv_req_count_ = i;
307 }
308
309 void local_comm_copy_into_recv_buffers_(
310 const GridDataType& data,
312 {
313 util::Timer timer( "ShellBoundaryCommPlan::local_comm" );
314
315 const auto& domain = *domain_;
316
317 for ( const auto& p : local_pairs_ )
318 {
319 if ( !domain.subdomains().contains( p.neighbor_subdomain ) )
320 Kokkos::abort( "Subdomain not found locally - but it should be there..." );
321
322 const auto neighbor_subdomain_id = std::get< 0 >( domain.subdomains().at( p.neighbor_subdomain ) );
323
324 if ( p.boundary_type == 0 )
325 {
326 auto& recv_buf = boundary_recv_buffers.buffer_vertex(
327 p.local_subdomain,
328 static_cast< grid::BoundaryVertex >( p.local_subdomain_boundary ),
329 p.neighbor_subdomain,
330 static_cast< grid::BoundaryVertex >( p.neighbor_subdomain_boundary ) );
331
332 copy_to_buffer<VecDim>(
333 recv_buf,
334 data,
335 neighbor_subdomain_id,
336 static_cast< grid::BoundaryVertex >( p.neighbor_subdomain_boundary ) );
337 }
338 else if ( p.boundary_type == 1 )
339 {
340 auto& recv_buf = boundary_recv_buffers.buffer_edge(
341 p.local_subdomain,
342 static_cast< grid::BoundaryEdge >( p.local_subdomain_boundary ),
343 p.neighbor_subdomain,
344 static_cast< grid::BoundaryEdge >( p.neighbor_subdomain_boundary ) );
345
346 copy_to_buffer<VecDim>(
347 recv_buf,
348 data,
349 neighbor_subdomain_id,
350 static_cast< grid::BoundaryEdge >( p.neighbor_subdomain_boundary ) );
351 }
352 else if ( p.boundary_type == 2 )
353 {
354 auto& recv_buf = boundary_recv_buffers.buffer_face(
355 p.local_subdomain,
356 static_cast< grid::BoundaryFace >( p.local_subdomain_boundary ),
357 p.neighbor_subdomain,
358 static_cast< grid::BoundaryFace >( p.neighbor_subdomain_boundary ) );
359
360 copy_to_buffer<VecDim>(
361 recv_buf,
362 data,
363 neighbor_subdomain_id,
364 static_cast< grid::BoundaryFace >( p.neighbor_subdomain_boundary ) );
365 }
366 else
367 {
368 Kokkos::abort( "Unknown boundary type" );
369 }
370 }
371 }
372
373 void pack_remote_sends_( const GridDataType& data ) const
374 {
375 util::Timer timer( "ShellBoundaryCommPlan::pack_remote" );
376
377 const auto& domain = *domain_;
378
379 for ( const auto& [rank, chunks] : send_chunks_by_rank_ )
380 {
381 auto& rank_buf = send_rank_buffers_.at( rank );
382
383 for ( const auto& ch : chunks )
384 {
385 const auto& p = ch.pair;
386 ScalarType* base_ptr = rank_buf.data() + ch.offset;
387
388 if ( p.boundary_type == 0 )
389 {
390 using BufT = grid::Grid0DDataVec< ScalarType, VecDim >;
391 auto unmanaged = detail::make_unmanaged_like< BufT >( base_ptr );
392
393 copy_to_buffer<VecDim>(
394 unmanaged,
395 data,
396 p.local_subdomain_id,
397 static_cast< grid::BoundaryVertex >( p.local_subdomain_boundary ) );
398 }
399 else if ( p.boundary_type == 1 )
400 {
401 using BufT = grid::Grid1DDataVec< ScalarType, VecDim >;
402 const auto local_edge_boundary = static_cast< grid::BoundaryEdge >( p.local_subdomain_boundary );
403 const int n_nodes = grid::is_edge_boundary_radial( local_edge_boundary ) ?
404 domain.domain_info().subdomain_num_nodes_radially() :
405 domain.domain_info().subdomain_num_nodes_per_side_laterally();
406
407 auto unmanaged = detail::make_unmanaged_like< BufT >( base_ptr, n_nodes );
408 copy_to_buffer<VecDim>( unmanaged, data, p.local_subdomain_id, local_edge_boundary );
409 }
410 else if ( p.boundary_type == 2 )
411 {
412 using BufT = grid::Grid2DDataVec< ScalarType, VecDim >;
413 const auto local_face_boundary = static_cast< grid::BoundaryFace >( p.local_subdomain_boundary );
414 const int ni = domain.domain_info().subdomain_num_nodes_per_side_laterally();
415 const int nj = grid::is_face_boundary_normal_to_radial_direction( local_face_boundary ) ?
416 domain.domain_info().subdomain_num_nodes_per_side_laterally() :
417 domain.domain_info().subdomain_num_nodes_radially();
418
419 auto unmanaged = detail::make_unmanaged_like< BufT >( base_ptr, ni, nj );
420 copy_to_buffer<VecDim>( unmanaged, data, p.local_subdomain_id, local_face_boundary );
421 }
422 else
423 {
424 Kokkos::abort( "Unknown boundary type" );
425 }
426 }
427 }
428
429 Kokkos::fence( "pack_rank_send_buffers" );
430 }
431
432 void post_isends_() const
433 {
434 util::Timer timer( "ShellBoundaryCommPlan::post_isends" );
435
436 int i = 0;
437 for ( const auto& [rank, buf] : send_rank_buffers_ )
438 {
439 const int total_sz = static_cast< int >( buf.extent( 0 ) );
440 MPI_Isend(
441 buf.data(),
442 total_sz,
443 mpi::mpi_datatype< ScalarType >(),
444 rank,
446 MPI_COMM_WORLD,
447 &data_send_requests_[i] );
448 ++i;
449 }
450 send_req_count_ = i;
451 }
452
453 void wait_all_() const
454 {
455 util::Timer timer( "ShellBoundaryCommPlan::waitall" );
456
457 if ( send_req_count_ > 0 )
458 MPI_Waitall( send_req_count_, data_send_requests_.data(), MPI_STATUSES_IGNORE );
459 if ( recv_req_count_ > 0 )
460 MPI_Waitall( recv_req_count_, data_recv_requests_.data(), MPI_STATUSES_IGNORE );
461 }
462
463 void scatter_recvs_into_boundary_buffers_(
465 {
466 util::Timer timer( "ShellBoundaryCommPlan::scatter_recvs" );
467
468 const auto& domain = *domain_;
469
470 for ( const auto& [rank, chunks] : recv_chunks_by_rank_ )
471 {
472 auto& rank_buf = recv_rank_buffers_.at( rank );
473
474 for ( const auto& ch : chunks )
475 {
476 const auto& p = ch.pair;
477 ScalarType* base_ptr = rank_buf.data() + ch.offset;
478
479 if ( p.boundary_type == 0 )
480 {
481 using BufT = grid::Grid0DDataVec< ScalarType, VecDim >;
482 auto unmanaged = detail::make_unmanaged_like< BufT >( base_ptr );
483
484 auto& recv_buf = boundary_recv_buffers.buffer_vertex(
485 p.local_subdomain,
486 static_cast< grid::BoundaryVertex >( p.local_subdomain_boundary ),
487 p.neighbor_subdomain,
488 static_cast< grid::BoundaryVertex >( p.neighbor_subdomain_boundary ) );
489
490 Kokkos::deep_copy( recv_buf, unmanaged );
491 }
492 else if ( p.boundary_type == 1 )
493 {
494 using BufT = grid::Grid1DDataVec< ScalarType, VecDim >;
495
496 const auto local_edge_boundary = static_cast< grid::BoundaryEdge >( p.local_subdomain_boundary );
497 const int n_nodes = grid::is_edge_boundary_radial( local_edge_boundary ) ?
498 domain.domain_info().subdomain_num_nodes_radially() :
499 domain.domain_info().subdomain_num_nodes_per_side_laterally();
500
501 auto unmanaged = detail::make_unmanaged_like< BufT >( base_ptr, n_nodes );
502
503 auto& recv_buf = boundary_recv_buffers.buffer_edge(
504 p.local_subdomain,
505 static_cast< grid::BoundaryEdge >( p.local_subdomain_boundary ),
506 p.neighbor_subdomain,
507 static_cast< grid::BoundaryEdge >( p.neighbor_subdomain_boundary ) );
508
509 Kokkos::deep_copy( recv_buf, unmanaged );
510 }
511 else if ( p.boundary_type == 2 )
512 {
513 using BufT = grid::Grid2DDataVec< ScalarType, VecDim >;
514
515 const auto local_face_boundary = static_cast< grid::BoundaryFace >( p.local_subdomain_boundary );
516 const int ni = domain.domain_info().subdomain_num_nodes_per_side_laterally();
517 const int nj = grid::is_face_boundary_normal_to_radial_direction( local_face_boundary ) ?
518 domain.domain_info().subdomain_num_nodes_per_side_laterally() :
519 domain.domain_info().subdomain_num_nodes_radially();
520
521 auto unmanaged = detail::make_unmanaged_like< BufT >( base_ptr, ni, nj );
522
523 auto& recv_buf = boundary_recv_buffers.buffer_face(
524 p.local_subdomain,
525 static_cast< grid::BoundaryFace >( p.local_subdomain_boundary ),
526 p.neighbor_subdomain,
527 static_cast< grid::BoundaryFace >( p.neighbor_subdomain_boundary ) );
528
529 Kokkos::deep_copy( recv_buf, unmanaged );
530 }
531 else
532 {
533 Kokkos::abort( "Unknown boundary type" );
534 }
535 }
536 }
537
538 Kokkos::fence( "scatter_rank_recv_buffers" );
539 }
540
541 void unpack_and_reduce_(
542 const GridDataType& data,
544 CommunicationReduction reduction ) const
545 {
546 util::Timer timer( "ShellBoundaryCommPlan::unpack_and_reduce" );
547
548 const auto& domain = *domain_;
549
550 for ( const auto& [local_subdomain_info, idx_and_neighborhood] : domain.subdomains() )
551 {
552 const auto& [local_subdomain_id, neighborhood] = idx_and_neighborhood;
553
554 for ( const auto& [local_vertex_boundary, neighbors] : neighborhood.neighborhood_vertex() )
555 {
556 for ( const auto& neighbor : neighbors )
557 {
558 const auto& [neighbor_subdomain_info, neighbor_local_boundary, neighbor_rank] = neighbor;
559
560 auto recv_buffer = boundary_recv_buffers.buffer_vertex(
561 local_subdomain_info, local_vertex_boundary, neighbor_subdomain_info, neighbor_local_boundary );
562
564 recv_buffer, data, local_subdomain_id, local_vertex_boundary, reduction );
565 }
566 }
567
568 for ( const auto& [local_edge_boundary, neighbors] : neighborhood.neighborhood_edge() )
569 {
570 for ( const auto& neighbor : neighbors )
571 {
572 const auto& [neighbor_subdomain_info, neighbor_local_boundary, boundary_direction, neighbor_rank] =
573 neighbor;
574
575 auto recv_buffer = boundary_recv_buffers.buffer_edge(
576 local_subdomain_info, local_edge_boundary, neighbor_subdomain_info, neighbor_local_boundary );
577
579 recv_buffer, data, local_subdomain_id, local_edge_boundary, boundary_direction, reduction );
580 }
581 }
582
583 for ( const auto& [local_face_boundary, neighbor] : neighborhood.neighborhood_face() )
584 {
585 const auto& [neighbor_subdomain_info, neighbor_local_boundary, boundary_directions, neighbor_rank] =
586 neighbor;
587
588 auto recv_buffer = boundary_recv_buffers.buffer_face(
589 local_subdomain_info, local_face_boundary, neighbor_subdomain_info, neighbor_local_boundary );
590
592 recv_buffer, data, local_subdomain_id, local_face_boundary, boundary_directions, reduction );
593 }
594 }
595
596 Kokkos::fence();
597 }
598
599 private:
600 const grid::shell::DistributedDomain* domain_ = nullptr;
601 bool enable_local_comm_ = true;
602
603 // Precomputed full list
604 std::vector< SendRecvPair > send_recv_pairs_;
605
606 // Precomputed local-only subset
607 std::vector< SendRecvPair > local_pairs_;
608
609 // Precomputed rank aggregation layouts
610 std::map< mpi::MPIRank, std::vector< ChunkInfo > > send_chunks_by_rank_;
611 std::map< mpi::MPIRank, std::vector< ChunkInfo > > recv_chunks_by_rank_;
612 std::map< mpi::MPIRank, int > send_total_by_rank_;
613 std::map< mpi::MPIRank, int > recv_total_by_rank_;
614
615 // Reused rank buffers
616 mutable std::map< mpi::MPIRank, rank_buffer_view > send_rank_buffers_;
617 mutable std::map< mpi::MPIRank, rank_buffer_view > recv_rank_buffers_;
618
619 // Reused request storage
620 mutable std::vector< MPI_Request > data_send_requests_;
621 mutable std::vector< MPI_Request > data_recv_requests_;
622 mutable int send_req_count_ = 0;
623 mutable int recv_req_count_ = 0;
624};
625
626// --------------------------------------------------------------------------------------
627// Unified one-call routine (plan is built once, then just executed each call)
628// --------------------------------------------------------------------------------------
629template < typename GridDataType >
632 const GridDataType& data,
633 SubdomainNeighborhoodSendRecvBuffer< typename GridDataType::value_type,
634 grid::grid_data_vec_dim< GridDataType >() >& recv_buffers,
636{
637 plan.exchange_and_reduce( data, recv_buffers, reduction );
638}
639
640} // namespace terra::communication::shell
Definition communication_plan.hpp:33
ShellBoundaryCommPlan(const grid::shell::DistributedDomain &domain, bool enable_local_comm=true)
Definition communication_plan.hpp:40
void exchange_and_reduce(const GridDataType &data, SubdomainNeighborhoodSendRecvBuffer< ScalarType, VecDim > &boundary_recv_buffers, CommunicationReduction reduction=CommunicationReduction::SUM) const
Definition communication_plan.hpp:48
typename GridDataType::memory_space memory_space
Definition communication_plan.hpp:37
void rebuild()
Definition communication_plan.hpp:71
typename GridDataType::value_type ScalarType
Definition communication_plan.hpp:35
Kokkos::View< ScalarType *, memory_space > rank_buffer_view
Definition communication_plan.hpp:38
static constexpr int VecDim
Definition communication_plan.hpp:36
const grid::Grid0DDataVec< ScalarType, VecDim > & buffer_vertex(const grid::shell::SubdomainInfo &local_subdomain, const grid::BoundaryVertex local_boundary_vertex, const grid::shell::SubdomainInfo &neighbor_subdomain, const grid::BoundaryVertex neighbor_boundary_vertex) const
Const reference to the view that is a buffer for a vertex of a subdomain.
Definition communication.hpp:73
const grid::Grid2DDataVec< ScalarType, VecDim > & buffer_face(const grid::shell::SubdomainInfo &local_subdomain, const grid::BoundaryFace local_boundary_face, const grid::shell::SubdomainInfo &neighbor_subdomain, const grid::BoundaryFace neighbor_boundary_face) const
Const reference to the view that is a buffer for a face of a subdomain.
Definition communication.hpp:110
const grid::Grid1DDataVec< ScalarType, VecDim > & buffer_edge(const grid::shell::SubdomainInfo &local_subdomain, const grid::BoundaryEdge local_boundary_edge, const grid::shell::SubdomainInfo &neighbor_subdomain, const grid::BoundaryEdge neighbor_boundary_edge) const
Const reference to the view that is a buffer for an edge of a subdomain.
Definition communication.hpp:92
Parallel data structure organizing the thick spherical shell metadata for distributed (MPI parallel) ...
Definition spherical_shell.hpp:2498
const DomainInfo & domain_info() const
Returns a const reference.
Definition spherical_shell.hpp:2577
int subdomain_num_nodes_radially() const
Equivalent to calling subdomain_num_nodes_radially( subdomain_refinement_level() )
Definition spherical_shell.hpp:861
(Sortable) Globally unique identifier for a single subdomain of a diamond.
Definition spherical_shell.hpp:595
Timer supporting RAII scope or manual stop.
Definition timer.hpp:270
local_subdomain_id
Definition EpsilonDivDiv_kernel_gen.py:23
Definition communication.hpp:14
void send_recv_with_plan(const ShellBoundaryCommPlan< GridDataType > &plan, const GridDataType &data, SubdomainNeighborhoodSendRecvBuffer< typename GridDataType::value_type, grid::grid_data_vec_dim< GridDataType >() > &recv_buffers, CommunicationReduction reduction=CommunicationReduction::SUM)
Definition communication_plan.hpp:630
constexpr int MPI_TAG_BOUNDARY_DATA
Definition communication.hpp:16
void copy_from_buffer_rotate_and_reduce(const grid::Grid0DDataVec< ScalarType, VecDim > &buffer, const ViewType &data, const int local_subdomain_id, const grid::BoundaryVertex boundary_vertex, const CommunicationReduction reduction)
Definition buffer_copy_kernels.hpp:404
CommunicationReduction
Communication reduction modes.
Definition buffer_copy_kernels.hpp:9
@ SUM
Sums up the node values during receive.
constexpr bool is_edge_boundary_radial(const BoundaryEdge id)
Definition grid_types.hpp:219
BoundaryVertex
Enum for identification of the 8 boundary vertices of a subdomain.
Definition grid_types.hpp:121
constexpr bool is_face_boundary_normal_to_radial_direction(const BoundaryFace id)
Definition grid_types.hpp:225
BoundaryFace
Enum for identification of the 6 boundary faces of a subdomain.
Definition grid_types.hpp:170
BoundaryEdge
Enum for identification of the 12 boundary edges of a subdomain.
Definition grid_types.hpp:142
int MPIRank
Definition mpi.hpp:8
MPIRank rank()
Definition mpi.hpp:10