2 // ************************************************************************
3 // This file is property of and copyright by the ALICE HLT Project *
4 // ALICE Experiment at CERN, All rights reserved. *
5 // See cxx source for full Copyright notice *
7 //*************************************************************************
9 #ifndef ALIHLTTPCCAPROCESS_H
10 #define ALIHLTTPCCAPROCESS_H
14 * Definitions needed for AliHLTTPCCATracker
18 #include "AliHLTTPCCADef.h"
19 #include "AliHLTTPCCATrackParam.h"
21 class AliHLTTPCCATracker;
23 #if defined(HLTCA_GPUCODE)
25 template<class TProcess>
26 GPUg() void AliHLTTPCCAProcess(int iSlice)
28 AliHLTTPCCATracker &tracker = ( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker )[iSlice];
29 GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
31 for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
33 TProcess::Thread( gridDim.x, blockDim.x, blockIdx.x, threadIdx.x, iSync, smem, tracker );
37 template<class TProcess>
38 GPUg() void AliHLTTPCCAProcessMulti(int firstSlice, int nSliceCount)
40 const int iSlice = nSliceCount * (blockIdx.x + (gridDim.x % nSliceCount != 0 && nSliceCount * (blockIdx.x + 1) % gridDim.x != 0)) / gridDim.x;
41 const int nSliceBlockOffset = gridDim.x * iSlice / nSliceCount;
42 const int sliceBlockId = blockIdx.x - nSliceBlockOffset;
43 const int sliceGridDim = gridDim.x * (iSlice + 1) / nSliceCount - gridDim.x * (iSlice) / nSliceCount;
44 AliHLTTPCCATracker &tracker = ( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker )[firstSlice + iSlice];
45 GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
47 for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
49 TProcess::Thread( sliceGridDim, blockDim.x, sliceBlockId, threadIdx.x, iSync, smem, tracker );
55 template<class TProcess>
56 GPUg() void AliHLTTPCCAProcess( int nBlocks, int nThreads, AliHLTTPCCATracker &tracker )
58 for ( int iB = 0; iB < nBlocks; iB++ ) {
59 typename TProcess::AliHLTTPCCASharedMemory smem;
60 for ( int iS = 0; iS <= TProcess::NThreadSyncPoints(); iS++ )
61 for ( int iT = 0; iT < nThreads; iT++ ) {
62 TProcess::Thread( nBlocks, nThreads, iB, iT, iS, smem, tracker );
67 #endif //HLTCA_GPUCODE
71 #if defined(HLTCA_GPUCODE)
73 template<typename TProcess>
74 GPUg() void AliHLTTPCCAProcess1()
76 AliHLTTPCCATracker &tracker = *( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker );
77 AliHLTTPCCATrackParam tParam;
79 GPUshared() typename TProcess::AliHLTTPCCASharedMemory sMem;
81 typename TProcess::AliHLTTPCCAThreadMemory rMem;
83 for ( int iSync = 0; iSync <= TProcess::NThreadSyncPoints(); iSync++ ) {
85 TProcess::Thread( gridDim.x, blockDim.x, blockIdx.x, threadIdx.x, iSync,
86 sMem, rMem, tracker, tParam );
92 template<typename TProcess>
93 GPUg() void AliHLTTPCCAProcess1( int nBlocks, int nThreads, AliHLTTPCCATracker &tracker )
95 for ( int iB = 0; iB < nBlocks; iB++ ) {
96 typename TProcess::AliHLTTPCCASharedMemory smem;
97 typename TProcess::AliHLTTPCCAThreadMemory *rMem = new typename TProcess::AliHLTTPCCAThreadMemory[nThreads];
98 AliHLTTPCCATrackParam *tParam = new AliHLTTPCCATrackParam[ nThreads ];
99 for ( int iS = 0; iS <= TProcess::NThreadSyncPoints(); iS++ ) {
100 for ( int iT = 0; iT < nThreads; iT++ )
101 TProcess::Thread( nBlocks, nThreads, iB, iT, iS, smem, rMem[iT], tracker, tParam[iT] );
108 #endif //HLTCA_GPUCODE
110 #endif //ALIHLTTPCCAPROCESS_H