16#ifndef CABANA_GRID_HALO_HPP
17#define CABANA_GRID_HALO_HPP
23#include <Cabana_Core_Config.hpp>
26#include <Kokkos_Core.hpp>
27#include <Kokkos_Profiling_ScopedRegion.hpp>
46template <std::
size_t NumSpaceDim>
61 setNeighbors(
const std::vector<std::array<int, num_space_dim>>& neighbors )
63 _neighbors = neighbors;
73 std::vector<std::array<int, num_space_dim>> _neighbors;
78template <std::
size_t NumSpaceDim>
84class NodeHaloPattern<3> :
public HaloPattern<3>
90 std::vector<std::array<int, 3>> neighbors;
91 neighbors.reserve( 26 );
92 for (
int i = -1; i < 2; ++i )
93 for (
int j = -1; j < 2; ++j )
94 for (
int k = -1; k < 2; ++k )
95 if ( !( i == 0 && j == 0 && k == 0 ) )
96 neighbors.push_back( { i, j, k } );
104class NodeHaloPattern<2> :
public HaloPattern<2>
110 std::vector<std::array<int, 2>> neighbors;
111 neighbors.reserve( 8 );
112 for (
int i = -1; i < 2; ++i )
113 for (
int j = -1; j < 2; ++j )
114 if ( !( i == 0 && j == 0 ) )
115 neighbors.push_back( { i, j } );
122template <std::
size_t NumSpaceDim>
128class FaceHaloPattern<3> :
public HaloPattern<3>
134 std::vector<std::array<int, 3>> neighbors = {
135 { -1, 0, 0 }, { 1, 0, 0 }, { 0, -1, 0 },
136 { 0, 1, 0 }, { 0, 0, -1 }, { 0, 0, 1 } };
144class FaceHaloPattern<2> :
public HaloPattern<2>
150 std::vector<std::array<int, 2>> neighbors = {
151 { -1, 0 }, { 1, 0 }, { 0, -1 }, { 0, 1 } };
159namespace ScatterReduce
202template <
class MemorySpace>
219 template <
class Pattern,
class... ArrayTypes>
221 const ArrayTypes&... arrays )
224 const std::size_t num_space_dim = Pattern::num_space_dim;
230 auto neighbor_id = [](
const std::array<int, num_space_dim>& ijk )
233 for ( std::size_t d = 1; d < num_space_dim; ++d )
234 id += num_space_dim *
id + ijk[d];
240 auto flip_id = [=](
const std::array<int, num_space_dim>& ijk )
242 std::array<int, num_space_dim> flip_ijk;
243 for ( std::size_t d = 0; d < num_space_dim; ++d )
244 flip_ijk[d] = -ijk[d];
251 auto neighbors = pattern.getNeighbors();
252 for (
const auto& n : neighbors )
255 int rank = local_grid->neighborRank( n );
284 template <
class Array_t>
285 MPI_Comm
getComm(
const Array_t& array )
const
287 return array.layout()->localGrid()->globalGrid().comm();
291 template <
class Array_t,
class... ArrayTypes>
292 MPI_Comm
getComm(
const Array_t& array,
const ArrayTypes&... arrays )
const
299 MPI_Comm_compare( comm,
getComm( arrays... ), &result );
300 if ( result != MPI_IDENT && result != MPI_CONGRUENT )
301 throw std::runtime_error(
"Cabana::Grid::Halo::getComm: Arrays "
302 "have different communicators" );
309 template <
class Array_t>
312 return array.layout()->localGrid();
317 template <
class Array_t,
class... ArrayTypes>
324 if ( local_grid->haloCellWidth() !=
325 array.layout()->localGrid()->haloCellWidth() )
327 throw std::runtime_error(
"Cabana::Grid::Halo::getlocalGrid: "
328 "Arrays have different halo widths" );
335 template <
class DecompositionTag, std::size_t NumSpaceDim,
339 const std::array<int, NumSpaceDim>& nid,
340 std::vector<Kokkos::View<char*, memory_space>>& buffers,
341 std::vector<Kokkos::View<int**, memory_space>>& steering,
342 const ArrayTypes&... arrays )
345 const std::size_t num_array =
sizeof...( ArrayTypes );
348 std::array<std::size_t, num_array> value_byte_sizes = {
349 sizeof(
typename ArrayTypes::value_type )... };
353 std::array<IndexSpace<NumSpaceDim + 1>, num_array> spaces = {
354 ( arrays.layout()->sharedIndexSpace( decomposition_tag, nid,
359 int buffer_bytes = 0;
360 int buffer_num_element = 0;
361 for ( std::size_t a = 0; a < num_array; ++a )
363 buffer_bytes += value_byte_sizes[a] * spaces[a].size();
364 buffer_num_element += spaces[a].size();
370 Kokkos::View<char*, memory_space>(
"halo_buffer", buffer_bytes ) );
373 steering.push_back( Kokkos::View<int**, memory_space>(
374 "steering", buffer_num_element, 3 + NumSpaceDim ) );
378 buffer_num_element, steering );
382 template <std::
size_t NumArray>
385 const std::array<std::size_t, NumArray>& value_byte_sizes,
386 const int buffer_bytes,
const int buffer_num_element,
387 std::vector<Kokkos::View<int**, memory_space>>& steering )
393 Kokkos::create_mirror_view( Kokkos::HostSpace(), steering.back() );
394 int elem_counter = 0;
395 int byte_counter = 0;
396 for ( std::size_t a = 0; a < NumArray; ++a )
398 for (
int i = spaces[a].min( 0 ); i < spaces[a].max( 0 ); ++i )
400 for (
int j = spaces[a].min( 1 ); j < spaces[a].max( 1 ); ++j )
402 for (
int k = spaces[a].min( 2 ); k < spaces[a].max( 2 );
405 for (
int l = spaces[a].min( 3 );
406 l < spaces[a].max( 3 ); ++l )
409 host_steering( elem_counter, 0 ) = byte_counter;
412 host_steering( elem_counter, 1 ) = a;
415 host_steering( elem_counter, 2 ) = i;
416 host_steering( elem_counter, 3 ) = j;
417 host_steering( elem_counter, 4 ) = k;
418 host_steering( elem_counter, 5 ) = l;
424 byte_counter += value_byte_sizes[a];
432 if ( byte_counter != buffer_bytes )
433 throw std::logic_error(
"Cabana::Grid::Halo::buildSteeringVector: "
434 "Steering vector contains different number "
435 "of bytes than buffer" );
436 if ( elem_counter != buffer_num_element )
437 throw std::logic_error(
"Cabana::Grid::Halo::buildSteeringVector: "
438 "Steering vector contains different number "
439 "of elements than buffer" );
442 Kokkos::deep_copy( steering.back(), host_steering );
446 template <std::
size_t NumArray>
449 const std::array<std::size_t, NumArray>& value_byte_sizes,
450 const int buffer_bytes,
const int buffer_num_element,
451 std::vector<Kokkos::View<int**, memory_space>>& steering )
457 Kokkos::create_mirror_view( Kokkos::HostSpace(), steering.back() );
458 int elem_counter = 0;
459 int byte_counter = 0;
460 for ( std::size_t a = 0; a < NumArray; ++a )
462 for (
int i = spaces[a].min( 0 ); i < spaces[a].max( 0 ); ++i )
464 for (
int j = spaces[a].min( 1 ); j < spaces[a].max( 1 ); ++j )
466 for (
int l = spaces[a].min( 2 ); l < spaces[a].max( 2 );
470 host_steering( elem_counter, 0 ) = byte_counter;
473 host_steering( elem_counter, 1 ) = a;
476 host_steering( elem_counter, 2 ) = i;
477 host_steering( elem_counter, 3 ) = j;
478 host_steering( elem_counter, 4 ) = l;
484 byte_counter += value_byte_sizes[a];
491 if ( byte_counter != buffer_bytes )
492 throw std::logic_error(
"Cabana::Grid::Halo::buildSteeringVector: "
493 "Steering vector contains different number "
494 "of bytes than buffer" );
495 if ( elem_counter != buffer_num_element )
496 throw std::logic_error(
"Cabana::Grid::Halo::buildSteeringVector: "
497 "Steering vector contains different number "
498 "of elements than buffer" );
501 Kokkos::deep_copy( steering.back(), host_steering );
506 template <
class ArrayView>
507 KOKKOS_INLINE_FUNCTION
static std::enable_if_t<4 == ArrayView::rank, void>
509 const Kokkos::View<int**, memory_space>& steering,
510 const int element_idx,
const ArrayView& array_view )
512 const char* elem_ptr =
reinterpret_cast<const char*
>( &array_view(
513 steering( element_idx, 2 ), steering( element_idx, 3 ),
514 steering( element_idx, 4 ), steering( element_idx, 5 ) ) );
515 for ( std::size_t b = 0; b <
sizeof(
typename ArrayView::value_type );
518 buffer( steering( element_idx, 0 ) + b ) = *( elem_ptr + b );
524 template <
class ArrayView>
525 KOKKOS_INLINE_FUNCTION
static std::enable_if_t<3 == ArrayView::rank, void>
527 const Kokkos::View<int**, memory_space>& steering,
528 const int element_idx,
const ArrayView& array_view )
530 const char* elem_ptr =
reinterpret_cast<const char*
>(
531 &array_view( steering( element_idx, 2 ), steering( element_idx, 3 ),
532 steering( element_idx, 4 ) ) );
533 for ( std::size_t b = 0; b <
sizeof(
typename ArrayView::value_type );
536 buffer( steering( element_idx, 0 ) + b ) = *( elem_ptr + b );
541 template <
class... ArrayViews>
542 KOKKOS_INLINE_FUNCTION
static void
543 packArray(
const Kokkos::View<char*, memory_space>& buffer,
544 const Kokkos::View<int**, memory_space>& steering,
545 const int element_idx,
546 const std::integral_constant<std::size_t, 0>,
550 if ( 0 == steering( element_idx, 1 ) )
556 template <std::size_t N,
class... ArrayViews>
557 KOKKOS_INLINE_FUNCTION
static void
558 packArray(
const Kokkos::View<char*, memory_space>& buffer,
559 const Kokkos::View<int**, memory_space>& steering,
560 const int element_idx,
561 const std::integral_constant<std::size_t, N>,
565 if ( N == steering( element_idx, 1 ) )
573 packArray( buffer, steering, element_idx,
574 std::integral_constant<std::size_t, N - 1>(), array_views );
578 template <
class ExecutionSpace,
class... ArrayViews>
580 const Kokkos::View<char*, memory_space>& buffer,
581 const Kokkos::View<int**, memory_space>& steering,
582 ArrayViews... array_views )
const
585 Kokkos::parallel_for(
586 "Cabana::Grid::Halo::pack_buffer",
587 Kokkos::RangePolicy<ExecutionSpace>( exec_space, 0,
588 steering.extent( 0 ) ),
589 KOKKOS_LAMBDA(
const int i ) {
592 std::integral_constant<std::size_t,
593 sizeof...( ArrayViews ) - 1>(),
601 KOKKOS_INLINE_FUNCTION
static void
604 array_val += buffer_val;
609 KOKKOS_INLINE_FUNCTION
static void
612 if ( buffer_val < array_val )
613 array_val = buffer_val;
618 KOKKOS_INLINE_FUNCTION
static void
621 if ( buffer_val > array_val )
622 array_val = buffer_val;
627 KOKKOS_INLINE_FUNCTION
static void
630 array_val = buffer_val;
635 template <
class ReduceOp,
class ArrayView>
636 KOKKOS_INLINE_FUNCTION
static std::enable_if_t<4 == ArrayView::rank, void>
638 const Kokkos::View<char*, memory_space>& buffer,
639 const Kokkos::View<int**, memory_space>& steering,
640 const int element_idx,
const ArrayView& array_view )
642 typename ArrayView::value_type elem;
643 char* elem_ptr =
reinterpret_cast<char*
>( &elem );
644 for ( std::size_t b = 0; b <
sizeof(
typename ArrayView::value_type );
647 *( elem_ptr + b ) = buffer( steering( element_idx, 0 ) + b );
650 array_view( steering( element_idx, 2 ),
651 steering( element_idx, 3 ),
652 steering( element_idx, 4 ),
653 steering( element_idx, 5 ) ) );
658 template <
class ReduceOp,
class ArrayView>
659 KOKKOS_INLINE_FUNCTION
static std::enable_if_t<3 == ArrayView::rank, void>
661 const Kokkos::View<char*, memory_space>& buffer,
662 const Kokkos::View<int**, memory_space>& steering,
663 const int element_idx,
const ArrayView& array_view )
665 typename ArrayView::value_type elem;
666 char* elem_ptr =
reinterpret_cast<char*
>( &elem );
667 for ( std::size_t b = 0; b <
sizeof(
typename ArrayView::value_type );
670 *( elem_ptr + b ) = buffer( steering( element_idx, 0 ) + b );
673 array_view( steering( element_idx, 2 ),
674 steering( element_idx, 3 ),
675 steering( element_idx, 4 ) ) );
679 template <
class ReduceOp,
class... ArrayViews>
680 KOKKOS_INLINE_FUNCTION
static void
682 const Kokkos::View<char*, memory_space>& buffer,
683 const Kokkos::View<int**, memory_space>& steering,
684 const int element_idx,
685 const std::integral_constant<std::size_t, 0>,
689 if ( 0 == steering( element_idx, 1 ) )
695 template <
class ReduceOp, std::size_t N,
class... ArrayViews>
696 KOKKOS_INLINE_FUNCTION
static void
698 const Kokkos::View<char*, memory_space>& buffer,
699 const Kokkos::View<int**, memory_space>& steering,
700 const int element_idx,
701 const std::integral_constant<std::size_t, N>,
705 if ( N == steering( element_idx, 1 ) )
713 unpackArray( reduce_op, buffer, steering, element_idx,
714 std::integral_constant<std::size_t, N - 1>(),
719 template <
class ExecutionSpace,
class ReduceOp,
class... ArrayViews>
721 const ExecutionSpace& exec_space,
722 const Kokkos::View<char*, memory_space>& buffer,
723 const Kokkos::View<int**, memory_space>& steering,
724 ArrayViews... array_views )
const
727 Kokkos::parallel_for(
728 "Cabana::Grid::Halo::unpack_buffer",
729 Kokkos::RangePolicy<ExecutionSpace>( exec_space, 0,
730 steering.extent( 0 ) ),
731 KOKKOS_LAMBDA(
const int i ) {
733 reduce_op, buffer, steering, i,
734 std::integral_constant<std::size_t,
735 sizeof...( ArrayViews ) - 1>(),
767template <
class ArrayT,
class... Types>
771 using type =
typename ArrayT::memory_space;
775template <
class MemorySpace,
class CommSpaceType = Mpi>
782#ifdef Cabana_ENABLE_MPI
799template <
class CommSpaceType =
Mpi,
class Pattern,
class... ArrayTypes>
801 const ArrayTypes&... arrays )
804 return std::make_shared<Halo<memory_space, CommSpaceType>>( pattern, width,
auto createHalo(const Pattern &pattern, const int width, const ArrayTypes &... arrays)
Halo creation function.
Definition Cabana_Grid_HaloBase.hpp:800
Multi-node grid scatter/gather implemented with vanilla MPI.
Pack variadic template parameters for device capture.
Definition Cabana_Grid_HaloBase.hpp:123
std::vector< int > _neighbor_ranks
The ranks we will send/receive from.
Definition Cabana_Grid_HaloBase.hpp:742
std::vector< Kokkos::View< char *, memory_space > > _ghosted_buffers
For each neighbor, send/receive buffers for data we ghost.
Definition Cabana_Grid_HaloBase.hpp:754
void packBuffer(const ExecutionSpace &exec_space, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, ArrayViews... array_views) const
Pack arrays into a buffer.
Definition Cabana_Grid_HaloBase.hpp:579
static KOKKOS_INLINE_FUNCTION void packArray(const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const std::integral_constant< std::size_t, N >, const Cabana::ParameterPack< ArrayViews... > &array_views)
Pack an array into a buffer.
Definition Cabana_Grid_HaloBase.hpp:558
void buildSteeringVector(const std::array< IndexSpace< 4 >, NumArray > &spaces, const std::array< std::size_t, NumArray > &value_byte_sizes, const int buffer_bytes, const int buffer_num_element, std::vector< Kokkos::View< int **, memory_space > > &steering)
Build 3d steering vector.
Definition Cabana_Grid_HaloBase.hpp:383
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_HaloBase.hpp:610
std::vector< int > _receive_tags
The tag we use for receiving from each neighbor.
Definition Cabana_Grid_HaloBase.hpp:748
MPI_Comm getComm(const Array_t &array, const ArrayTypes &... arrays) const
Get the communicator and check to make sure all are the same.
Definition Cabana_Grid_HaloBase.hpp:292
static KOKKOS_INLINE_FUNCTION void unpackArray(const ReduceOp reduce_op, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const std::integral_constant< std::size_t, N >, const Cabana::ParameterPack< ArrayViews... > &array_views)
Unpack an array from a buffer.
Definition Cabana_Grid_HaloBase.hpp:697
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_HaloBase.hpp:619
HaloBase(const Pattern &pattern, const int width, const ArrayTypes &... arrays)
Constructor.
Definition Cabana_Grid_HaloBase.hpp:220
static KOKKOS_INLINE_FUNCTION void packArray(const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const std::integral_constant< std::size_t, 0 >, const Cabana::ParameterPack< ArrayViews... > &array_views)
Pack an array into a buffer.
Definition Cabana_Grid_HaloBase.hpp:543
static KOKKOS_INLINE_FUNCTION std::enable_if_t< 4==ArrayView::rank, void > unpackElement(const ReduceOp &reduce_op, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const ArrayView &array_view)
Definition Cabana_Grid_HaloBase.hpp:637
std::vector< Kokkos::View< int **, memory_space > > _owned_steering
For each neighbor, steering vector for the owned buffer.
Definition Cabana_Grid_HaloBase.hpp:757
void buildSteeringVector(const std::array< IndexSpace< 3 >, NumArray > &spaces, const std::array< std::size_t, NumArray > &value_byte_sizes, const int buffer_bytes, const int buffer_num_element, std::vector< Kokkos::View< int **, memory_space > > &steering)
Build 2d steering vector.
Definition Cabana_Grid_HaloBase.hpp:447
void buildCommData(DecompositionTag decomposition_tag, const int width, const std::array< int, NumSpaceDim > &nid, std::vector< Kokkos::View< char *, memory_space > > &buffers, std::vector< Kokkos::View< int **, memory_space > > &steering, const ArrayTypes &... arrays)
Build communication data.
Definition Cabana_Grid_HaloBase.hpp:338
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_HaloBase.hpp:602
static KOKKOS_INLINE_FUNCTION std::enable_if_t< 3==ArrayView::rank, void > packElement(const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const ArrayView &array_view)
Definition Cabana_Grid_HaloBase.hpp:526
MemorySpace memory_space
Memory space.
Definition Cabana_Grid_HaloBase.hpp:207
static KOKKOS_INLINE_FUNCTION std::enable_if_t< 3==ArrayView::rank, void > unpackElement(const ReduceOp &reduce_op, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const ArrayView &array_view)
Definition Cabana_Grid_HaloBase.hpp:660
auto getLocalGrid(const Array_t &array, const ArrayTypes &... arrays)
Definition Cabana_Grid_HaloBase.hpp:318
std::vector< Kokkos::View< char *, memory_space > > _owned_buffers
For each neighbor, send/receive buffers for data we own.
Definition Cabana_Grid_HaloBase.hpp:751
static KOKKOS_INLINE_FUNCTION void unpackArray(const ReduceOp &reduce_op, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const std::integral_constant< std::size_t, 0 >, const Cabana::ParameterPack< ArrayViews... > &array_views)
Unpack an array from a buffer.
Definition Cabana_Grid_HaloBase.hpp:681
std::vector< Kokkos::View< int **, memory_space > > _ghosted_steering
For each neighbor, steering vector for the ghosted buffer.
Definition Cabana_Grid_HaloBase.hpp:760
void unpackBuffer(const ReduceOp &reduce_op, const ExecutionSpace &exec_space, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, ArrayViews... array_views) const
Unpack arrays from a buffer.
Definition Cabana_Grid_HaloBase.hpp:720
MPI_Comm getComm(const Array_t &array) const
Get the communicator.
Definition Cabana_Grid_HaloBase.hpp:285
static KOKKOS_INLINE_FUNCTION std::enable_if_t< 4==ArrayView::rank, void > packElement(const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const ArrayView &array_view)
Definition Cabana_Grid_HaloBase.hpp:508
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_HaloBase.hpp:628
std::vector< int > _send_tags
The tag we use for sending to each neighbor.
Definition Cabana_Grid_HaloBase.hpp:745
auto getLocalGrid(const Array_t &array)
Definition Cabana_Grid_HaloBase.hpp:310
Base halo exchange pattern class.
Definition Cabana_Grid_HaloBase.hpp:48
std::vector< std::array< int, num_space_dim > > getNeighbors() const
Get the neighbors that are in the halo pattern.
Definition Cabana_Grid_HaloBase.hpp:67
void setNeighbors(const std::vector< std::array< int, num_space_dim > > &neighbors)
Assign the neighbors that are in the halo pattern.
Definition Cabana_Grid_HaloBase.hpp:61
static constexpr std::size_t num_space_dim
Spatial dimension.
Definition Cabana_Grid_HaloBase.hpp:51
Definition Cabana_Grid_HaloBase.hpp:776
Structured index space.
Definition Cabana_Grid_IndexSpace.hpp:37
Definition Cabana_Grid_HaloBase.hpp:79
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
ParameterPack< Types... > makeParameterPack(const Types &... ts)
Create a parameter pack.
Definition Cabana_ParameterPack.hpp:189
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
Grid and particle-grid data structures and algorithms.
Infer array memory space.
Definition Cabana_Grid_HaloBase.hpp:769
typename ArrayT::memory_space type
Memory space.
Definition Cabana_Grid_HaloBase.hpp:771
Ghosted decomposition tag.
Definition Cabana_Grid_Types.hpp:197
Owned decomposition tag.
Definition Cabana_Grid_Types.hpp:190
Definition Cabana_Grid_HaloBase.hpp:176
Definition Cabana_Grid_HaloBase.hpp:170
Definition Cabana_Grid_HaloBase.hpp:184
Sum values from neighboring ranks into this rank's data.
Definition Cabana_Grid_HaloBase.hpp:164
Vanilla MPI backend tag - default.
Definition Cabana_Tags.hpp:28
Definition Cabana_ParameterPack.hpp:86