-
Notifications
You must be signed in to change notification settings - Fork 10
/
Copy pathexampleArrayOfArrays.cpp
397 lines (334 loc) · 13.8 KB
/
exampleArrayOfArrays.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
/*
* Copyright (c) 2020, Lawrence Livermore National Security, LLC and LvArray contributors.
* All rights reserved.
* See the LICENSE file for details.
* SPDX-License-Identifier: (BSD-3-Clause)
*/
// Source includes
#include "ArrayOfArrays.hpp"
#include "MallocBuffer.hpp"
#if defined(LVARRAY_USE_CHAI)
#include "ChaiBuffer.hpp"
#endif
#include "sortedArrayManipulation.hpp"
// TPL includes
#include <RAJA/RAJA.hpp>
#include <gtest/gtest.h>
// system includes
#include <string>
#include <random>
// You can't define a nvcc extended host-device or device lambda in a TEST statement.
#define CUDA_TEST( X, Y ) \
static void cuda_test_ ## X ## _ ## Y(); \
TEST( X, Y ) { cuda_test_ ## X ## _ ## Y(); } \
static void cuda_test_ ## X ## _ ## Y()
// Sphinx start after examples
TEST( ArrayOfArrays, construction )
{
// Create an empty ArrayOfArrays.
LvArray::ArrayOfArrays< std::string, std::ptrdiff_t, LvArray::MallocBuffer > arrayOfArrays;
EXPECT_EQ( arrayOfArrays.size(), 0 );
// Append an array of length 2.
arrayOfArrays.appendArray( 2 );
EXPECT_EQ( arrayOfArrays.size(), 1 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 2 );
EXPECT_EQ( arrayOfArrays.capacityOfArray( 0 ), 2 );
arrayOfArrays( 0, 0 ) = "First array, first entry.";
arrayOfArrays( 0, 1 ) = "First array, second entry.";
// Append another array of length 3.
arrayOfArrays.appendArray( 3 );
EXPECT_EQ( arrayOfArrays.size(), 2 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 1 ), 3 );
EXPECT_EQ( arrayOfArrays.capacityOfArray( 1 ), 3 );
arrayOfArrays( 1, 0 ) = "Second array, first entry.";
arrayOfArrays( 1, 2 ) = "Second array, third entry.";
EXPECT_EQ( arrayOfArrays[ 0 ][ 1 ], "First array, second entry." );
EXPECT_EQ( arrayOfArrays[ 1 ][ 2 ], "Second array, third entry." );
// Values are default initialized.
EXPECT_EQ( arrayOfArrays[ 1 ][ 1 ], std::string() );
}
TEST( ArrayOfArrays, modification )
{
LvArray::ArrayOfArrays< std::string, std::ptrdiff_t, LvArray::MallocBuffer > arrayOfArrays;
// Append an array.
arrayOfArrays.appendArray( 3 );
arrayOfArrays( 0, 0 ) = "First array, first entry.";
arrayOfArrays( 0, 1 ) = "First array, second entry.";
arrayOfArrays( 0, 2 ) = "First array, third entry.";
// Insert a new array at the beginning.
std::array< std::string, 2 > newFirstArray = { "New first array, first entry.",
"New first array, second entry." };
arrayOfArrays.insertArray( 0, newFirstArray.begin(), newFirstArray.end() );
EXPECT_EQ( arrayOfArrays.size(), 2 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 2 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 1 ), 3 );
EXPECT_EQ( arrayOfArrays( 0, 1 ), "New first array, second entry." );
EXPECT_EQ( arrayOfArrays[ 1 ][ 1 ], "First array, second entry." );
// Erase the values from what is now the second array.
arrayOfArrays.clearArray( 1 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 1 ), 0 );
// Append a value to the end of each array.
arrayOfArrays.emplaceBack( 1, "Second array, first entry." );
arrayOfArrays.emplaceBack( 0, "New first array, third entry." );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 3 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 1 ), 1 );
EXPECT_EQ( arrayOfArrays[ 1 ][ 0 ], "Second array, first entry." );
EXPECT_EQ( arrayOfArrays( 0, 2 ), "New first array, third entry." );
}
// Sphinx end before examples
#if defined(LVARRAY_USE_OPENMP)
// Sphinx start after view
TEST( ArrayOfArrays, view )
{
// Create an ArrayOfArrays with 10 inner arrays each with capacity 9.
LvArray::ArrayOfArrays< int, std::ptrdiff_t, LvArray::MallocBuffer > arrayOfArrays( 10, 9 );
{
// Then create a view.
LvArray::ArrayOfArraysView< int,
std::ptrdiff_t const,
false,
LvArray::MallocBuffer > const view = arrayOfArrays.toView();
EXPECT_EQ( view.size(), 10 );
EXPECT_EQ( view.capacityOfArray( 7 ), 9 );
// Modify every inner array in parallel
RAJA::forall< RAJA::omp_parallel_for_exec >(
RAJA::TypedRangeSegment< std::ptrdiff_t >( 0, view.size() ),
[view] ( std::ptrdiff_t const i )
{
for( std::ptrdiff_t j = 0; j < i; ++j )
{
view.emplaceBack( i, 10 * i + j );
}
}
);
EXPECT_EQ( view.sizeOfArray( 9 ), view.capacityOfArray( 9 ) );
// The last array is at capacity. Therefore any attempt at increasing its size through
// a view is invalid and may lead to undefined behavior!
// view.emplaceBack( 9, 0 );
}
{
// Create a view which cannot modify the sizes of the inner arrays.
LvArray::ArrayOfArraysView< int,
std::ptrdiff_t const,
true,
LvArray::MallocBuffer > const viewConstSizes = arrayOfArrays.toViewConstSizes();
for( std::ptrdiff_t i = 0; i < viewConstSizes.size(); ++i )
{
for( int & value : viewConstSizes[ i ] )
{
value *= 2;
}
// This would be a compilation error
// viewConstSizes.emplaceBack( i, 0 );
}
}
{
// Create a view which has read only access.
LvArray::ArrayOfArraysView< int const,
std::ptrdiff_t const,
true,
LvArray::MallocBuffer > const viewConst = arrayOfArrays.toViewConst();
for( std::ptrdiff_t i = 0; i < viewConst.size(); ++i )
{
for( std::ptrdiff_t j = 0; j < viewConst.sizeOfArray( i ); ++j )
{
EXPECT_EQ( viewConst( i, j ), 2 * ( 10 * i + j ) );
// This would be a compilation error
// viewConst( i, j ) = 0;
}
}
}
// Verify that all the modifications are present in the parent ArrayOfArrays.
EXPECT_EQ( arrayOfArrays.size(), 10 );
for( std::ptrdiff_t i = 0; i < arrayOfArrays.size(); ++i )
{
EXPECT_EQ( arrayOfArrays.sizeOfArray( i ), i );
for( std::ptrdiff_t j = 0; j < arrayOfArrays.sizeOfArray( i ); ++j )
{
EXPECT_EQ( arrayOfArrays( i, j ), 2 * ( 10 * i + j ) );
}
}
}
// Sphinx end before view
// Sphinx start after atomic
TEST( ArrayOfArrays, atomic )
{
// Create an ArrayOfArrays with 1 array with a capacity of 100.
LvArray::ArrayOfArrays< int, std::ptrdiff_t, LvArray::MallocBuffer > arrayOfArrays( 1, 100 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 0 );
// Then create a view.
LvArray::ArrayOfArraysView< int,
std::ptrdiff_t const,
false,
LvArray::MallocBuffer > const view = arrayOfArrays.toView();
// Append to the single inner array in parallel
RAJA::forall< RAJA::omp_parallel_for_exec >(
RAJA::TypedRangeSegment< std::ptrdiff_t >( 0, view.capacityOfArray( 0 ) ),
[view] ( std::ptrdiff_t const i )
{
view.emplaceBackAtomic< RAJA::builtin_atomic >( 0, i );
}
);
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 100 );
// Sort the entries in the array since they were appended in an arbitrary order
LvArray::sortedArrayManipulation::makeSorted( arrayOfArrays[ 0 ].begin(), arrayOfArrays[ 0 ].end() );
for( std::ptrdiff_t i = 0; i < arrayOfArrays.sizeOfArray( 0 ); ++i )
{
EXPECT_EQ( arrayOfArrays( 0, i ), i );
}
}
// Sphinx end before atomic
#endif
// Sphinx start after compress
TEST( ArrayOfArrays, compress )
{
// Create an ArrayOfArrays with three inner arrays each with capacity 5.
LvArray::ArrayOfArrays< int, std::ptrdiff_t, LvArray::MallocBuffer > arrayOfArrays( 3, 5 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 1 ), 0 );
EXPECT_EQ( arrayOfArrays.capacityOfArray( 1 ), 5 );
// Append values to the inner arrays.
std::array< int, 5 > values = { 0, 1, 2, 3, 4 };
arrayOfArrays.appendToArray( 0, values.begin(), values.begin() + 3 );
arrayOfArrays.appendToArray( 1, values.begin(), values.begin() + 4 );
arrayOfArrays.appendToArray( 2, values.begin(), values.begin() + 5 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 3 );
EXPECT_EQ( arrayOfArrays.capacityOfArray( 0 ), 5 );
// Since the first array has some extra capacity it is not contiguous in memory
// with the second array.
EXPECT_NE( &arrayOfArrays( 0, 2 ) + 1, &arrayOfArrays( 1, 0 ) );
// Compress the ArrayOfArrays. This doesn't change the sizes or values of the
// inner arrays, only their capacities.
arrayOfArrays.compress();
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 3 );
EXPECT_EQ( arrayOfArrays.capacityOfArray( 0 ), 3 );
// The values are preserved.
EXPECT_EQ( arrayOfArrays( 2, 4 ), 4 );
// And the inner arrays are now contiguous.
EXPECT_EQ( &arrayOfArrays( 0, 2 ) + 1, &arrayOfArrays( 1, 0 ) );
}
// Sphinx end before compress
// Sphinx start after resizeFromCapacities
TEST( ArrayOfArrays, resizeFromCapacities )
{
// Create an ArrayOfArrays with two inner arrays
LvArray::ArrayOfArrays< int, std::ptrdiff_t, LvArray::MallocBuffer > arrayOfArrays( 2 );
// Append values to the inner arrays.
std::array< int, 5 > values = { 0, 1, 2, 3, 4 };
arrayOfArrays.appendToArray( 0, values.begin(), values.begin() + 3 );
arrayOfArrays.appendToArray( 1, values.begin(), values.begin() + 4 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 3 );
// Resize the ArrayOfArrays from a new list of capacities.
std::array< std::ptrdiff_t, 3 > newCapacities = { 3, 5, 2 };
arrayOfArrays.resizeFromCapacities< RAJA::loop_exec >( newCapacities.size(), newCapacities.data() );
// This will clear any existing arrays.
EXPECT_EQ( arrayOfArrays.size(), 3 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 1 ), 0 );
EXPECT_EQ( arrayOfArrays.capacityOfArray( 1 ), newCapacities[ 1 ] );
}
// Sphinx end before resizeFromCapacities
#if defined(LVARRAY_USE_CUDA) && defined(LVARRAY_USE_CHAI)
// Sphinx start after ChaiBuffer
CUDA_TEST( ArrayOfArrays, ChaiBuffer )
{
LvArray::ArrayOfArrays< int, std::ptrdiff_t, LvArray::ChaiBuffer > arrayOfArrays( 10, 9 );
{
// Create a view.
LvArray::ArrayOfArraysView< int,
std::ptrdiff_t const,
false,
LvArray::ChaiBuffer > const view = arrayOfArrays.toView();
// Capture the view on device. This will copy the values, sizes and offsets.
// The values and sizes will be touched.
RAJA::forall< RAJA::cuda_exec< 32 > >(
RAJA::TypedRangeSegment< std::ptrdiff_t >( 0, view.size() ),
[view] __device__ ( std::ptrdiff_t const i )
{
for( std::ptrdiff_t j = 0; j < i; ++j )
{
view.emplace( i, 0, 10 * i + j );
}
}
);
}
{
// Create a view which cannot modify the sizes of the inner arrays.
LvArray::ArrayOfArraysView< int,
std::ptrdiff_t const,
true,
LvArray::ChaiBuffer > const viewConstSizes = arrayOfArrays.toViewConstSizes();
// Capture the view on the host. This will copy back the values and sizes since they were previously touched
// on device. It will only touch the values on host.
RAJA::forall< RAJA::loop_exec >(
RAJA::TypedRangeSegment< std::ptrdiff_t >( 0, viewConstSizes.size() ),
[viewConstSizes] ( std::ptrdiff_t const i )
{
for( int & value : viewConstSizes[ i ] )
{
value *= 2;
}
}
);
}
{
// Create a view which has read only access.
LvArray::ArrayOfArraysView< int const,
std::ptrdiff_t const,
true,
LvArray::ChaiBuffer > const viewConst = arrayOfArrays.toViewConst();
// Capture the view on device. Since the values were previously touched on host it will copy them over.
// Both the sizes and offsets are current on device so they are not copied over. Nothing is touched.
RAJA::forall< RAJA::loop_exec >(
RAJA::TypedRangeSegment< std::ptrdiff_t >( 0, viewConst.size() ),
[viewConst] ( std::ptrdiff_t const i )
{
for( std::ptrdiff_t j = 0; j < viewConst.sizeOfArray( i ); ++j )
{
LVARRAY_ERROR_IF_NE( viewConst( i, j ), 2 * ( 10 * i + i - j - 1 ) );
}
}
);
}
// This won't copy any data since everything is current on host. It will however touch the values,
// sizes and offsets.
arrayOfArrays.move( LvArray::MemorySpace::CPU );
// Verify that all the modifications are present in the parent ArrayOfArrays.
EXPECT_EQ( arrayOfArrays.size(), 10 );
for( std::ptrdiff_t i = 0; i < arrayOfArrays.size(); ++i )
{
EXPECT_EQ( arrayOfArrays.sizeOfArray( i ), i );
for( std::ptrdiff_t j = 0; j < arrayOfArrays.sizeOfArray( i ); ++j )
{
EXPECT_EQ( arrayOfArrays( i, j ), 2 * ( 10 * i + i - j - 1 ) );
}
}
}
// Sphinx end before ChaiBuffer
#endif
// Sphinx start after bounds check
TEST( ArrayOfArrays, boundsCheck )
{
#if defined(LVARRAY_BOUNDS_CHECK)
LvArray::ArrayOfArrays< int, std::ptrdiff_t, LvArray::MallocBuffer > arrayOfArrays;
// Append an array.
std::array< int, 5 > values = { 0, 1, 2, 3, 4 };
arrayOfArrays.appendArray( values.begin(), values.end() );
EXPECT_EQ( arrayOfArrays.size(), 1 );
EXPECT_EQ( arrayOfArrays.sizeOfArray( 0 ), 5 );
// Out of bounds access aborts the program.
EXPECT_DEATH_IF_SUPPORTED( arrayOfArrays( 0, -1 ), "" );
EXPECT_DEATH_IF_SUPPORTED( arrayOfArrays[ 0 ][ 6 ], "" );
EXPECT_DEATH_IF_SUPPORTED( arrayOfArrays[ 1 ][ 5 ], "" );
EXPECT_DEATH_IF_SUPPORTED( arrayOfArrays.capacityOfArray( 5 ), "" );
EXPECT_DEATH_IF_SUPPORTED( arrayOfArrays.insertArray( 5, values.begin(), values.end() ), "" );
EXPECT_DEATH_IF_SUPPORTED( arrayOfArrays.emplace( 0, 44, 4 ), "" );
EXPECT_DEATH_IF_SUPPORTED( arrayOfArrays.emplace( 1, 44, 4 ), "" );
#endif
}
// Sphinx end before bounds check
// This is the default gtest main method. It is included for ease of debugging.
int main( int argc, char * * argv )
{
::testing::InitGoogleTest( &argc, argv );
int const result = RUN_ALL_TESTS();
return result;
}