]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/tracking-ca/AliHLTTPCCAProcess.h
OpenCL version of the HLT tracker added (the new code is not used in standard compila...
[u/mrichter/AliRoot.git] / HLT / TPCLib / tracking-ca / AliHLTTPCCAProcess.h
1 //-*- Mode: C++ -*-
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                               *
6 //                                                                        *
7 //*************************************************************************
8
9 #ifndef ALIHLTTPCCAPROCESS_H
10 #define ALIHLTTPCCAPROCESS_H
11
12
13 /**
14  * Definitions needed for AliHLTTPCCATracker
15  *
16  */
17
18 #include "AliHLTTPCCADef.h"
19 #include "AliHLTTPCCATrackParam.h"
20
21 MEM_CLASS_PRE() class AliHLTTPCCATracker;
22
23 #if defined(HLTCA_GPUCODE)
24
25 #ifdef __OPENCL__
26
27 #else //__OPENCL__
28
29 template<class TProcess>
30 GPUg() void AliHLTTPCCAProcess(int iSlice)
31 {
32   AliHLTTPCCATracker &tracker = ( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker )[iSlice];
33   GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
34
35   for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
36     GPUsync();
37     TProcess::Thread( get_num_groups(0), get_local_size(0), get_group_id(0), get_local_id(0), iSync, smem, tracker  );
38   }
39 }
40
41 template <class TProcess>
42 GPUg() void AliHLTTPCCAProcessMultiA(int firstSlice, int nSliceCount, int nVirtualBlocks)
43 {
44         if (get_group_id(0) >= nSliceCount) return;
45         AliHLTTPCCATracker &tracker = ( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker )[firstSlice + get_group_id(0)];
46
47         GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
48
49         for (int i = 0;i < nVirtualBlocks;i++)
50         {
51                 for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
52                         GPUsync();
53                         TProcess::Thread( nVirtualBlocks, get_local_size(0), i, get_local_id(0), iSync, smem, tracker  );
54                 }               
55         }
56 }
57
58 template<class TProcess>
59 GPUg() void AliHLTTPCCAProcessMulti(int firstSlice, int nSliceCount)
60 {
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;
67
68   for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
69     GPUsync();
70     TProcess::Thread( sliceGridDim, get_local_size(0), sliceBlockId, get_local_id(0), iSync, smem, tracker  );
71   }
72 }
73
74 template<typename TProcess>
75 GPUg() void AliHLTTPCCAProcess1()
76 {
77   AliHLTTPCCATracker &tracker = *( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker );
78   AliHLTTPCCATrackParam tParam;
79
80   GPUshared() typename TProcess::AliHLTTPCCASharedMemory sMem;
81
82   typename TProcess::AliHLTTPCCAThreadMemory rMem;
83
84   for ( int iSync = 0; iSync <= TProcess::NThreadSyncPoints(); iSync++ ) {
85     GPUsync();
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  );
88   }
89 }
90
91 #endif //__OPENCL__
92
93 #else //HLTCA_GPUCODE
94
95 template<class TProcess>
96 GPUg() void AliHLTTPCCAProcess( int nBlocks, int nThreads, AliHLTTPCCATracker &tracker )
97 {
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  );
103       }
104   }
105 }
106
107 template<typename TProcess>
108 GPUg() void AliHLTTPCCAProcess1( int nBlocks, int nThreads, AliHLTTPCCATracker &tracker )
109 {
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]  );
117     }
118     delete[] rMem;
119     delete[] tParam;
120   }
121 }
122
123 #endif //HLTCA_GPUCODE
124
125 #endif //ALIHLTTPCCAPROCESS_H