ATLAS Offline Software
Loading...
Searching...
No Matches
LinearTransformTaskExampleAlg.cxx
Go to the documentation of this file.
1//
2// Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3//
4
5// Local include(s).
7#include "CudaMultiplyTask.h"
8
9// AthCUDA include(s).
11
12namespace AthCUDAExamples {
13
15
16 // Retrieve the necessary component(s).
17 ATH_CHECK( m_kernelRunnerSvc.retrieve() );
18
19 // Return gracefully.
20 return StatusCode::SUCCESS;
21 }
22
23 StatusCode
24 LinearTransformTaskExampleAlg::execute( const EventContext& ) const {
25
26 // Create a dummy array variable that will be multiplied by some amount.
27 static const std::size_t ARRAY_SIZE = 10000;
28 std::vector< float > dummyArray;
29 dummyArray.reserve( ARRAY_SIZE );
30 static const float ARRAY_ELEMENT = 3.141592f;
31 for( std::size_t i = 0; i < ARRAY_SIZE; ++i ) {
32 dummyArray.push_back( ARRAY_ELEMENT );
33 }
34
35 // Create and launch the calculation using a task object.
37 static const float MULTIPLIER = 1.23f;
38 auto task = make_CudaMultiplyTask( status, dummyArray, MULTIPLIER );
39 ATH_CHECK( m_kernelRunnerSvc->execute( std::move( task ) ) );
40
41 // Now wait for the task to finish.
42 if( status.wait() != 0 ) {
43 ATH_MSG_ERROR( "Something went wrong in the execution of the CUDA "
44 "task" );
45 return StatusCode::FAILURE;
46 }
47
48 // Check if the operation succeeded.
49 static const float EXPECTED_RESULT = ARRAY_ELEMENT * MULTIPLIER;
50 for( std::size_t i = 0; i < ARRAY_SIZE; ++i ) {
51 if( std::abs( dummyArray[ i ] - EXPECTED_RESULT ) > 0.001 ) {
52 ATH_MSG_ERROR( "The CUDA transformation failed to run" );
53 return StatusCode::FAILURE;
54 }
55 }
56
57 // Return gracefully.
58 return StatusCode::SUCCESS;
59 }
60
61} // namespace AthCUDAExamples
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
virtual StatusCode execute(const EventContext &ctx) const override
Function executing the algorithm.
virtual StatusCode initialize() override
Function initialising the algorithm.
ServiceHandle< AthCUDA::IKernelRunnerSvc > m_kernelRunnerSvc
The kernel runner service to use.
Helper object used for synchronising GPU kernel tasks.
std::unique_ptr< AthCUDA::IKernelTask > make_CudaMultiplyTask(AthCUDA::KernelStatus &status, std::vector< float > &array, float multiplier)
Function setting up the task performing a simple linear transformation.