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 MEM_CLASS_PRE() class AliHLTTPCCATracker;
23 #if defined(HLTCA_GPUCODE)
29 template<class TProcess>
30 GPUg() void AliHLTTPCCAProcess(int iSlice)
32 AliHLTTPCCATracker &tracker = ( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker )[iSlice];
33 GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
35 for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
37 TProcess::Thread( get_num_groups(0), get_local_size(0), get_group_id(0), get_local_id(0), iSync, smem, tracker );
41 template <class TProcess>
42 GPUg() void AliHLTTPCCAProcessMultiA(int firstSlice, int nSliceCount, int nVirtualBlocks)
44 if (get_group_id(0) >= nSliceCount) return;
45 AliHLTTPCCATracker &tracker = ( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker )[firstSlice + get_group_id(0)];
47 GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
49 for (int i = 0;i < nVirtualBlocks;i++)
51 for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
53 TProcess::Thread( nVirtualBlocks, get_local_size(0), i, get_local_id(0), iSync, smem, tracker );
58 template<class TProcess>
59 GPUg() void AliHLTTPCCAProcessMulti(int firstSlice, int nSliceCount)
61 const int iSlice = nSliceCount * (get_group_id(0) + (get_num_groups(0) % nSliceCount != 0 && nSliceCount * (get_group_id(0) + 1) % get_num_groups(0) != 0)) / get_num_groups(0);
62 const int nSliceBlockOffset = get_num_groups(0) * iSlice / nSliceCount;
63 const int sliceBlockId = get_group_id(0) - nSliceBlockOffset;
64 const int sliceGridDim = get_num_groups(0) * (iSlice + 1) / nSliceCount - get_num_groups(0) * (iSlice) / nSliceCount;
65 AliHLTTPCCATracker &tracker = ( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker )[firstSlice + iSlice];
66 GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
68 for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
70 TProcess::Thread( sliceGridDim, get_local_size(0), sliceBlockId, get_local_id(0), iSync, smem, tracker );
74 template<typename TProcess>
75 GPUg() void AliHLTTPCCAProcess1()
77 AliHLTTPCCATracker &tracker = *( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker );
78 AliHLTTPCCATrackParam tParam;
80 GPUshared() typename TProcess::AliHLTTPCCASharedMemory sMem;
82 typename TProcess::AliHLTTPCCAThreadMemory rMem;
84 for ( int iSync = 0; iSync <= TProcess::NThreadSyncPoints(); iSync++ ) {
86 TProcess::Thread( get_num_groups(0), get_local_size(0), get_group_id(0), get_local_id(0), iSync,
87 sMem, rMem, tracker, tParam );
95 template<class TProcess>
96 GPUg() void AliHLTTPCCAProcess( int nBlocks, int nThreads, AliHLTTPCCATracker &tracker )
98 for ( int iB = 0; iB < nBlocks; iB++ ) {
99 typename TProcess::AliHLTTPCCASharedMemory smem;
100 for ( int iS = 0; iS <= TProcess::NThreadSyncPoints(); iS++ )
101 for ( int iT = 0; iT < nThreads; iT++ ) {
102 TProcess::Thread( nBlocks, nThreads, iB, iT, iS, smem, tracker );
107 template<typename TProcess>
108 GPUg() void AliHLTTPCCAProcess1( int nBlocks, int nThreads, AliHLTTPCCATracker &tracker )
110 for ( int iB = 0; iB < nBlocks; iB++ ) {
111 typename TProcess::AliHLTTPCCASharedMemory smem;
112 typename TProcess::AliHLTTPCCAThreadMemory *rMem = new typename TProcess::AliHLTTPCCAThreadMemory[nThreads];
113 AliHLTTPCCATrackParam *tParam = new AliHLTTPCCATrackParam[ nThreads ];
114 for ( int iS = 0; iS <= TProcess::NThreadSyncPoints(); iS++ ) {
115 for ( int iT = 0; iT < nThreads; iT++ )
116 TProcess::Thread( nBlocks, nThreads, iB, iT, iS, smem, rMem[iT], tracker, tParam[iT] );
123 #endif //HLTCA_GPUCODE
125 #endif //ALIHLTTPCCAPROCESS_H