]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/tracking-ca/AliHLTTPCCAProcess.h
bug fix: reconstruction crash when the output buffer size exceed
[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 class AliHLTTPCCATracker;
22
23 #if defined(HLTCA_GPUCODE)
24
25 template<class TProcess>
26 GPUg() void AliHLTTPCCAProcess(int iSlice)
27 {
28   AliHLTTPCCATracker &tracker = ( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker )[iSlice];
29   GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
30
31   for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
32     __syncthreads();
33     TProcess::Thread( gridDim.x, blockDim.x, blockIdx.x, threadIdx.x, iSync, smem, tracker  );
34   }
35 }
36
37 template<class TProcess>
38 GPUg() void AliHLTTPCCAProcessMulti(int firstSlice, int nSliceCount)
39 {
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;
46
47   for( int iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
48     __syncthreads();
49     TProcess::Thread( sliceGridDim, blockDim.x, sliceBlockId, threadIdx.x, iSync, smem, tracker  );
50   }
51 }
52
53 #else
54
55 template<class TProcess>
56 GPUg() void AliHLTTPCCAProcess( int nBlocks, int nThreads, AliHLTTPCCATracker &tracker )
57 {
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  );
63       }
64   }
65 }
66
67 #endif //HLTCA_GPUCODE
68
69
70
71 #if defined(HLTCA_GPUCODE)
72
73 template<typename TProcess>
74 GPUg() void AliHLTTPCCAProcess1()
75 {
76   AliHLTTPCCATracker &tracker = *( ( AliHLTTPCCATracker* ) gAliHLTTPCCATracker );
77   AliHLTTPCCATrackParam tParam;
78
79   GPUshared() typename TProcess::AliHLTTPCCASharedMemory sMem;
80
81   typename TProcess::AliHLTTPCCAThreadMemory rMem;
82
83   for ( int iSync = 0; iSync <= TProcess::NThreadSyncPoints(); iSync++ ) {
84     GPUsync();
85     TProcess::Thread( gridDim.x, blockDim.x, blockIdx.x, threadIdx.x, iSync,
86                       sMem, rMem, tracker, tParam  );
87   }
88 }
89
90 #else
91
92 template<typename TProcess>
93 GPUg() void AliHLTTPCCAProcess1( int nBlocks, int nThreads, AliHLTTPCCATracker &tracker )
94 {
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]  );
102     }
103     delete[] rMem;
104     delete[] tParam;
105   }
106 }
107
108 #endif //HLTCA_GPUCODE
109
110 #endif //ALIHLTTPCCAPROCESS_H