]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/tracking-ca/AliHLTTPCCAProcess.h
Update of the HLT CA tracker
[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 #ifndef ALIHLTTPCCAPROCESS_H
8 #define ALIHLTTPCCAPROCESS_H
9
10
11 /**
12  * Definitions needed for AliHLTTPCCATracker
13  *
14  */
15
16 #include "AliHLTTPCCADef.h"
17 #include "AliHLTTPCCATrackParam.h"
18
19 class AliHLTTPCCATracker;
20
21 #if defined(HLTCA_GPUCODE)
22
23 template<class TProcess>
24 GPUg() void AliHLTTPCCAProcess()
25 {
26   AliHLTTPCCATracker &tracker = *((AliHLTTPCCATracker*) cTracker);
27
28   GPUshared() typename TProcess::AliHLTTPCCASharedMemory smem;
29
30   TProcess::Thread( gridDim.x, blockDim.x, blockIdx.x, threadIdx.x, 0, smem, tracker  );
31
32 #define GPUPROCESS(iSync) \
33   if( TProcess::NThreadSyncPoints()>=iSync ){ \
34     GPUsync(); \
35     TProcess::Thread( gridDim.x, blockDim.x, blockIdx.x, threadIdx.x, iSync, smem, tracker  ); \
36   }
37
38   GPUPROCESS(1)
39   GPUPROCESS(2)
40   GPUPROCESS(3)
41     
42     //for( Int_t iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
43     //__syncthreads();
44     //TProcess::ThreadGPU( gridDim.x, blockDim.x, blockIdx.x, threadIdx.x, iSync, smem, tracker  ); 
45     //}
46     
47 #undef GPUPROCESS
48 }
49
50 #else
51
52 template<class TProcess>
53 GPUg() void AliHLTTPCCAProcess( Int_t nBlocks, Int_t nThreads, AliHLTTPCCATracker &tracker )
54 {  
55   for( Int_t iB=0; iB<nBlocks; iB++ ){
56     typename TProcess::AliHLTTPCCASharedMemory smem;
57     for( Int_t iS=0; iS<=TProcess::NThreadSyncPoints(); iS++)
58       for( Int_t iT=0; iT<nThreads; iT++ ){     
59         TProcess::Thread( nBlocks, nThreads, iB, iT, iS, smem, tracker  );
60       }
61   }
62 }
63
64 #endif
65
66
67
68 #if defined(HLTCA_GPUCODE)
69
70 template<typename TProcess>
71 GPUg() void AliHLTTPCCAProcess1()
72 {
73   AliHLTTPCCATracker &tracker = *((AliHLTTPCCATracker*) cTracker);
74   AliHLTTPCCATrackParam tParam; 
75   
76   GPUshared() typename TProcess::AliHLTTPCCASharedMemory sMem;  
77   
78   typename TProcess::AliHLTTPCCAThreadMemory rMem;
79
80   for( Int_t iSync=0; iSync<=TProcess::NThreadSyncPoints(); iSync++){
81     GPUsync(); 
82     TProcess::Thread( gridDim.x, blockDim.x, blockIdx.x, threadIdx.x, iSync, 
83                       sMem, rMem, tracker, tParam  ); 
84   }  
85 }
86
87 #else
88
89 template<typename TProcess>
90 GPUg() void AliHLTTPCCAProcess1( Int_t nBlocks, Int_t nThreads, AliHLTTPCCATracker &tracker )
91 {
92   for( Int_t iB=0; iB<nBlocks; iB++ ){
93     typename TProcess::AliHLTTPCCASharedMemory smem;
94     typename TProcess::AliHLTTPCCAThreadMemory rMem[nThreads];
95     AliHLTTPCCATrackParam tParam[nThreads];
96     for( Int_t iS=0; iS<=TProcess::NThreadSyncPoints(); iS++){
97       for( Int_t iT=0; iT<nThreads; iT++ )
98         TProcess::Thread( nBlocks, nThreads, iB, iT, iS, smem, rMem[iT], tracker, tParam[iT]  );
99     }
100   }
101 }
102
103 #endif
104
105 #endif