3. LvArray::Array
¶
The LvArray::Array
implements a multidimensional array of values. Unlike Kokkos mdspan or the RAJA View the LvArray::Array
owns the allocation it is associated with and supports slicing with operator[]
in addition to the standard operator()
. Further more a one dimensional LvArray::Array
supports operations such as emplace_back
and erase
with functionality similar to std::vector
.
3.1. Template arguments¶
The LvArray::Array
requires five template arguments.
T
: The type of values stored in the array.NDIM
: The number of dimensionality of the array or the number of indices required to access a value.PERMUTATION
: Acamp::idx_seq
which describes the mapping from the multidimensional index space to a linear index. Must be of lengthNDIM
and contain all the values between0
andNDIM - 1
. The way to read a permutation is that the indices go from the slowest on the left to the fastest on the right. Equivalently the left most index has the largest stride whereas the right most index has unit stride.INDEX_TYPE
: An integral type used in index calculations, the suggested type isstd::ptrdiff_t
.BUFFER_TYPE
: A template template parameter specifying the buffer type used for allocation and de-allocation, theLvArray::Array
contains aBUFFER_TYPE< T >
.
Note
LvArray
uses the same permutation conventions as the RAJA::View
, in fact it is recommended to use the types defined in RAJA/Permutations.hpp.
3.2. Creating and accessing a LvArray::Array
¶
The LvArray::Array
has two primary constructors a default constructor which creates an empty array and a constructor that takes the size of each dimension. When using the sized constructor the values of the array are default initialized.
TEST( Array, constructors )
{
{
// Create an empty 2D array of integers.
LvArray::Array< int,
2,
camp::idx_seq< 0, 1 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array;
EXPECT_TRUE( array.empty() );
EXPECT_EQ( array.size(), 0 );
EXPECT_EQ( array.size( 0 ), 0 );
EXPECT_EQ( array.size( 1 ), 0 );
}
{
// Create a 3D array of std::string of size 3 x 4 x 5.
LvArray::Array< std::string,
3,
camp::idx_seq< 0, 1, 2 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array( 3, 4, 5 );
EXPECT_FALSE( array.empty() );
EXPECT_EQ( array.size(), 3 * 4 * 5 );
EXPECT_EQ( array.size( 0 ), 3 );
EXPECT_EQ( array.size( 1 ), 4 );
EXPECT_EQ( array.size( 2 ), 5 );
// The values are default initialized.
std::string const * const values = array.data();
for( std::ptrdiff_t i = 0; i < array.size(); ++i )
{
EXPECT_EQ( values[ i ], std::string() );
}
}
}
[Source: examples/exampleArray.cpp]
LvArray::Array
supports two indexing methods operator()
which takes all of the indices to a value and operator[]
which takes a single index at a time and can be chained together.
TEST( Array, accessors )
{
// Create a 2D array of integers.
LvArray::Array< int,
2,
camp::idx_seq< 0, 1 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array( 3, 4 );
// Access using operator().
for( std::ptrdiff_t i = 0; i < array.size( 0 ); ++i )
{
for( std::ptrdiff_t j = 0; j < array.size( 1 ); ++j )
{
array( i, j ) = array.size( 1 ) * i + j;
}
}
// Access using operator[].
for( std::ptrdiff_t i = 0; i < array.size( 0 ); ++i )
{
for( std::ptrdiff_t j = 0; j < array.size( 1 ); ++j )
{
EXPECT_EQ( array[ i ][ j ], array.size( 1 ) * i + j );
}
}
}
[Source: examples/exampleArray.cpp]
The two indexing methods work consistently regardless of how the data is layed out in memory.
TEST( Array, permutations )
{
{
// Create a 3D array of doubles in the standard layout.
LvArray::Array< int,
3,
camp::idx_seq< 0, 1, 2 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array( 3, 4, 5 );
// Index 0 has the largest stride while index 2 has unit stride.
EXPECT_EQ( array.strides()[ 0 ], array.size( 2 ) * array.size( 1 ) );
EXPECT_EQ( array.strides()[ 1 ], array.size( 2 ) );
EXPECT_EQ( array.strides()[ 2 ], 1 );
int const * const pointer = array.data();
for( std::ptrdiff_t i = 0; i < array.size( 0 ); ++i )
{
for( std::ptrdiff_t j = 0; j < array.size( 1 ); ++j )
{
for( std::ptrdiff_t k = 0; k < array.size( 2 ); ++k )
{
std::ptrdiff_t const offset = array.size( 2 ) * array.size( 1 ) * i +
array.size( 2 ) * j + k;
EXPECT_EQ( &array( i, j, k ), pointer + offset );
}
}
}
}
{
// Create a 3D array of doubles in a flipped layout.
LvArray::Array< int,
3,
camp::idx_seq< 2, 1, 0 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array( 3, 4, 5 );
// Index 0 has the unit stride while index 2 has the largest stride.
EXPECT_EQ( array.strides()[ 0 ], 1 );
EXPECT_EQ( array.strides()[ 1 ], array.size( 0 ) );
EXPECT_EQ( array.strides()[ 2 ], array.size( 0 ) * array.size( 1 ) );
int const * const pointer = array.data();
for( std::ptrdiff_t i = 0; i < array.size( 0 ); ++i )
{
for( std::ptrdiff_t j = 0; j < array.size( 1 ); ++j )
{
for( std::ptrdiff_t k = 0; k < array.size( 2 ); ++k )
{
std::ptrdiff_t const offset = i + array.size( 0 ) * j +
array.size( 0 ) * array.size( 1 ) * k;
EXPECT_EQ( &array[ i ][ j ][ k ], pointer + offset );
}
}
}
}
}
[Source: examples/exampleArray.cpp]
3.3. Resizing a LvArray::Array
¶
LvArray::Array
supports a multitude of resizing options. You can resize all of the dimensions at once or specify a set of dimensions to resize. These methods default initialize newly created values and destroy any values no longer in the Array.
TEST( Array, resize )
{
LvArray::Array< int,
3,
camp::idx_seq< 0, 1, 2 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array;
// Resize using a pointer
std::ptrdiff_t const sizes[ 3 ] = { 2, 5, 6 };
array.resize( 3, sizes );
EXPECT_EQ( array.size(), 2 * 5 * 6 );
EXPECT_EQ( array.size( 0 ), 2 );
EXPECT_EQ( array.size( 1 ), 5 );
EXPECT_EQ( array.size( 2 ), 6 );
// Resizing using a variadic parameter pack.
array.resize( 3, 4, 2 );
EXPECT_EQ( array.size(), 3 * 4 * 2 );
EXPECT_EQ( array.size( 0 ), 3 );
EXPECT_EQ( array.size( 1 ), 4 );
EXPECT_EQ( array.size( 2 ), 2 );
// Resize the second and third dimensions
array.resizeDimension< 1, 2 >( 3, 6 );
EXPECT_EQ( array.size(), 3 * 3 * 6 );
EXPECT_EQ( array.size( 0 ), 3 );
EXPECT_EQ( array.size( 1 ), 3 );
EXPECT_EQ( array.size( 2 ), 6 );
}
[Source: examples/exampleArray.cpp]
LvArray::Array
also has a method resizeWithoutInitializationOrDestruction
that is only enabled if the value type of the array T
is trivially destructible. This method does not initialize new values or destroy old values and as such it can be much faster for large allocations of trivial types.
It is important to note that unless the array being resized is one dimensional the resize methods above do not preserve the values in the array. That is if you have a two dimensional array A
of size and you resize it to using any of the methods above then you cannot rely on A( i, j )
having the same value it did before the resize.
There is also a method resize
which takes a single parameter and will resize the dimension given by getSingleParameterResizeIndex
. Unlike the previous methods this will preserve the values in the array. By default the first dimension is resized but you can choose the dimension with setSingleParameterResizeIndex
.
TEST( Array, resizeSingleDimension )
{
LvArray::Array< int,
2,
camp::idx_seq< 1, 0 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array( 5, 6 );
for( std::ptrdiff_t i = 0; i < array.size( 0 ); ++i )
{
for( std::ptrdiff_t j = 0; j < array.size( 1 ); ++j )
{
array( i, j ) = 6 * i + j;
}
}
// Grow the first dimension from 5 to 8.
array.resize( 8 );
for( std::ptrdiff_t i = 0; i < array.size( 0 ); ++i )
{
for( std::ptrdiff_t j = 0; j < array.size( 1 ); ++j )
{
if( i < 5 )
{
EXPECT_EQ( array( i, j ), 6 * i + j );
}
else
{
EXPECT_EQ( array( i, j ), 0 );
}
}
}
// Shrink the second dimension from 6 to 3;
array.setSingleParameterResizeIndex( 1 );
array.resize( 3 );
for( std::ptrdiff_t i = 0; i < array.size( 0 ); ++i )
{
for( std::ptrdiff_t j = 0; j < array.size( 1 ); ++j )
{
if( i < 5 )
{
EXPECT_EQ( array( i, j ), 6 * i + j );
}
else
{
EXPECT_EQ( array( i, j ), 0 );
}
}
}
}
[Source: examples/exampleArray.cpp]
The single dimension resize should only be used when it is necessary to preserve the values as it is a much more complicated operation than the multi-dimension resize methods.
3.4. The one dimensional LvArray::Array
¶
The one dimensional LvArray::Array
supports a couple methods that are not available to multidimensional arrays. These methods are emplace_back
, emplace
, insert
, pop_back
and erase
. They all behave exactly like their std::vector
counter part, the only difference being that emplace
, insert
and erase
take an integer specifying the position to perform the operation instead of an iterator.
3.5. Lambda capture and LvArray::ArrayView
¶
LvArray::ArrayView
is the parent class and the view class of LvArray::Array
. It shares the same template parameters as LvArray::Array
except that the PERMUTATION
type is replaced by an integer corresponding to the unit stride dimension (the last entry in the permutation).
There are multiple ways to create an LvArray::ArrayView
from an LvArray::Array
. Because it is the parent class you can just create a reference to a LvArray::ArrayView
from an existing LvArray::Array
or by using the method toView
. LvArray::Array
also has a user defined conversion operator to a LvArray::ArrayView
with a const value type or you can use the toViewConst
method.
The LvArray::ArrayView
has a copy constructor, move constructor, copy assignment operator and move assignment operator all of which perform the same operation on the BUFFER_TYPE
. So ArrayView( ArrayView const & )
calls BUFFER_TYPE< T >( BUFFER_TYPE< T > const & )
and ArrayView::operator=( ArrayView && )
calls BUFFER_TYPE< T >::operator=( BUFFER_TYPE< T > && )
. With the exception of the StackBuffer
all of these operations are shallow copies.
TEST( Array, arrayView )
{
LvArray::Array< int,
2,
camp::idx_seq< 1, 0 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array( 5, 6 );
// Create a view.
LvArray::ArrayView< int,
2,
0,
std::ptrdiff_t,
LvArray::MallocBuffer > const view = array;
EXPECT_EQ( view.data(), array.data() );
// Create a view with const values.
LvArray::ArrayView< int const,
2,
0,
std::ptrdiff_t,
LvArray::MallocBuffer > const viewConst = array.toViewConst();
EXPECT_EQ( viewConst.data(), array.data() );
// Copy a view.
LvArray::ArrayView< int,
2,
0,
std::ptrdiff_t,
LvArray::MallocBuffer > const viewCopy = view;
EXPECT_EQ( viewCopy.data(), array.data() );
}
[Source: examples/exampleArray.cpp]
An LvArray::ArrayView
is almost always constructed from an existing LvArray::Array
but it does has a default constructor which constructs an uninitialized ArrayView
. The only valid use of an uninitialized LvArray::ArrayView
is to assign to it from an existing LvArray::ArrayView
. This behavior is primarily used when putting an LvArray::ArrayView
into a container.
LvArray::ArrayView
supports the subset of operations on LvArray::Array
that do not require reallocation or resizing. Every method in LvArray::ArrayView
is const
except the assignment operators. Most methods of are callable on device.
3.6. ArraySlice¶
Up until now all the examples using operator[]
have consumed all of the available indices, i.e. a two dimensional array has always had two immediate applications of operator[]
to access a value. But what is the result of a single application of operator[]
to a two dimensional array? Well it’s a one dimensional LvArray::ArraySlice
! LvArray::ArraySlice
shares the same template parameters as LvArray::ArrayView
except that it doesn’t need the BUFFER_TYPE
parameter. LvArray::ArrayView
is a feather weight object that consists only of three pointers: a pointer to the values, a pointer to the dimensions of the array and a pointer to the strides of the dimensions. LvArray::ArraySlice
has shallow copy semantics but no assignment operators and every method is const
.
Unlike LvArray::Array
and LvArray::ArrayView
the LvArray::ArraySlice
can refer to a non-contiguous set of values. As such it does not have a method data
but instead has dataIfContiguous
that performs a runtime check and aborts execution if the LvArray::ArraySlice
does not refer to a contiguous set of values. Every method is const
and callable on device.
TEST( Array, arraySlice )
{
{
LvArray::Array< int,
2,
camp::idx_seq< 0, 1 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array( 5, 6 );
// The unit stride dimension of array is 1 so when we slice off
// the first dimension the unit stride dimension of the slice is 0.
LvArray::ArraySlice< int,
1,
0,
std::ptrdiff_t > const slice = array[ 2 ];
EXPECT_TRUE( slice.isContiguous() );
EXPECT_EQ( slice.size(), 6 );
EXPECT_EQ( slice.size( 0 ), 6 );
slice[ 3 ] = 1;
}
{
LvArray::Array< int,
3,
camp::idx_seq< 2, 1, 0 >,
std::ptrdiff_t,
LvArray::MallocBuffer > array( 3, 5, 6 );
// The unit stride dimension of array is 0 so when we slice off
// the first dimension the unit stride dimension of the slice is -1.
LvArray::ArraySlice< int,
2,
-1,
std::ptrdiff_t > const slice = array[ 2 ];
EXPECT_FALSE( slice.isContiguous() );
EXPECT_EQ( slice.size(), 5 * 6 );
EXPECT_EQ( slice.size( 0 ), 5 );
EXPECT_EQ( slice.size( 1 ), 6 );
slice( 3, 4 ) = 1;
}
}
[Source: examples/exampleArray.cpp]
Note
A LvArray::ArraySlice
should not be captured in a device kernel, even if the slice comes from an array that has be moved to the device. This is because LvArray::ArraySlice
only contains a pointer to the dimension sizes and strides. Therefore when the slice is memcpy
’d to device the size and stride pointers will still point to host memory. This does not apply if you construct the slice manually and pass it device pointers but this is not a common use case.
3.7. Usage with LvArray::ChaiBuffer
¶
When using the LvArray::ChaiBuffer
as the buffer type the LvArray::Array
can exist in multiple memory spaces. It can be explicitly moved between spaces with the method move
and when the RAJA execution context is set the LvArray::ArrayView
copy constructor will ensure that the newly constructed view’s allocation is in the associated memory space.
CUDA_TEST( Array, chaiBuffer )
{
LvArray::Array< int,
2,
camp::idx_seq< 1, 0 >,
std::ptrdiff_t,
LvArray::ChaiBuffer > array( 5, 6 );
// Move the array to the device.
array.move( LvArray::MemorySpace::cuda );
int * const devicePointer = array.data();
RAJA::forall< RAJA::cuda_exec< 32 > >(
RAJA::TypedRangeSegment< std::ptrdiff_t >( 0, array.size() ),
[devicePointer] __device__ ( std::ptrdiff_t const i )
{
devicePointer[ i ] = i;
}
);
LvArray::ArrayView< int,
2,
0,
std::ptrdiff_t,
LvArray::ChaiBuffer > const & view = array;
// Capture the view in a host kernel which moves the data back to the host.
RAJA::forall< RAJA::loop_exec >(
RAJA::TypedRangeSegment< std::ptrdiff_t >( 0, view.size() ),
[view] ( std::ptrdiff_t const i )
{
EXPECT_EQ( view.data()[ i ], i );
}
);
}
[Source: examples/exampleArray.cpp]
Like the LvArray::ChaiBuffer
a new allocation is created the first time the array is moved to a new space. Every time the array is moved between spaces it is touched in the new space unless the value type T
is const
or the array was moved with the move
method and the optional touch
parameter was set to false
. The data is only copied between memory spaces if the last space the array was touched in is different from the space to move to. The name associated with the Array can be set with setName
.
TEST( Array, setName )
{
LvArray::Array< int,
2,
camp::idx_seq< 1, 0 >,
std::ptrdiff_t,
LvArray::ChaiBuffer > array( 1024, 1024 );
// Move the array to the device.
array.move( LvArray::MemorySpace::cuda );
// Provide a name and move the array to the host.
array.setName( "my_array" );
array.move( LvArray::MemorySpace::host );
}
[Source: examples/exampleArray.cpp]
Output
Moved 4.0 MB to the DEVICE: LvArray::Array<int, 2, camp::int_seq<long, 1l, 0l>, long, LvArray::ChaiBuffer>
Moved 4.0 MB to the HOST : LvArray::Array<int, 2, camp::int_seq<long, 1l, 0l>, long, LvArray::ChaiBuffer> my_array
3.8. Interacting with an LvArray::Array
of arbitrary dimension.¶
Often it can be useful to write a template function that can operate on an LvArray::Array
with an arbitrary number of dimension. If the order the data is accessed in is not important you can use data()
or the iterator interface begin
and end
. For example you could write a function to sum up the values in an array as follows
template< int NDIM, int USD >
int sum( LvArray::ArraySlice< int const, NDIM, USD, std::ptrdiff_t > const slice )
{
int value = 0;
for( int const val : slice )
{
value += val;
}
return value;
}
[Source: examples/exampleArray.cpp]
However this same approach might not work as intended for floating point arrays because the order of additions and hence the answer would change with the data layout. To calculate the sum consistently you need to iterate over the multidimensional index space of the array. For an array with a fixed number of dimensions you could just nest the appropriate number of for loops. LvArray::forValuesInSlice
defined in sliceHelpers.hpp
solves this problem, it iterates over the values in a LvArray::ArraySlice
in a consistent order regardless of the permutation. Using this sum
could be written as
template< int NDIM, int USD >
double sum( LvArray::ArraySlice< double const, NDIM, USD, std::ptrdiff_t > const slice )
{
double value = 0;
LvArray::forValuesInSlice( slice, [&value] ( double const val )
{
value += val;
} );
return value;
}
[Source: examples/exampleArray.cpp]
sliceHelpers.hpp
also provides a function forAllValuesInSliceWithIndices
which calls the user provided callback with the indicies to each value in addition to the value. This can be used to easily copy data between arrays with different layouts.
template< int NDIM, int DST_USD, int SRC_USD >
void copy( LvArray::ArraySlice< int, NDIM, DST_USD, std::ptrdiff_t > const dst,
LvArray::ArraySlice< int const, NDIM, SRC_USD, std::ptrdiff_t > const src )
{
for( int dim = 0; dim < NDIM; ++dim )
{
LVARRAY_ERROR_IF_NE( dst.size( dim ), src.size( dim ) );
}
LvArray::forValuesInSliceWithIndices( dst,
[src] ( int & val, auto const ... indices )
{
val = src( indices ... );
}
);
}
[Source: examples/exampleArray.cpp]
Finally you can write a recursive function that operates on an LvArray::ArraySlice
. An example of this is the stream output method for LvArray::Array
.
template< typename T, int NDIM, int USD, typename INDEX_TYPE >
std::ostream & operator<<( std::ostream & stream,
::LvArray::ArraySlice< T, NDIM, USD, INDEX_TYPE > const slice )
{
stream << "{ ";
if( slice.size( 0 ) > 0 )
{
stream << slice[ 0 ];
}
for( INDEX_TYPE i = 1; i < slice.size( 0 ); ++i )
{
stream << ", " << slice[ i ];
}
stream << " }";
return stream;
}
[Source: src/output.hpp]
3.9. Usage with LVARRAY_BOUNDS_CHECK
¶
When LVARRAY_BOUNDS_CHECK
is defined all accesses via operator[]
and operator()
is checked to make sure the indices are valid. If invalid indices are detected an error message is printed to standard out and the program is aborted. It should be noted that access via operator()
is able to provide a more useful error message upon an invalid access because it has access to all of the indices whereas operator[]
only has access to a single index at a time. size( int dim )
, linearIndex
, emplace
and insert
will also check that their arguments are in bounds.
TEST( Array, boundsCheck )
{
#if defined(ARRAY_USE_BOUNDS_CHECK)
LvArray::Array< int, 3, camp::idx_seq< 0, 1, 2 >, std::ptrdiff_t, LvArray::MallocBuffer > x( 3, 4, 5 );
// Out of bounds access aborts the program.
EXPECT_DEATH_IF_SUPPORTED( x( 2, 3, 4 ), "" );
EXPECT_DEATH_IF_SUPPORTED( x( -1, 4, 6 ), "" );
EXPECT_DEATH_IF_SUPPORTED( x[ 0 ][ 10 ][ 2 ], "" );
// Out of bounds emplace
LvArray::Array< int, 1, camp::idx_seq< 0 >, std::ptrdiff_t, LvArray::MallocBuffer > x( 10 );
EXPECT_DEATH_IF_SUPPORTED( x.emplace( -1, 5 ) );
#endif
}
[Source: examples/exampleArray.cpp]
3.10. Guidelines¶
In general you should opt to pass around the most restrictive array possible. If a function takes in an array and only needs to read and write the values then it should take in an LvArray::ArraySlice
. If it needs to read the values and captures the array in a kernel then it should take in an LvArray::ArrayView< T const, ... >
. Only when a function needs to resize or reallocate the array should it accept an LvArray::Array
.
Our benchmarks have shown no significant performance difference between using operator()
and nested applications of operator[]
, if you do see a difference let us know!