Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana_Grid_SparseHalo.hpp
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
12#ifndef CABANA_GRID_SPARSEHALO_HPP
13#define CABANA_GRID_SPARSEHALO_HPP
14
16#include <Cabana_SoA.hpp>
17#include <Cabana_Tuple.hpp>
18
19#include <Cabana_Grid_Halo.hpp> // to get the pattern and tags defined here
22#include <Cabana_Grid_Types.hpp>
23
24#include <Kokkos_Core.hpp>
25
26#include <numeric>
27#include <type_traits>
28
29namespace Cabana
30{
31namespace Grid
32{
33namespace Experimental
34{
35//---------------------------------------------------------------------------//
36// Halo
37// ---------------------------------------------------------------------------//
48template <class MemorySpace, class DataTypes, class EntityType,
49 std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim,
50 typename Value = int, typename Key = uint64_t>
52{
53 public:
55 static constexpr std::size_t num_space_dim = NumSpaceDim;
56
58 using memory_space = MemorySpace;
60 using entity_type = EntityType;
61
64
66 using value_type = Value;
68 using key_type = Key;
71 {
72 invalid_key = ~static_cast<key_type>( 0 )
73 };
74
76 using aosoa_member_types = DataTypes;
79
81 template <std::size_t M>
85 static constexpr std::size_t member_num = aosoa_member_types::size;
87 static constexpr unsigned long long cell_bits_per_tile_dim =
88 cellBitsPerTileDim;
90 static constexpr unsigned long long cell_num_per_tile =
91 1 << ( cell_bits_per_tile_dim * 3 );
92
94 using buffer_view = Kokkos::View<tuple_type*, memory_space>;
97 using steering_view = Kokkos::View<key_type*, memory_space>;
102 enum Index
103 {
104 own = 0,
105 ghost = 1,
106 total = 2
107 };
108
111 using counting_view = Kokkos::View<int[2], memory_space,
112 Kokkos::MemoryTraits<Kokkos::Atomic>>;
113
114 //---------------------------------------------------------------------------//
121 template <class SparseArrayType>
123 const std::shared_ptr<SparseArrayType>& sparse_array )
124 : _pattern( pattern )
125 {
126 // Function to get the local id of the neighbor.
127 auto neighbor_id = []( const std::array<int, num_space_dim>& ijk )
128 {
129 int id = ijk[0];
130 for ( std::size_t d = 1; d < num_space_dim; ++d )
131 id += num_space_dim * id + ijk[d];
132 return id;
133 };
134
135 // Neighbor id flip function. This lets us compute what neighbor we
136 // are relative to a given neighbor.
137 auto flip_id = [=]( const std::array<int, num_space_dim>& ijk )
138 {
139 std::array<int, num_space_dim> flip_ijk;
140 for ( std::size_t d = 0; d < num_space_dim; ++d )
141 flip_ijk[d] = -ijk[d];
142 return flip_ijk;
143 };
144
145 // compute total number of bytes in each SoA structure of the sparse
146 // grid (the spasre grid is in AoSoA manner)
147 auto soa_byte_array = compute_member_size_list();
148 for ( std::size_t i = 0; i < member_num; ++i )
149 _soa_member_bytes[i] = soa_byte_array[i];
150 _soa_total_bytes = std::max(
151 std::accumulate( soa_byte_array.begin(), soa_byte_array.end(), 0 ),
152 static_cast<int>( sizeof( tuple_type ) ) );
153
154 // Get the local grid the array uses.
155 auto local_grid = sparse_array->layout().localGrid();
156
157 // linear MPI rank ID of the current working rank
158 _self_rank =
159 local_grid->neighborRank( std::array<int, 3>( { 0, 0, 0 } ) );
160
161 // set the linear neighbor rank ID
162 // set up correspondence between sending and receiving buffers
163 auto neighbors = _pattern.getNeighbors();
164 for ( const auto& n : neighbors )
165 {
166 // neighbor rank linear ID
167 int rank = local_grid->neighborRank( n );
168
169 // if neighbor is valid
170 if ( rank >= 0 )
171 {
172 // add neighbor to the list
173 _neighbor_ranks.push_back( rank );
174 // set the tag to use for sending data to this neighbor
175 // the corresponding neighbor will have the same receiving tag
176 _send_tags.push_back( neighbor_id( n ) );
177 // set the tag to use for receiving data from this neighbor
178 // the corresponding neighbor will have the same sending tag
179 _receive_tags.push_back( neighbor_id( flip_id( n ) ) );
180
181 // build communication data for owned entries
182 buildCommData( Own(), local_grid, n, _owned_buffers,
183 _owned_tile_steering, _owned_tile_spaces );
184 // build communication data for ghosted entries
185 buildCommData( Ghost(), local_grid, n, _ghosted_buffers,
186 _ghosted_tile_steering, _ghosted_tile_spaces );
187
188 auto& own_index_space = _owned_tile_spaces.back();
189 auto& ghost_index_space = _ghosted_tile_spaces.back();
190 int tmp_steering_size =
191 own_index_space.sizeTile() > ghost_index_space.sizeTile()
192 ? own_index_space.sizeTile()
193 : ghost_index_space.sizeTile();
194 _tmp_tile_steering.push_back(
195 steering_view( "tmp_tile_steering", tmp_steering_size ) );
196
197 _valid_counting.push_back(
198 counting_view( "halo_valid_counting" ) );
199 _neighbor_counting.push_back(
200 counting_view( "halo_neighbor_valid_counting" ) );
201 Kokkos::deep_copy( _valid_counting.back(), 0 );
202 Kokkos::deep_copy( _neighbor_counting.back(), 0 );
203
204 _valid_neighbor_ids.emplace_back( n );
205 }
206 }
207 }
208
210 template <class SparseArrayType>
211 MPI_Comm getComm( const SparseArrayType sparse_array ) const
212 {
213 return sparse_array.layout().localGrid()->globalGrid().comm();
214 }
215
227 template <class DecompositionTag, class LocalGridType>
228 void buildCommData( DecompositionTag decomposition_tag,
229 const std::shared_ptr<LocalGridType>& local_grid,
230 const std::array<int, num_space_dim>& nid,
231 std::vector<buffer_view>& buffers,
232 std::vector<steering_view>& steering,
233 std::vector<tile_index_space>& spaces )
234 {
235 // get the halo sparse tile index space sharsed with the neighbor
236 spaces.push_back(
237 local_grid->template sharedTileIndexSpace<cell_bits_per_tile_dim>(
238 decomposition_tag, entity_type(), nid ) );
239 auto& index_space = spaces.back();
240
241 // allocate the buffer to store shared data with given neighbor
242 buffers.push_back(
243 buffer_view( "halo_buffer", index_space.sizeCell() ) );
244 // allocate the steering to guide the communication with the given
245 // neighbor
246 steering.push_back(
247 steering_view( "halo_tile_steering", index_space.sizeTile() ) );
248 // clear steering (init)
249 Kokkos::deep_copy( steering.back(), invalid_key );
250 }
251
252 //---------------------------------------------------------------------------//
258 template <class LocalGridType>
259 void updateTileSpace( const std::shared_ptr<LocalGridType>& local_grid )
260 {
261 // clear index space array first
262 _owned_tile_spaces.clear();
263 _ghosted_tile_spaces.clear();
264
265 // loop over all neighbors and update the shared tile index space
266 auto neighbors = _pattern.getNeighbors();
267 for ( std::size_t i = 0; i < _valid_neighbor_ids.size(); ++i )
268 {
269 // get neighbor relative id
270 auto& n = _valid_neighbor_ids[i];
271 // get neighbor linear MPI rank ID
272 int rank = local_grid->neighborRank( n );
273 // check if neighbor rank is valid
274 // the neighbor id should always be valid (as all should be
275 // well-prepared during construction/initialization)
276 if ( rank == _neighbor_ranks[i] )
277 {
278 // get shared tile index space from local grid
279 _owned_tile_spaces.push_back(
280 local_grid
281 ->template sharedTileIndexSpace<cell_bits_per_tile_dim>(
282 Own(), entity_type(), n ) );
283 _ghosted_tile_spaces.push_back(
284 local_grid
285 ->template sharedTileIndexSpace<cell_bits_per_tile_dim>(
286 Ghost(), entity_type(), n ) );
287
288 // reference to tile index spaces
289 auto& own_index_space = _owned_tile_spaces.back();
290 auto& ghost_index_space = _ghosted_tile_spaces.back();
291
292 // number of tiles inside each shared tile space
293 int own_tile_size = own_index_space.sizeTile();
294 int ghost_tile_size = ghost_index_space.sizeTile();
295 int tmp_steering_size = own_tile_size > ghost_tile_size
296 ? own_tile_size
297 : ghost_tile_size;
298
299 // check if the data steering and buffers require resize
300 // since the partition may changes during simulation process
301 // the size of shared spaces may also change accordingly
302 if ( own_tile_size > _owned_tile_steering[i].extent( 0 ) )
303 {
304 Kokkos::resize( _owned_tile_steering[i], own_tile_size );
305 Kokkos::resize( _owned_buffers[i],
306 own_index_space.sizeCell() );
307 }
308 if ( ghost_tile_size > _ghosted_tile_steering[i].extent( 0 ) )
309 {
310 Kokkos::resize( _ghosted_tile_steering[i],
311 ghost_tile_size );
312 Kokkos::resize( _ghosted_buffers[i],
313 ghost_index_space.sizeCell() );
314 }
315 if ( tmp_steering_size > _tmp_tile_steering[i].extent( 0 ) )
316 Kokkos::resize( _tmp_tile_steering[i], tmp_steering_size );
317 }
318 else
319 std::runtime_error(
320 "Cabana::Grid::Experimental::SparseHalo::updateTileSpace: "
321 "Neighbor rank doesn't match id" );
322 }
323 }
324
325 //---------------------------------------------------------------------------//
334 template <class ExecSpace, class SparseMapType>
335 void register_halo( SparseMapType& map )
336 {
337 // return if there's no valid neighbors
338 int num_n = _neighbor_ranks.size();
339 if ( 0 == num_n )
340 return;
341
342 // loop for all valid neighbors
343 for ( int nid = 0; nid < num_n; nid++ )
344 {
345 auto& owned_space = _owned_tile_spaces[nid];
346 auto& owned_steering = _owned_tile_steering[nid];
347
348 auto& ghosted_space = _ghosted_tile_spaces[nid];
349 auto& ghosted_steering = _ghosted_tile_steering[nid];
350
351 auto& counting = _valid_counting[nid];
352
353 // check the sparse map, if there are valid tiles in corresponding
354 // shared tile index space, add the tile key into the steering view,
355 // this steering keys will later be used for halo data collection
356
357 Kokkos::parallel_for(
358 Kokkos::RangePolicy<ExecSpace>( 0, map.capacity() ),
359 KOKKOS_LAMBDA( const int index ) {
360 if ( map.valid_at( index ) )
361 {
362 auto tile_key = map.key_at( index );
363 int ti, tj, tk;
364 map.key2ijk( tile_key, ti, tj, tk );
365 if ( owned_space.tileInRange( ti, tj, tk ) )
366 {
367 owned_steering( counting( Index::own )++ ) =
368 tile_key;
369 }
370 else if ( ghosted_space.tileInRange( ti, tj, tk ) )
371 {
372 ghosted_steering( counting( Index::ghost )++ ) =
373 tile_key;
374 }
375 }
376 } );
377 Kokkos::fence();
378 }
379 }
380
381 //---------------------------------------------------------------------------//
389 void clear( MPI_Comm comm )
390 {
391 // clear counting
392 for ( std::size_t i = 0; i < _valid_counting.size(); ++i )
393 Kokkos::deep_copy( _valid_counting[i], 0 );
394 for ( std::size_t i = 0; i < _neighbor_counting.size(); ++i )
395 Kokkos::deep_copy( _neighbor_counting[i], 0 );
396 // clear steering
397 for ( std::size_t i = 0; i < _owned_tile_steering.size(); ++i )
398 Kokkos::deep_copy( _owned_tile_steering[i], invalid_key );
399 for ( std::size_t i = 0; i < _ghosted_tile_steering.size(); ++i )
400 Kokkos::deep_copy( _ghosted_tile_steering[i], invalid_key );
401 for ( std::size_t i = 0; i < _tmp_tile_steering.size(); ++i )
402 Kokkos::deep_copy( _tmp_tile_steering[i], invalid_key );
403 // sync
404 MPI_Barrier( comm );
405 }
406
407 //---------------------------------------------------------------------------//
417 MPI_Comm comm, const bool is_neighbor_counting_collected = false ) const
418 {
419 // the valid halo size is already counted, no need to recount
420 if ( is_neighbor_counting_collected )
421 return;
422
423 // number of neighbors
424 int num_n = _neighbor_ranks.size();
425
426 // MPI request to collect counting information in shared owned and
427 // shared ghosted space
428 std::vector<MPI_Request> counting_requests( 2 * num_n,
429 MPI_REQUEST_NULL );
430 const int mpi_tag_counting = 1234;
431
432 // receive from all neighbors
433 for ( int nid = 0; nid < num_n; ++nid )
434 {
435 MPI_Irecv( _neighbor_counting[nid].data(),
436 Index::total * sizeof( int ), MPI_BYTE,
437 _neighbor_ranks[nid],
438 mpi_tag_counting + _receive_tags[nid], comm,
439 &counting_requests[nid] );
440 }
441 // send to all valid neighbors
442 for ( int nid = 0; nid < num_n; ++nid )
443 {
444 MPI_Isend( _valid_counting[nid].data(),
445 Index::total * sizeof( int ), MPI_BYTE,
446 _neighbor_ranks[nid], mpi_tag_counting + _send_tags[nid],
447 comm, &counting_requests[nid + num_n] );
448 }
449
450 // wait until all counting data sending finished
451 const int ec = MPI_Waitall( num_n, counting_requests.data() + num_n,
452 MPI_STATUSES_IGNORE );
453
454 // check if the counting communication succeed
455 if ( MPI_SUCCESS != ec )
456 throw std::logic_error(
457 "Cabana::Grid::Experimental::SparseHalo::"
458 "collectNeighborCounting: counting sending failed." );
459 MPI_Barrier( comm );
460 }
461
472 MPI_Comm comm, std::vector<int>& valid_sends,
473 std::vector<int>& valid_recvs,
474 const bool is_neighbor_counting_collected = false ) const
475 {
476 // collect neighbor counting if needed
477 collectNeighborCounting( comm, is_neighbor_counting_collected );
478
479 // loop over all valid neighbors to check if there's data communication
480 // needed
481 for ( std::size_t nid = 0; nid < _neighbor_ranks.size(); ++nid )
482 {
483 // reference to counting in owned and ghosted share spaces
484 auto h_counting = Kokkos::create_mirror_view_and_copy(
485 Kokkos::HostSpace(), _valid_counting[nid] );
486 auto h_neighbor_counting = Kokkos::create_mirror_view_and_copy(
487 Kokkos::HostSpace(), _neighbor_counting[nid] );
488
489 // if the current ghosted space and the neighbor's owned share
490 // space is non-emty, we need to send data to the neighbor
491 if ( !( h_counting( Index::ghost ) == 0 ||
492 h_neighbor_counting( Index::own ) == 0 ) )
493 {
494 valid_sends.push_back( nid );
495 }
496 // if the current owned share space and the neighbor's ghosted share
497 // space is non-emty, the neighbor will send data to us and we need
498 // to receive data accordingly
499 if ( !( h_counting( Index::own ) == 0 ||
500 h_neighbor_counting( Index::ghost ) == 0 ) )
501 {
502 valid_recvs.push_back( nid );
503 }
504 }
505 }
506
517 MPI_Comm comm, std::vector<int>& valid_sends,
518 std::vector<int>& valid_recvs,
519 const bool is_neighbor_counting_collected = false ) const
520 {
521 // collect neighbor counting if needed
522 collectNeighborCounting( comm, is_neighbor_counting_collected );
523
524 // loop over all valid neighbors to check if there's data communication
525 // needed
526 for ( std::size_t nid = 0; nid < _neighbor_ranks.size(); ++nid )
527 {
528 auto h_counting = Kokkos::create_mirror_view_and_copy(
529 Kokkos::HostSpace(), _valid_counting[nid] );
530 auto h_neighbor_counting = Kokkos::create_mirror_view_and_copy(
531 Kokkos::HostSpace(), _neighbor_counting[nid] );
532
533 // if the current owned space and the neighbor's ghosted share
534 // space is non-emty, we need to send data to the neighbor
535 if ( !( h_counting( Index::own ) == 0 ||
536 h_neighbor_counting( Index::ghost ) == 0 ) )
537 {
538 valid_sends.push_back( nid );
539 }
540
541 // if the current ghosted space and the neighbor's owned share
542 // space is non-emty, we need to send data to the neighbor
543 if ( !( h_counting( Index::ghost ) == 0 ||
544 h_neighbor_counting( Index::own ) == 0 ) )
545 {
546 valid_recvs.push_back( nid );
547 }
548 }
549 }
550
562 template <class ExecSpace, class SparseArrayType>
563 void gather( const ExecSpace& exec_space, SparseArrayType& sparse_array,
564 const bool is_neighbor_counting_collected = false ) const
565 {
566 // return if no valid neighbor
567 if ( 0 == _neighbor_ranks.size() )
568 return;
569
570 // Get the MPI communicator.
571 auto comm = getComm( sparse_array );
572
573 const auto& map = sparse_array.layout().sparseMap();
574
575 // communicate "counting" among neighbors, to decide if the grid data
576 // communication is needed
577 std::vector<int> valid_sends;
578 std::vector<int> valid_recvs;
579 gatherValidSendAndRecvRanks( comm, valid_sends, valid_recvs,
580 is_neighbor_counting_collected );
581 MPI_Barrier( comm );
582
583 // ------------------------------------------------------------------
584 // communicate steering (array keys) for all valid sends and receives
585 std::vector<MPI_Request> steering_requests(
586 valid_recvs.size() + valid_sends.size(), MPI_REQUEST_NULL );
587 const int mpi_tag_steering = 3214;
588
589 // get the steering keys from valid neighbors to get all grids that
590 // we need to receive
591 // loop over all neighbors that will send data to the current rank
592 for ( std::size_t i = 0; i < valid_recvs.size(); ++i )
593 {
594 int nid = valid_recvs[i];
595 Kokkos::View<int[2], Kokkos::HostSpace> h_neighbor_counting(
596 "tmp_host_neighbor_counting" );
597 Kokkos::deep_copy( h_neighbor_counting, _neighbor_counting[nid] );
598
599 MPI_Irecv( _tmp_tile_steering[nid].data(),
600 h_neighbor_counting( Index::own ) * sizeof( key_type ),
601 MPI_BYTE, _neighbor_ranks[nid],
602 mpi_tag_steering + _receive_tags[nid], comm,
603 &steering_requests[i] );
604 }
605
606 // send the steering keys to valid neighbors
607 // loop over all neighbors that requires our owned data
608 for ( std::size_t i = 0; i < valid_sends.size(); ++i )
609 {
610 int nid = valid_sends[i];
611 Kokkos::View<int[2], Kokkos::HostSpace> h_counting(
612 "tmp_host_counting" );
613 Kokkos::deep_copy( h_counting, _valid_counting[nid] );
614
615 MPI_Isend( _owned_tile_steering[nid].data(),
616 h_counting( Index::own ) * sizeof( key_type ), MPI_BYTE,
617 _neighbor_ranks[nid], mpi_tag_steering + _send_tags[nid],
618 comm, &steering_requests[i + valid_recvs.size()] );
619 }
620
621 // wait for all sending work finish
622 const int ec_ss = MPI_Waitall(
623 valid_sends.size(), steering_requests.data() + valid_recvs.size(),
624 MPI_STATUSES_IGNORE );
625 if ( MPI_SUCCESS != ec_ss )
626 throw std::logic_error( "Cabana::Grid::Experimental::SparseHalo::"
627 "gather: steering sending failed." );
628 MPI_Barrier( comm );
629
630 // ------------------------------------------------------------------
631 // communicate sparse array data
632 // Pick a tag to use for communication. This object has its own
633 // communication space so any tag will do.
634 std::vector<MPI_Request> requests(
635 valid_recvs.size() + valid_sends.size(), MPI_REQUEST_NULL );
636 const int mpi_tag = 2345;
637
638 // post receives
639 for ( std::size_t i = 0; i < valid_recvs.size(); ++i )
640 {
641 int nid = valid_recvs[i];
642 Kokkos::View<int[2], Kokkos::HostSpace> h_neighbor_counting(
643 "tmp_host_neighbor_counting" );
644 ;
645 Kokkos::deep_copy( h_neighbor_counting, _neighbor_counting[nid] );
646
647 MPI_Irecv( _ghosted_buffers[nid].data(),
648 h_neighbor_counting( Index::own ) * cell_num_per_tile *
649 _soa_total_bytes,
650 MPI_BYTE, _neighbor_ranks[nid],
651 mpi_tag + _receive_tags[nid], comm, &requests[i] );
652 }
653
654 // pack send buffers and post sends
655 for ( std::size_t i = 0; i < valid_sends.size(); ++i )
656 {
657 int nid = valid_sends[i];
658 Kokkos::View<int[2], Kokkos::HostSpace> h_counting(
659 "tmp_host_counting" );
660 ;
661 Kokkos::deep_copy( h_counting, _valid_counting[nid] );
662
663 packBuffer( exec_space, _owned_buffers[nid],
664 _owned_tile_steering[nid], sparse_array,
665 h_counting( Index::own ) );
666 Kokkos::fence();
667
668 MPI_Isend(
669 _owned_buffers[nid].data(),
670 h_counting( Index::own ) * cell_num_per_tile * _soa_total_bytes,
671 MPI_BYTE, _neighbor_ranks[nid], mpi_tag + _send_tags[nid], comm,
672 &requests[i + valid_recvs.size()] );
673 }
674
675 // unpack receive buffers
676 for ( std::size_t i = 0; i < valid_recvs.size(); ++i )
677 {
678 // get the next buffer to unpack
679 int unpack_index = MPI_UNDEFINED;
680 MPI_Waitany( valid_recvs.size(), requests.data(), &unpack_index,
681 MPI_STATUS_IGNORE );
682
683 // in theory we should receive enough buffers to unpack
684 // if not there could be some problems
685 if ( MPI_UNDEFINED == unpack_index )
686 std::runtime_error(
687 std::string( "Cabana::Grid::Experimental::SparseHalo::"
688 "gather: data receiving failed, "
689 "get only " ) +
690 std::to_string( i ) + ", need " +
691 std::to_string( valid_recvs.size() ) );
692 // otherwise unpack the next buffer
693 else
694 {
695 int nid = valid_recvs[unpack_index];
696 auto h_neighbor_counting = Kokkos::create_mirror_view_and_copy(
697 Kokkos::HostSpace(), _neighbor_counting[nid] );
699 _ghosted_buffers[nid], _tmp_tile_steering[nid],
700 sparse_array, map,
701 h_neighbor_counting( Index::own ) );
702 Kokkos::fence();
703 }
704 }
705
706 // wait to finish all send requests
707 const int ec_data = MPI_Waitall( valid_sends.size(),
708 requests.data() + valid_recvs.size(),
709 MPI_STATUSES_IGNORE );
710 if ( MPI_SUCCESS != ec_data )
711 throw std::logic_error(
712 "sparse_halo_gather: data sending failed." );
713
714 // reinit steerings for next round of communication
715 for ( std::size_t i = 0; i < _tmp_tile_steering.size(); ++i )
716 Kokkos::deep_copy( _tmp_tile_steering[i], invalid_key );
717 MPI_Barrier( comm );
718 }
719
734 template <class ExecSpace, class ReduceOp, class SparseArrayType>
735 void scatter( const ExecSpace& exec_space, const ReduceOp& reduce_op,
736 SparseArrayType& sparse_array,
737 const bool is_neighbor_counting_collected = false ) const
738 {
739 // return if no valid neighbor
740 if ( 0 == _neighbor_ranks.size() )
741 return;
742
743 // Get the MPI communicator.
744 auto comm = getComm( sparse_array );
745
746 const auto& map = sparse_array.layout().sparseMap();
747
748 // communicate "counting" among neighbors, to decide if the grid data
749 // transfer is needed
750 std::vector<int> valid_sends;
751 std::vector<int> valid_recvs;
752 scatterValidSendAndRecvRanks( comm, valid_sends, valid_recvs,
753 is_neighbor_counting_collected );
754 MPI_Barrier( comm );
755
756 // ------------------------------------------------------------------
757 // communicate steering (array keys) for all valid sends and receives
758 std::vector<MPI_Request> steering_requests(
759 valid_recvs.size() + valid_sends.size(), MPI_REQUEST_NULL );
760 const int mpi_tag_steering = 214;
761
762 // get the steering keys from valid neighbors to know all grids that
763 // we need to receive
764 // loop over all neighbors that will send data to the current rank
765 for ( std::size_t i = 0; i < valid_recvs.size(); ++i )
766 {
767 int nid = valid_recvs[i];
768 Kokkos::View<int[2], Kokkos::HostSpace> h_neighbor_counting(
769 "tmp_host_neighbor_counting" );
770 ;
771 Kokkos::deep_copy( h_neighbor_counting, _neighbor_counting[nid] );
772
773 MPI_Irecv( _tmp_tile_steering[nid].data(),
774 h_neighbor_counting( Index::ghost ) * sizeof( key_type ),
775 MPI_BYTE, _neighbor_ranks[nid],
776 mpi_tag_steering + _receive_tags[nid], comm,
777 &steering_requests[i] );
778 }
779
780 // send the steering keys to valid neighbors
781 // loop over all neighbors that requires our owned data
782 for ( std::size_t i = 0; i < valid_sends.size(); ++i )
783 {
784 int nid = valid_sends[i];
785 Kokkos::View<int[2], Kokkos::HostSpace> h_counting(
786 "tmp_host_counting" );
787 ;
788 Kokkos::deep_copy( h_counting, _valid_counting[nid] );
789
790 MPI_Isend( _ghosted_tile_steering[nid].data(),
791 h_counting( Index::ghost ) * sizeof( key_type ),
792 MPI_BYTE, _neighbor_ranks[nid],
793 mpi_tag_steering + _send_tags[nid], comm,
794 &steering_requests[i + valid_recvs.size()] );
795 }
796
797 // wait for all sending work finish
798 const int ec_ss = MPI_Waitall(
799 valid_sends.size(), steering_requests.data() + valid_recvs.size(),
800 MPI_STATUSES_IGNORE );
801 if ( MPI_SUCCESS != ec_ss )
802 throw std::logic_error( "Cabana::Grid::Experimental::SparseHalo::"
803 "scatter: steering sending failed." );
804 MPI_Barrier( comm );
805
806 // ------------------------------------------------------------------
807 // communicate sparse array data
808 // Pick a tag to use for communication. This object has its own
809 // communication space so any tag will do.
810 std::vector<MPI_Request> requests(
811 valid_recvs.size() + valid_sends.size(), MPI_REQUEST_NULL );
812 const int mpi_tag = 345;
813
814 // post receives
815 for ( std::size_t i = 0; i < valid_recvs.size(); ++i )
816 {
817 int nid = valid_recvs[i];
818 Kokkos::View<int[2], Kokkos::HostSpace> h_neighbor_counting(
819 "tmp_host_neighbor_counting" );
820 ;
821 Kokkos::deep_copy( h_neighbor_counting, _neighbor_counting[nid] );
822
823 MPI_Irecv( _owned_buffers[nid].data(),
824 h_neighbor_counting( Index::ghost ) * cell_num_per_tile *
825 _soa_total_bytes,
826 MPI_BYTE, _neighbor_ranks[nid],
827 mpi_tag + _receive_tags[nid], comm, &requests[i] );
828 }
829
830 // pack send buffers and post sends
831 for ( std::size_t i = 0; i < valid_sends.size(); ++i )
832 {
833 int nid = valid_sends[i];
834 Kokkos::View<int[2], Kokkos::HostSpace> h_counting(
835 "tmp_host_counting" );
836 ;
837 Kokkos::deep_copy( h_counting, _valid_counting[nid] );
838 packBuffer( exec_space, _ghosted_buffers[nid],
839 _ghosted_tile_steering[nid], sparse_array,
840 h_counting( Index::ghost ) );
841 Kokkos::fence();
842
843 MPI_Isend( _ghosted_buffers[nid].data(),
844 h_counting( Index::ghost ) * cell_num_per_tile *
845 _soa_total_bytes,
846 MPI_BYTE, _neighbor_ranks[nid],
847 mpi_tag + _send_tags[nid], comm,
848 &requests[i + valid_recvs.size()] );
849 }
850
851 // unpack receive buffers
852 for ( std::size_t i = 0; i < valid_recvs.size(); ++i )
853 {
854 // get the next buffer to unpack
855 int unpack_index = MPI_UNDEFINED;
856 MPI_Waitany( valid_recvs.size(), requests.data(), &unpack_index,
857 MPI_STATUS_IGNORE );
858
859 // in theory we should receive enough buffers to unpack
860 // if not there could be some problems
861 if ( MPI_UNDEFINED == unpack_index )
862 std::runtime_error(
863 std::string( "sparse_halo_scatter: data receiving failed, "
864 "get only " ) +
865 std::to_string( i ) + ", need " +
866 std::to_string( valid_recvs.size() ) );
867 // otherwise unpack the next buffer with the given reduce operator
868 else
869 {
870 int nid = valid_recvs[unpack_index];
871 auto h_neighbor_counting = Kokkos::create_mirror_view_and_copy(
872 Kokkos::HostSpace(), _neighbor_counting[nid] );
873 unpackBuffer( reduce_op, exec_space, _owned_buffers[nid],
874 _tmp_tile_steering[nid], sparse_array, map,
875 h_neighbor_counting( Index::ghost ) );
876 Kokkos::fence();
877 }
878 }
879
880 // wait to finish all send requests
881 const int ec_data = MPI_Waitall( valid_sends.size(),
882 requests.data() + valid_recvs.size(),
883 MPI_STATUSES_IGNORE );
884 if ( MPI_SUCCESS != ec_data )
885 throw std::logic_error( "Cabana::Grid::Experimental::SparseHalo::"
886 "scatter: data sending failed." );
887
888 // reinit steerings for next round of communication
889 for ( std::size_t i = 0; i < _tmp_tile_steering.size(); ++i )
890 Kokkos::deep_copy( _tmp_tile_steering[i], invalid_key );
891 MPI_Barrier( comm );
892 }
893
894 //---------------------------------------------------------------------------//
905 template <class ExecSpace, class SparseArrayType>
906 void packBuffer( const ExecSpace& exec_space, const buffer_view& buffer,
907 const steering_view& tile_steering,
908 SparseArrayType& sparse_array, const int count ) const
909 {
910 Kokkos::parallel_for(
911 "Cabana::Grid::Experimental::SparseHalo::packBuffer",
912 Kokkos::RangePolicy<ExecSpace>( exec_space, 0, count ),
913 KOKKOS_LAMBDA( const int i ) {
914 if ( tile_steering( i ) != invalid_key )
915 {
916 const int buffer_idx = i * cell_num_per_tile;
917 auto tile_key = tile_steering( i );
918 for ( int lcid = 0; lcid < (int)cell_num_per_tile; ++lcid )
919 {
920 buffer( buffer_idx + lcid ) =
921 sparse_array.getTuple( tile_key, lcid );
922 }
923 }
924 } );
925 }
926
927 //---------------------------------------------------------------------------//
929 template <class T>
930 KOKKOS_INLINE_FUNCTION static void
931 unpackOp( ScatterReduce::Sum, const T& buffer_val, T& array_val )
932 {
933 array_val += buffer_val;
934 }
935
937 template <class T>
938 KOKKOS_INLINE_FUNCTION static void
939 unpackOp( ScatterReduce::Min, const T& buffer_val, T& array_val )
940 {
941 if ( buffer_val < array_val )
942 array_val = buffer_val;
943 }
944
946 template <class T>
947 KOKKOS_INLINE_FUNCTION static void
948 unpackOp( ScatterReduce::Max, const T& buffer_val, T& array_val )
949 {
950 if ( buffer_val > array_val )
951 array_val = buffer_val;
952 }
953
955 template <class T>
956 KOKKOS_INLINE_FUNCTION static void
957 unpackOp( ScatterReduce::Replace, const T& buffer_val, T& array_val )
958 {
959 array_val = buffer_val;
960 }
961
962 //---------------------------------------------------------------------------//
976 template <class ReduceOp, std::size_t N, std::size_t M, class SoAType>
977 KOKKOS_FORCEINLINE_FUNCTION static std::enable_if_t<3 == M, void>
978 unpackTupleMember( const ReduceOp& reduce_op, const tuple_type& src_tuple,
979 SoAType& dst_soa, const int soa_idx,
980 const Kokkos::Array<std::size_t, M>& extents,
981 const std::integral_constant<std::size_t, N>,
982 const std::integral_constant<std::size_t, M> )
983 {
984 for ( std::size_t d0 = 0; d0 < extents[0]; ++d0 )
985 for ( int d1 = 0; d1 < extents[1]; ++d1 )
986 for ( int d2 = 0; d2 < extents[2]; ++d2 )
987 {
988 unpackOp( reduce_op,
989 Cabana::get<N>( src_tuple, d0, d1, d2 ),
990 Cabana::get<N>( dst_soa, soa_idx, d0, d1, d2 ) );
991 }
992 }
993
1007 template <class ReduceOp, std::size_t N, std::size_t M, class SoAType>
1008 KOKKOS_FORCEINLINE_FUNCTION static std::enable_if_t<2 == M, void>
1009 unpackTupleMember( const ReduceOp& reduce_op, const tuple_type& src_tuple,
1010 SoAType& dst_soa, const int soa_idx,
1011 const Kokkos::Array<std::size_t, M>& extents,
1012 const std::integral_constant<std::size_t, N>,
1013 const std::integral_constant<std::size_t, M> )
1014 {
1015 for ( std::size_t d0 = 0; d0 < extents[0]; ++d0 )
1016 for ( int d1 = 0; d1 < extents[1]; ++d1 )
1017 {
1018 unpackOp( reduce_op, Cabana::get<N>( src_tuple, d0, d1 ),
1019 Cabana::get<N>( dst_soa, soa_idx, d0, d1 ) );
1020 }
1021 }
1022
1036 template <class ReduceOp, std::size_t N, std::size_t M, class SoAType>
1037 KOKKOS_FORCEINLINE_FUNCTION static std::enable_if_t<1 == M, void>
1038 unpackTupleMember( const ReduceOp& reduce_op, const tuple_type& src_tuple,
1039 SoAType& dst_soa, const int soa_idx,
1040 const Kokkos::Array<std::size_t, M>& extents,
1041 const std::integral_constant<std::size_t, N>,
1042 const std::integral_constant<std::size_t, M> )
1043 {
1044 for ( std::size_t d0 = 0; d0 < extents[0]; ++d0 )
1045 {
1046 unpackOp( reduce_op, Cabana::get<N>( src_tuple, d0 ),
1047 Cabana::get<N>( dst_soa, soa_idx, d0 ) );
1048 }
1049 }
1050
1063 template <class ReduceOp, std::size_t N, std::size_t M, class SoAType>
1064 KOKKOS_FORCEINLINE_FUNCTION static std::enable_if_t<0 == M, void>
1065 unpackTupleMember( const ReduceOp& reduce_op, const tuple_type& src_tuple,
1066 SoAType& dst_soa, const int soa_idx,
1067 const Kokkos::Array<std::size_t, M>&,
1068 const std::integral_constant<std::size_t, N>,
1069 const std::integral_constant<std::size_t, M> )
1070 {
1071 unpackOp( reduce_op, Cabana::get<N>( src_tuple ),
1072 Cabana::get<N>( dst_soa, soa_idx ) );
1073 }
1074
1084 template <class ReduceOp, class SoAType>
1085 KOKKOS_FORCEINLINE_FUNCTION static void
1086 unpackTuple( const ReduceOp& reduce_op, const tuple_type& src_tuple,
1087 SoAType& dst_soa, const int soa_idx,
1088 const std::integral_constant<std::size_t, 0> )
1089 {
1090 using current_type = member_data_type<0>;
1091 auto extents = compute_member_extents<current_type>();
1093 reduce_op, src_tuple, dst_soa, soa_idx, extents,
1094 std::integral_constant<std::size_t, 0>(),
1095 std::integral_constant<std::size_t,
1096 std::rank<current_type>::value>() );
1097 }
1098
1109 template <class ReduceOp, std::size_t N, class SoAType>
1110 KOKKOS_FORCEINLINE_FUNCTION static void
1111 unpackTuple( const ReduceOp& reduce_op, const tuple_type& src_tuple,
1112 SoAType& dst_soa, const int soa_idx,
1113 const std::integral_constant<std::size_t, N> )
1114 {
1115 using current_type = member_data_type<N>;
1116 auto extents = compute_member_extents<current_type>();
1118 reduce_op, src_tuple, dst_soa, soa_idx, extents,
1119 std::integral_constant<std::size_t, N>(),
1120 std::integral_constant<std::size_t,
1121 std::rank<current_type>::value>() );
1122
1123 if ( N > 1 )
1124 {
1125 // recurcively unpack the next tuple element
1126 unpackTuple( reduce_op, src_tuple, dst_soa, soa_idx,
1127 std::integral_constant<std::size_t, N - 1>() );
1128 }
1129 else
1130 {
1131 unpackTuple( reduce_op, src_tuple, dst_soa, soa_idx,
1132 std::integral_constant<std::size_t, 0>() );
1133 }
1134 }
1135
1150 template <class ReduceOp, class ExecSpace, class SparseArrayType,
1151 class SparseMapType>
1152 void unpackBuffer( const ReduceOp& reduce_op, const ExecSpace& exec_space,
1153 const buffer_view& buffer,
1154 const steering_view& tile_steering,
1155 const SparseArrayType& sparse_array, SparseMapType& map,
1156 const int count ) const
1157 {
1158 Kokkos::parallel_for(
1159 "Cabana::Grid::Experimental::SparseHalo::unpackBuffer",
1160 Kokkos::RangePolicy<ExecSpace>( exec_space, 0, count ),
1161 KOKKOS_LAMBDA( const int i ) {
1162 if ( tile_steering( i ) != invalid_key )
1163 {
1164 auto tile_key = tile_steering( i );
1165 if ( map.isValidKey( tile_key ) )
1166 {
1167 int ti, tj, tk;
1168 map.key2ijk( tile_key, ti, tj, tk );
1169
1170 auto tile_id = map.queryTileFromTileKey( tile_key );
1171 const int buffer_idx = i * cell_num_per_tile;
1172 for ( int lcid = 0; lcid < (int)cell_num_per_tile;
1173 ++lcid )
1174 {
1175 auto& tuple = buffer( buffer_idx + lcid );
1176 auto& data_access =
1177 sparse_array.accessTile( tile_id );
1179 reduce_op, tuple, data_access, lcid,
1180 std::integral_constant<std::size_t,
1181 member_num - 1>() );
1182 }
1183 }
1184 }
1185 } );
1186 }
1187
1188 private:
1189 // These functions may be useful if added to Cabana_MemberType.hpp
1190 // compute member size
1191 template <std::size_t M>
1192 static constexpr std::size_t compute_member_size()
1193 {
1194 return sizeof( member_data_type<M> );
1195 }
1196
1197 template <typename Sequence>
1198 struct compute_member_size_list_impl;
1199
1200 template <std::size_t... Is>
1201 struct compute_member_size_list_impl<std::index_sequence<Is...>>
1202 {
1203 std::array<std::size_t, member_num> operator()()
1204 {
1205 return { compute_member_size<Is>()... };
1206 }
1207 };
1208
1209 template <std::size_t N = member_num,
1210 typename Indices = std::make_index_sequence<N>>
1211 std::array<std::size_t, member_num> compute_member_size_list()
1212 {
1213 compute_member_size_list_impl<Indices> op;
1214 return op();
1215 }
1216
1217 // member extent
1218 template <typename Type, std::size_t M>
1219 KOKKOS_FORCEINLINE_FUNCTION static constexpr std::size_t
1220 compute_one_member_extent()
1221 {
1222 return std::extent<Type, M>::value;
1223 }
1224
1225 template <class Type, std::size_t M, typename Sequence>
1226 struct compute_member_extents_impl;
1227
1228 template <class Type, std::size_t M, std::size_t... Is>
1229 struct compute_member_extents_impl<Type, M, std::index_sequence<Is...>>
1230 {
1231 KOKKOS_FORCEINLINE_FUNCTION
1232 Kokkos::Array<std::size_t, M> operator()()
1233 {
1234 return { compute_one_member_extent<Type, Is>()... };
1235 }
1236 };
1237
1238 template <class Type, std::size_t M = std::rank<Type>::value,
1239 typename Indices = std::make_index_sequence<M>>
1240 KOKKOS_FORCEINLINE_FUNCTION static Kokkos::Array<std::size_t, M>
1241 compute_member_extents()
1242 {
1243 compute_member_extents_impl<Type, M, Indices> op;
1244 return op();
1245 }
1246
1247 private:
1248 // Current MPI linear rank ID
1249 int _self_rank;
1250 // Hallo pattern
1251 halo_pattern_type _pattern;
1252
1253 // neighbor rank linear MPI rank IDs
1254 std::vector<int> _neighbor_ranks;
1255 // valid neigber rank indices; valid means require data communication
1256 std::vector<std::array<int, num_space_dim>> _valid_neighbor_ids;
1257 // sending tags
1258 std::vector<int> _send_tags;
1259 // receiving tags
1260 std::vector<int> _receive_tags;
1261
1262 // owned view buffers
1263 std::vector<buffer_view> _owned_buffers;
1264 // ghosted view buffers
1265 std::vector<buffer_view> _ghosted_buffers;
1266
1267 // owned tile key steerings
1268 std::vector<steering_view> _owned_tile_steering;
1269 // key steering buffers (used to store valid keys get from neighbors)
1270 std::vector<steering_view> _tmp_tile_steering;
1271 // ghosted tile key steerings
1272 std::vector<steering_view> _ghosted_tile_steering;
1273
1274 // valid halo grid counting on current rank (each element map to a neighbor)
1275 std::vector<counting_view> _valid_counting;
1276 // valid halo grid counting on corresponding neighbor ranks
1277 std::vector<counting_view> _neighbor_counting;
1278
1279 // owned tile space
1280 std::vector<tile_index_space> _owned_tile_spaces;
1281 // ghosted tile space
1282 std::vector<tile_index_space> _ghosted_tile_spaces;
1283
1284 // SoA member bytes num
1285 Kokkos::Array<std::size_t, member_num> _soa_member_bytes;
1286 // SoA total bytes count
1287 std::size_t _soa_total_bytes;
1288};
1289
1290//---------------------------------------------------------------------------//
1291// Sparse halo creation.
1292//---------------------------------------------------------------------------//
1298template <class MemorySpace, unsigned long long cellBitsPerTileDim,
1299 class DataTypes, class EntityType, class MeshType,
1300 class SparseMapType, class Pattern, typename Value = int,
1301 typename Key = uint64_t>
1302auto createSparseHalo(
1303 const Pattern& pattern,
1304 const std::shared_ptr<SparseArray<DataTypes, MemorySpace, EntityType,
1305 MeshType, SparseMapType>>
1306 array )
1307{
1308 using array_type = SparseArray<DataTypes, MemorySpace, EntityType, MeshType,
1309 SparseMapType>;
1310 using memory_space = typename array_type::memory_space;
1311 static constexpr std::size_t num_space_dim = array_type::num_space_dim;
1312 return std::make_shared<
1313 SparseHalo<memory_space, DataTypes, EntityType, num_space_dim,
1314 cellBitsPerTileDim, Value, Key>>( pattern, array );
1315}
1316
1317} // namespace Experimental
1318} // namespace Grid
1319} // namespace Cabana
1320
1321#endif // CABANA_GRID_SPARSEHALO_HPP
Multi-node grid scatter/gather.
Sparse grid fields arrays using AoSoA.
Grid type tags.
AoSoA tuple member types.
Struct-of-Arrays for building AoSoA.
Tuple of single particle information to build AoSoA.
Sparse array of field data on the local sparse mesh; Array data is stored in AoSoA manner,...
Definition Cabana_Grid_SparseArray.hpp:316
Definition Cabana_Grid_SparseHalo.hpp:52
static KOKKOS_INLINE_FUNCTION void unpackOp(ScatterReduce::Min, const T &buffer_val, T &array_val)
Reduce an element into the buffer. Min reduction.
Definition Cabana_Grid_SparseHalo.hpp:939
void gather(const ExecSpace &exec_space, SparseArrayType &sparse_array, const bool is_neighbor_counting_collected=false) const
Gather data into our ghosted share space from their owners.
Definition Cabana_Grid_SparseHalo.hpp:563
DataTypes aosoa_member_types
data members in AoSoA structure
Definition Cabana_Grid_SparseHalo.hpp:76
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 2==M, void > unpackTupleMember(const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const Kokkos::Array< std::size_t, M > &extents, const std::integral_constant< std::size_t, N >, const std::integral_constant< std::size_t, M >)
Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 2)
Definition Cabana_Grid_SparseHalo.hpp:1009
Kokkos::View< key_type *, memory_space > steering_view
Definition Cabana_Grid_SparseHalo.hpp:97
void register_halo(SparseMapType &map)
register valid halos (according to grid activation status in sparse map) in the steerings
Definition Cabana_Grid_SparseHalo.hpp:335
Index
index (own or ghost)
Definition Cabana_Grid_SparseHalo.hpp:103
typename Cabana::MemberTypeAtIndex< M, aosoa_member_types >::type member_data_type
AoSoA member data type.
Definition Cabana_Grid_SparseHalo.hpp:82
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 3==M, void > unpackTupleMember(const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const Kokkos::Array< std::size_t, M > &extents, const std::integral_constant< std::size_t, N >, const std::integral_constant< std::size_t, M >)
Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 3)
Definition Cabana_Grid_SparseHalo.hpp:978
static KOKKOS_INLINE_FUNCTION void unpackOp(ScatterReduce::Replace, const T &buffer_val, T &array_val)
Reduce an element into the buffer. Replace reduction.
Definition Cabana_Grid_SparseHalo.hpp:957
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 0==M, void > unpackTupleMember(const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const Kokkos::Array< std::size_t, M > &, const std::integral_constant< std::size_t, N >, const std::integral_constant< std::size_t, M >)
Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 0)
Definition Cabana_Grid_SparseHalo.hpp:1065
void updateTileSpace(const std::shared_ptr< LocalGridType > &local_grid)
update tile index space according to current partition
Definition Cabana_Grid_SparseHalo.hpp:259
void scatter(const ExecSpace &exec_space, const ReduceOp &reduce_op, SparseArrayType &sparse_array, const bool is_neighbor_counting_collected=false) const
Scatter data from our ghosts to their owners using the given type of reduce operation.
Definition Cabana_Grid_SparseHalo.hpp:735
void gatherValidSendAndRecvRanks(MPI_Comm comm, std::vector< int > &valid_sends, std::vector< int > &valid_recvs, const bool is_neighbor_counting_collected=false) const
collect all valid ranks for sparse grid gather operations
Definition Cabana_Grid_SparseHalo.hpp:516
static constexpr std::size_t member_num
AoSoA member #.
Definition Cabana_Grid_SparseHalo.hpp:85
EntityType entity_type
entity type on sparse grid
Definition Cabana_Grid_SparseHalo.hpp:60
SparseHalo(const halo_pattern_type pattern, const std::shared_ptr< SparseArrayType > &sparse_array)
constructor
Definition Cabana_Grid_SparseHalo.hpp:122
KeyValue
invalid key in sparse map
Definition Cabana_Grid_SparseHalo.hpp:71
Key key_type
key type in sparse map
Definition Cabana_Grid_SparseHalo.hpp:68
static KOKKOS_FORCEINLINE_FUNCTION void unpackTuple(const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const std::integral_constant< std::size_t, 0 >)
Unpack a sparse arrays tuple for it's member with index 0.
Definition Cabana_Grid_SparseHalo.hpp:1086
NodeHaloPattern< NumSpaceDim > halo_pattern_type
sparse grid halo pattern (TODO currently reusing Halo's Node pattern)
Definition Cabana_Grid_SparseHalo.hpp:63
static KOKKOS_INLINE_FUNCTION void unpackOp(ScatterReduce::Sum, const T &buffer_val, T &array_val)
Reduce an element into the buffer. Sum reduction.
Definition Cabana_Grid_SparseHalo.hpp:931
MPI_Comm getComm(const SparseArrayType sparse_array) const
Get the communicator.
Definition Cabana_Grid_SparseHalo.hpp:211
void clear(MPI_Comm comm)
clear guiding information in sparse halo,
Definition Cabana_Grid_SparseHalo.hpp:389
void packBuffer(const ExecSpace &exec_space, const buffer_view &buffer, const steering_view &tile_steering, SparseArrayType &sparse_array, const int count) const
Pack sparse arrays at halo regions into a buffer.
Definition Cabana_Grid_SparseHalo.hpp:906
Kokkos::View< int[2], memory_space, Kokkos::MemoryTraits< Kokkos::Atomic > > counting_view
Definition Cabana_Grid_SparseHalo.hpp:111
Value value_type
value type of entities on sparse grid
Definition Cabana_Grid_SparseHalo.hpp:66
void collectNeighborCounting(MPI_Comm comm, const bool is_neighbor_counting_collected=false) const
neighbor tile counting, communication needed only if the counting is non-zero
Definition Cabana_Grid_SparseHalo.hpp:416
void scatterValidSendAndRecvRanks(MPI_Comm comm, std::vector< int > &valid_sends, std::vector< int > &valid_recvs, const bool is_neighbor_counting_collected=false) const
collect all valid ranks for sparse grid scatter operations
Definition Cabana_Grid_SparseHalo.hpp:471
static constexpr std::size_t num_space_dim
sparse array dimension number
Definition Cabana_Grid_SparseHalo.hpp:55
void unpackBuffer(const ReduceOp &reduce_op, const ExecSpace &exec_space, const buffer_view &buffer, const steering_view &tile_steering, const SparseArrayType &sparse_array, SparseMapType &map, const int count) const
Unpack a sparse array communication buffer.
Definition Cabana_Grid_SparseHalo.hpp:1152
void buildCommData(DecompositionTag decomposition_tag, const std::shared_ptr< LocalGridType > &local_grid, const std::array< int, num_space_dim > &nid, std::vector< buffer_view > &buffers, std::vector< steering_view > &steering, std::vector< tile_index_space > &spaces)
Build communication data.
Definition Cabana_Grid_SparseHalo.hpp:228
Cabana::Tuple< aosoa_member_types > tuple_type
AoSoA tuple type.
Definition Cabana_Grid_SparseHalo.hpp:78
Kokkos::View< tuple_type *, memory_space > buffer_view
communication data buffer view type
Definition Cabana_Grid_SparseHalo.hpp:94
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 1==M, void > unpackTupleMember(const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const Kokkos::Array< std::size_t, M > &extents, const std::integral_constant< std::size_t, N >, const std::integral_constant< std::size_t, M >)
Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 1)
Definition Cabana_Grid_SparseHalo.hpp:1038
static KOKKOS_FORCEINLINE_FUNCTION void unpackTuple(const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const std::integral_constant< std::size_t, N >)
Unpack a sparse arrays tuple for all members when element ID!=0.
Definition Cabana_Grid_SparseHalo.hpp:1111
TileIndexSpace< num_space_dim, cell_bits_per_tile_dim > tile_index_space
tile index space type TODO
Definition Cabana_Grid_SparseHalo.hpp:99
MemorySpace memory_space
memory space to store the sparse grid
Definition Cabana_Grid_SparseHalo.hpp:58
static KOKKOS_INLINE_FUNCTION void unpackOp(ScatterReduce::Max, const T &buffer_val, T &array_val)
Reduce an element into the buffer. Max reduction.
Definition Cabana_Grid_SparseHalo.hpp:948
static constexpr unsigned long long cell_num_per_tile
sparse grid hierarchy: cell # per dimension
Definition Cabana_Grid_SparseHalo.hpp:90
static constexpr unsigned long long cell_bits_per_tile_dim
sparse grid hierarchy: cell id bit# per dimension
Definition Cabana_Grid_SparseHalo.hpp:87
Definition Cabana_Grid_Halo.hpp:76
Index space with tile as unit; _min and _max forms the tile range. Note this is for sparse grid only,...
Definition Cabana_Grid_SparseIndexSpace.hpp:1137
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
KOKKOS_FORCEINLINE_FUNCTION std::enable_if< is_parameter_pack< ParameterPack_t >::value, typenameParameterPack_t::templatevalue_type< N > & >::type get(ParameterPack_t &pp)
Get an element from a parameter pack.
Definition Cabana_ParameterPack.hpp:129
Ghosted decomposition tag.
Definition Cabana_Grid_Types.hpp:197
Owned decomposition tag.
Definition Cabana_Grid_Types.hpp:190
Definition Cabana_Grid_Halo.hpp:173
Definition Cabana_Grid_Halo.hpp:167
Definition Cabana_Grid_Halo.hpp:181
Sum values from neighboring ranks into this rank's data.
Definition Cabana_Grid_Halo.hpp:161
Get the type of the member at a given index.
Definition Cabana_MemberTypes.hpp:75
Definition Cabana_Tuple.hpp:32