Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana_Grid_HaloBase.hpp
Go to the documentation of this file.
1/****************************************************************************
2 * Copyright (c) 2018-2023 by the Cabana authors *
3 * All rights reserved. *
4 * *
5 * This file is part of the Cabana library. Cabana is distributed under a *
6 * BSD 3-clause license. For the licensing terms see the LICENSE file in *
7 * the top-level directory. *
8 * *
9 * SPDX-License-Identifier: BSD-3-Clause *
10 ****************************************************************************/
11
16#ifndef CABANA_GRID_HALO_HPP
17#define CABANA_GRID_HALO_HPP
18
19#include <Cabana_Grid_Array.hpp>
21#include <Cabana_Tags.hpp>
22
23#include <Cabana_Core_Config.hpp>
25
26#include <Kokkos_Core.hpp>
27#include <Kokkos_Profiling_ScopedRegion.hpp>
28
29#include <mpi.h>
30
31#include <algorithm>
32#include <array>
33#include <cmath>
34#include <type_traits>
35#include <vector>
36
37namespace Cabana
38{
39namespace Grid
40{
41
42//---------------------------------------------------------------------------//
43// Halo exchange patterns.
44//---------------------------------------------------------------------------//
46template <std::size_t NumSpaceDim>
47class HaloPattern
48{
49 public:
51 static constexpr std::size_t num_space_dim = NumSpaceDim;
52
53 // Default constructor.
54 HaloPattern() {}
55
56 // Destructor
57 virtual ~HaloPattern() = default;
58
60 void
61 setNeighbors( const std::vector<std::array<int, num_space_dim>>& neighbors )
62 {
63 _neighbors = neighbors;
64 }
65
67 std::vector<std::array<int, num_space_dim>> getNeighbors() const
68 {
69 return _neighbors;
70 }
71
72 private:
73 std::vector<std::array<int, num_space_dim>> _neighbors;
74};
75
78template <std::size_t NumSpaceDim>
80
83template <>
84class NodeHaloPattern<3> : public HaloPattern<3>
85{
86 public:
87 NodeHaloPattern()
88 : HaloPattern<3>()
89 {
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 } );
97 this->setNeighbors( neighbors );
98 }
99};
100
103template <>
104class NodeHaloPattern<2> : public HaloPattern<2>
105{
106 public:
107 NodeHaloPattern()
108 : HaloPattern<2>()
109 {
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 } );
116 this->setNeighbors( neighbors );
117 }
118};
119
122template <std::size_t NumSpaceDim>
124
127template <>
128class FaceHaloPattern<3> : public HaloPattern<3>
129{
130 public:
131 FaceHaloPattern()
132 : HaloPattern<3>()
133 {
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 } };
137 this->setNeighbors( neighbors );
138 }
139};
140
143template <>
144class FaceHaloPattern<2> : public HaloPattern<2>
145{
146 public:
147 FaceHaloPattern()
148 : HaloPattern<2>()
149 {
150 std::vector<std::array<int, 2>> neighbors = {
151 { -1, 0 }, { 1, 0 }, { 0, -1 }, { 0, 1 } };
152 this->setNeighbors( neighbors );
153 }
154};
155
156//---------------------------------------------------------------------------//
157// Scatter reduction.
158//---------------------------------------------------------------------------//
159namespace ScatterReduce
160{
161
163struct Sum
164{
165};
166
169struct Min
170{
171};
172
175struct Max
176{
177};
178
184{
185};
186
187} // end namespace ScatterReduce
188
189//---------------------------------------------------------------------------//
190// Halo
191// ---------------------------------------------------------------------------//
202template <class MemorySpace>
204{
205 public:
207 using memory_space = MemorySpace;
208
209 protected:
219 template <class Pattern, class... ArrayTypes>
220 HaloBase( const Pattern& pattern, const int width,
221 const ArrayTypes&... arrays )
222 {
223 // Spatial dimension.
224 const std::size_t num_space_dim = Pattern::num_space_dim;
225
226 // Get the local grid.
227 auto local_grid = getLocalGrid( arrays... );
228
229 // Function to get the local id of the neighbor.
230 auto neighbor_id = []( const std::array<int, num_space_dim>& ijk )
231 {
232 int id = ijk[0];
233 for ( std::size_t d = 1; d < num_space_dim; ++d )
234 id += num_space_dim * id + ijk[d];
235 return id;
236 };
237
238 // Neighbor id flip function. This lets us compute what neighbor we
239 // are relative to a given neighbor.
240 auto flip_id = [=]( const std::array<int, num_space_dim>& ijk )
241 {
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];
245 return flip_ijk;
246 };
247
248 // Get the neighbor ranks we will exchange with in the halo and
249 // allocate buffers. If any of the exchanges are self sends mark these
250 // so we know which send buffers correspond to which receive buffers.
251 auto neighbors = pattern.getNeighbors();
252 for ( const auto& n : neighbors )
253 {
254 // Get the rank of the neighbor.
255 int rank = local_grid->neighborRank( n );
256
257 // If this is a valid rank add it as a neighbor.
258 if ( rank >= 0 )
259 {
260 // Add the rank.
261 _neighbor_ranks.push_back( rank );
262
263 // Set the tag we will use to send data to this neighbor. The
264 // receiving rank should have a matching tag.
265 _send_tags.push_back( neighbor_id( n ) );
266
267 // Set the tag we will use to receive data from this
268 // neighbor. The sending rank should have a matching tag.
269 _receive_tags.push_back( neighbor_id( flip_id( n ) ) );
270
271 // Create communication data for owned entities.
273 arrays... );
274
275 // Create communication data for ghosted entities.
277 _ghosted_steering, arrays... );
278 }
279 }
280 }
281
282 public:
284 template <class Array_t>
285 MPI_Comm getComm( const Array_t& array ) const
286 {
287 return array.layout()->localGrid()->globalGrid().comm();
288 }
289
291 template <class Array_t, class... ArrayTypes>
292 MPI_Comm getComm( const Array_t& array, const ArrayTypes&... arrays ) const
293 {
294 auto comm = getComm( array );
295
296 // Check that the communicator of this array is the same as the other
297 // arrays.
298 int result;
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" );
303
304 return comm;
305 }
306
309 template <class Array_t>
310 auto getLocalGrid( const Array_t& array )
311 {
312 return array.layout()->localGrid();
313 }
314
317 template <class Array_t, class... ArrayTypes>
318 auto getLocalGrid( const Array_t& array, const ArrayTypes&... arrays )
319 {
320 // Recurse.
321 auto local_grid = getLocalGrid( arrays... );
322
323 // Check that the halo sizes same.
324 if ( local_grid->haloCellWidth() !=
325 array.layout()->localGrid()->haloCellWidth() )
326 {
327 throw std::runtime_error( "Cabana::Grid::Halo::getlocalGrid: "
328 "Arrays have different halo widths" );
329 }
330
331 return local_grid;
332 }
333
335 template <class DecompositionTag, std::size_t NumSpaceDim,
336 class... ArrayTypes>
337 void
338 buildCommData( DecompositionTag decomposition_tag, const int width,
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 )
343 {
344 // Number of arrays.
345 const std::size_t num_array = sizeof...( ArrayTypes );
346
347 // Get the byte sizes of array value types.
348 std::array<std::size_t, num_array> value_byte_sizes = {
349 sizeof( typename ArrayTypes::value_type )... };
350
351 // Get the index spaces we share with this neighbor. We
352 // get a shared index space for each array.
353 std::array<IndexSpace<NumSpaceDim + 1>, num_array> spaces = {
354 ( arrays.layout()->sharedIndexSpace( decomposition_tag, nid,
355 width ) )... };
356
357 // Compute the buffer size of this neighbor and the
358 // number of elements in the buffer.
359 int buffer_bytes = 0;
360 int buffer_num_element = 0;
361 for ( std::size_t a = 0; a < num_array; ++a )
362 {
363 buffer_bytes += value_byte_sizes[a] * spaces[a].size();
364 buffer_num_element += spaces[a].size();
365 }
366
367 // Allocate the buffer of data that we share with this neighbor. All
368 // arrays will be packed into a single buffer.
369 buffers.push_back(
370 Kokkos::View<char*, memory_space>( "halo_buffer", buffer_bytes ) );
371
372 // Allocate the steering vector for building the buffer.
373 steering.push_back( Kokkos::View<int**, memory_space>(
374 "steering", buffer_num_element, 3 + NumSpaceDim ) );
375
376 // Build steering vector.
377 buildSteeringVector( spaces, value_byte_sizes, buffer_bytes,
378 buffer_num_element, steering );
379 }
380
382 template <std::size_t NumArray>
384 const std::array<IndexSpace<4>, NumArray>& spaces,
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 )
388 {
389 // Create the steering vector. For each element in the buffer it gives
390 // the starting byte location of the element, the array the element is
391 // in, and the ijkl structured index in the array of the element.
392 auto host_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 )
397 {
398 for ( int i = spaces[a].min( 0 ); i < spaces[a].max( 0 ); ++i )
399 {
400 for ( int j = spaces[a].min( 1 ); j < spaces[a].max( 1 ); ++j )
401 {
402 for ( int k = spaces[a].min( 2 ); k < spaces[a].max( 2 );
403 ++k )
404 {
405 for ( int l = spaces[a].min( 3 );
406 l < spaces[a].max( 3 ); ++l )
407 {
408 // Byte starting location in buffer.
409 host_steering( elem_counter, 0 ) = byte_counter;
410
411 // Array location of element.
412 host_steering( elem_counter, 1 ) = a;
413
414 // Structured index in array of element.
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;
419
420 // Update element id.
421 ++elem_counter;
422
423 // Update buffer position.
424 byte_counter += value_byte_sizes[a];
425 }
426 }
427 }
428 }
429 }
430
431 // Check that all elements and bytes are accounted for.
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" );
440
441 // Copy steering vector to device.
442 Kokkos::deep_copy( steering.back(), host_steering );
443 }
444
446 template <std::size_t NumArray>
448 const std::array<IndexSpace<3>, NumArray>& spaces,
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 )
452 {
453 // Create the steering vector. For each element in the buffer it gives
454 // the starting byte location of the element, the array the element is
455 // in, and the ijkl structured index in the array of the element.
456 auto host_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 )
461 {
462 for ( int i = spaces[a].min( 0 ); i < spaces[a].max( 0 ); ++i )
463 {
464 for ( int j = spaces[a].min( 1 ); j < spaces[a].max( 1 ); ++j )
465 {
466 for ( int l = spaces[a].min( 2 ); l < spaces[a].max( 2 );
467 ++l )
468 {
469 // Byte starting location in buffer.
470 host_steering( elem_counter, 0 ) = byte_counter;
471
472 // Array location of element.
473 host_steering( elem_counter, 1 ) = a;
474
475 // Structured index in array of element.
476 host_steering( elem_counter, 2 ) = i;
477 host_steering( elem_counter, 3 ) = j;
478 host_steering( elem_counter, 4 ) = l;
479
480 // Update element id.
481 ++elem_counter;
482
483 // Update buffer position.
484 byte_counter += value_byte_sizes[a];
485 }
486 }
487 }
488 }
489
490 // Check that all elements and bytes are accounted for.
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" );
499
500 // Copy steering vector to device.
501 Kokkos::deep_copy( steering.back(), host_steering );
502 }
503
506 template <class ArrayView>
507 KOKKOS_INLINE_FUNCTION static std::enable_if_t<4 == ArrayView::rank, void>
508 packElement( const Kokkos::View<char*, memory_space>& buffer,
509 const Kokkos::View<int**, memory_space>& steering,
510 const int element_idx, const ArrayView& array_view )
511 {
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 );
516 ++b )
517 {
518 buffer( steering( element_idx, 0 ) + b ) = *( elem_ptr + b );
519 }
520 }
521
524 template <class ArrayView>
525 KOKKOS_INLINE_FUNCTION static std::enable_if_t<3 == ArrayView::rank, void>
526 packElement( const Kokkos::View<char*, memory_space>& buffer,
527 const Kokkos::View<int**, memory_space>& steering,
528 const int element_idx, const ArrayView& array_view )
529 {
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 );
534 ++b )
535 {
536 buffer( steering( element_idx, 0 ) + b ) = *( elem_ptr + b );
537 }
538 }
539
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>,
547 const Cabana::ParameterPack<ArrayViews...>& array_views )
548 {
549 // If the pack element_idx is in the current array, pack it.
550 if ( 0 == steering( element_idx, 1 ) )
551 packElement( buffer, steering, element_idx,
552 Cabana::get<0>( array_views ) );
553 }
554
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>,
562 const Cabana::ParameterPack<ArrayViews...>& array_views )
563 {
564 // If the pack element_idx is in the current array, pack it.
565 if ( N == steering( element_idx, 1 ) )
566 {
567 packElement( buffer, steering, element_idx,
568 Cabana::get<N>( array_views ) );
569 return;
570 }
571
572 // Recurse.
573 packArray( buffer, steering, element_idx,
574 std::integral_constant<std::size_t, N - 1>(), array_views );
575 }
576
578 template <class ExecutionSpace, class... ArrayViews>
579 void packBuffer( const ExecutionSpace& exec_space,
580 const Kokkos::View<char*, memory_space>& buffer,
581 const Kokkos::View<int**, memory_space>& steering,
582 ArrayViews... array_views ) const
583 {
584 auto pp = Cabana::makeParameterPack( array_views... );
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 ) {
590 packArray(
591 buffer, steering, i,
592 std::integral_constant<std::size_t,
593 sizeof...( ArrayViews ) - 1>(),
594 pp );
595 } );
596 exec_space.fence();
597 }
598
600 template <class T>
601 KOKKOS_INLINE_FUNCTION static void
602 unpackOp( ScatterReduce::Sum, const T& buffer_val, T& array_val )
603 {
604 array_val += buffer_val;
605 }
606
608 template <class T>
609 KOKKOS_INLINE_FUNCTION static void
610 unpackOp( ScatterReduce::Min, const T& buffer_val, T& array_val )
611 {
612 if ( buffer_val < array_val )
613 array_val = buffer_val;
614 }
615
617 template <class T>
618 KOKKOS_INLINE_FUNCTION static void
619 unpackOp( ScatterReduce::Max, const T& buffer_val, T& array_val )
620 {
621 if ( buffer_val > array_val )
622 array_val = buffer_val;
623 }
624
626 template <class T>
627 KOKKOS_INLINE_FUNCTION static void
628 unpackOp( ScatterReduce::Replace, const T& buffer_val, T& array_val )
629 {
630 array_val = buffer_val;
631 }
632
635 template <class ReduceOp, class ArrayView>
636 KOKKOS_INLINE_FUNCTION static std::enable_if_t<4 == ArrayView::rank, void>
637 unpackElement( const ReduceOp& reduce_op,
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 )
641 {
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 );
645 ++b )
646 {
647 *( elem_ptr + b ) = buffer( steering( element_idx, 0 ) + b );
648 }
649 unpackOp( reduce_op, elem,
650 array_view( steering( element_idx, 2 ),
651 steering( element_idx, 3 ),
652 steering( element_idx, 4 ),
653 steering( element_idx, 5 ) ) );
654 }
655
658 template <class ReduceOp, class ArrayView>
659 KOKKOS_INLINE_FUNCTION static std::enable_if_t<3 == ArrayView::rank, void>
660 unpackElement( const ReduceOp& reduce_op,
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 )
664 {
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 );
668 ++b )
669 {
670 *( elem_ptr + b ) = buffer( steering( element_idx, 0 ) + b );
671 }
672 unpackOp( reduce_op, elem,
673 array_view( steering( element_idx, 2 ),
674 steering( element_idx, 3 ),
675 steering( element_idx, 4 ) ) );
676 }
677
679 template <class ReduceOp, class... ArrayViews>
680 KOKKOS_INLINE_FUNCTION static void
681 unpackArray( const ReduceOp& reduce_op,
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>,
686 const Cabana::ParameterPack<ArrayViews...>& array_views )
687 {
688 // If the unpack element_idx is in the current array, unpack it.
689 if ( 0 == steering( element_idx, 1 ) )
690 unpackElement( reduce_op, buffer, steering, element_idx,
691 Cabana::get<0>( array_views ) );
692 }
693
695 template <class ReduceOp, std::size_t N, class... ArrayViews>
696 KOKKOS_INLINE_FUNCTION static void
697 unpackArray( const ReduceOp reduce_op,
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>,
702 const Cabana::ParameterPack<ArrayViews...>& array_views )
703 {
704 // If the unpack element_idx is in the current array, unpack it.
705 if ( N == steering( element_idx, 1 ) )
706 {
707 unpackElement( reduce_op, buffer, steering, element_idx,
708 Cabana::get<N>( array_views ) );
709 return;
710 }
711
712 // Recurse.
713 unpackArray( reduce_op, buffer, steering, element_idx,
714 std::integral_constant<std::size_t, N - 1>(),
715 array_views );
716 }
717
719 template <class ExecutionSpace, class ReduceOp, class... ArrayViews>
720 void unpackBuffer( const ReduceOp& reduce_op,
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
725 {
726 auto pp = Cabana::makeParameterPack( array_views... );
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>(),
736 pp );
737 } );
738 }
739
740 protected:
742 std::vector<int> _neighbor_ranks;
743
745 std::vector<int> _send_tags;
746
748 std::vector<int> _receive_tags;
749
751 std::vector<Kokkos::View<char*, memory_space>> _owned_buffers;
752
754 std::vector<Kokkos::View<char*, memory_space>> _ghosted_buffers;
755
757 std::vector<Kokkos::View<int**, memory_space>> _owned_steering;
758
760 std::vector<Kokkos::View<int**, memory_space>> _ghosted_steering;
761};
762
763//---------------------------------------------------------------------------//
764// Creation function.
765//---------------------------------------------------------------------------//
767template <class ArrayT, class... Types>
769{
771 using type = typename ArrayT::memory_space;
772};
773
774// Forward declaration of the primary grid Halo template.
775template <class MemorySpace, class CommSpaceType = Mpi>
776class Halo;
777
778} // end namespace Grid
779} // end namespace Cabana
780
781// Include communication backends from what is enabled in CMake.
782#ifdef Cabana_ENABLE_MPI
784#endif // Enable MPI
785
786namespace Cabana
787{
788namespace Grid
789{
790
791//---------------------------------------------------------------------------//
799template <class CommSpaceType = Mpi, class Pattern, class... ArrayTypes>
800auto createHalo( const Pattern& pattern, const int width,
801 const ArrayTypes&... arrays )
802{
803 using memory_space = typename ArrayPackMemorySpace<ArrayTypes...>::type;
804 return std::make_shared<Halo<memory_space, CommSpaceType>>( pattern, width,
805 arrays... );
806}
807
808} // end namespace Grid
809} // end namespace Cabana
810
811//---------------------------------------------------------------------------//
812
813#endif // end CABANA_GRID_HALOBASE_HPP
Grid field arrays.
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.
Logical grid indexing.
Pack variadic template parameters for device capture.
Type tags used in Cabana.
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