1 // **************************************************************************
2 // This file is property of and copyright by the ALICE HLT Project *
3 // ALICE Experiment at CERN, All rights reserved. *
5 // Primary Authors: Sergey Gorbunov <sergey.gorbunov@kip.uni-heidelberg.de> *
6 // Ivan Kisel <kisel@kip.uni-heidelberg.de> *
7 // David Rohr <drohr@kip.uni-heidelberg.de> *
8 // for The ALICE HLT Project. *
10 // Permission to use, copy, modify and distribute this software and its *
11 // documentation strictly for non-commercial purposes is hereby granted *
12 // without fee, provided that the above copyright notice appears in all *
13 // copies and that both the copyright notice and this permission notice *
14 // appear in the supporting documentation. The authors make no claims *
15 // about the suitability of this software for any purpose. It is *
16 // provided "as is" without express or implied warranty. *
18 //***************************************************************************
20 #include "AliHLTTPCCAGPUTracker.h"
28 #include <sys/syscall.h>
31 #include "AliHLTTPCCADef.h"
32 #include "AliHLTTPCCAGPUConfig.h"
34 #include <sm_11_atomic_functions.h>
35 #include <sm_12_atomic_functions.h>
39 //Disable assertions since they produce errors in GPU Code
45 __constant__ float4 gAliHLTTPCCATracker[HLTCA_GPU_TRACKER_CONSTANT_MEM / sizeof( float4 )];
46 #ifdef HLTCA_GPU_TEXTURE_FETCH
47 texture<ushort2, 1, cudaReadModeElementType> gAliTexRefu2;
48 texture<unsigned short, 1, cudaReadModeElementType> gAliTexRefu;
49 texture<signed short, 1, cudaReadModeElementType> gAliTexRefs;
52 //Include CXX Files, GPUd() macro will then produce CUDA device code out of the tracker source code
53 #include "AliHLTTPCCATrackParam.cxx"
54 #include "AliHLTTPCCATrack.cxx"
56 #include "AliHLTTPCCAHitArea.cxx"
57 #include "AliHLTTPCCAGrid.cxx"
58 #include "AliHLTTPCCARow.cxx"
59 #include "AliHLTTPCCAParam.cxx"
60 #include "AliHLTTPCCATracker.cxx"
62 #include "AliHLTTPCCAProcess.h"
64 #include "AliHLTTPCCATrackletSelector.cxx"
65 #include "AliHLTTPCCANeighboursFinder.cxx"
66 #include "AliHLTTPCCANeighboursCleaner.cxx"
67 #include "AliHLTTPCCAStartHitsFinder.cxx"
68 #include "AliHLTTPCCAStartHitsSorter.cxx"
69 #include "AliHLTTPCCATrackletConstructor.cxx"
70 #include "AliHLTTPCCASliceOutput.cxx"
72 #include "MemoryAssignmentHelpers.h"
74 #ifndef HLTCA_STANDALONE
75 #include "AliHLTDefinitions.h"
76 #include "AliHLTSystem.h"
79 ClassImp( AliHLTTPCCAGPUTracker )
81 bool AliHLTTPCCAGPUTracker::fgGPUUsed = false;
83 int AliHLTTPCCAGPUTracker::InitGPU(int sliceCount, int forceDeviceID)
85 //Find best CUDA device, initialize and allocate memory
89 HLTWarning("CUDA already used by another AliHLTTPCCAGPUTracker running in same process");
93 fThreadId = GetThread();
95 cudaDeviceProp fCudaDeviceProp;
97 fGPUMemSize = HLTCA_GPU_ROWS_MEMORY + HLTCA_GPU_COMMON_MEMORY + sliceCount * (HLTCA_GPU_SLICE_DATA_MEMORY + HLTCA_GPU_GLOBAL_MEMORY);
99 #ifndef CUDA_DEVICE_EMULATION
100 int count, bestDevice = -1;
101 long long int bestDeviceSpeed = 0, deviceSpeed;
102 if (CudaFailedMsg(cudaGetDeviceCount(&count)))
104 HLTError("Error getting CUDA Device Count");
108 if (fDebugLevel >= 2) HLTInfo("Available CUDA devices:");
109 for (int i = 0;i < count;i++)
111 unsigned int free, total;
114 cuDeviceGet(&tmpDevice, i);
115 CUcontext tmpContext;
116 cuCtxCreate(&tmpContext, 0, tmpDevice);
117 if(cuMemGetInfo(&free, &total)) std::cout << "Error\n";
118 cuCtxDestroy(tmpContext);
119 CudaFailedMsg(cudaGetDeviceProperties(&fCudaDeviceProp, i));
121 int deviceOK = fCudaDeviceProp.major < 9 && !(fCudaDeviceProp.major < 1 || (fCudaDeviceProp.major == 1 && fCudaDeviceProp.minor < 2)) && free >= fGPUMemSize;
123 if (fDebugLevel >= 2) HLTInfo("%s%2d: %s (Rev: %d.%d - Mem Avail %d / %d)%s", deviceOK ? " " : "[", i, fCudaDeviceProp.name, fCudaDeviceProp.major, fCudaDeviceProp.minor, free, fCudaDeviceProp.totalGlobalMem, deviceOK ? "" : " ]");
124 deviceSpeed = (long long int) fCudaDeviceProp.multiProcessorCount * (long long int) fCudaDeviceProp.clockRate * (long long int) fCudaDeviceProp.warpSize * (long long int) free;
125 if (deviceOK && deviceSpeed > bestDeviceSpeed)
128 bestDeviceSpeed = deviceSpeed;
131 if (bestDevice == -1)
133 HLTWarning("No CUDA Device available, aborting CUDA Initialisation");
139 if (forceDeviceID == -1)
140 cudaDevice = bestDevice;
142 cudaDevice = forceDeviceID;
147 cudaGetDeviceProperties(&fCudaDeviceProp ,cudaDevice );
149 if (fDebugLevel >= 1)
151 HLTInfo("Using CUDA Device %s with Properties:", fCudaDeviceProp.name);
152 HLTInfo("totalGlobalMem = %d", fCudaDeviceProp.totalGlobalMem);
153 HLTInfo("sharedMemPerBlock = %d", fCudaDeviceProp.sharedMemPerBlock);
154 HLTInfo("regsPerBlock = %d", fCudaDeviceProp.regsPerBlock);
155 HLTInfo("warpSize = %d", fCudaDeviceProp.warpSize);
156 HLTInfo("memPitch = %d", fCudaDeviceProp.memPitch);
157 HLTInfo("maxThreadsPerBlock = %d", fCudaDeviceProp.maxThreadsPerBlock);
158 HLTInfo("maxThreadsDim = %d %d %d", fCudaDeviceProp.maxThreadsDim[0], fCudaDeviceProp.maxThreadsDim[1], fCudaDeviceProp.maxThreadsDim[2]);
159 HLTInfo("maxGridSize = %d %d %d", fCudaDeviceProp.maxGridSize[0], fCudaDeviceProp.maxGridSize[1], fCudaDeviceProp.maxGridSize[2]);
160 HLTInfo("totalConstMem = %d", fCudaDeviceProp.totalConstMem);
161 HLTInfo("major = %d", fCudaDeviceProp.major);
162 HLTInfo("minor = %d", fCudaDeviceProp.minor);
163 HLTInfo("clockRate %d= ", fCudaDeviceProp.clockRate);
164 HLTInfo("textureAlignment %d= ", fCudaDeviceProp.textureAlignment);
167 if (fCudaDeviceProp.major < 1 || (fCudaDeviceProp.major == 1 && fCudaDeviceProp.minor < 2))
169 HLTError( "Unsupported CUDA Device" );
174 if (sizeof(AliHLTTPCCATracker) * sliceCount > HLTCA_GPU_TRACKER_OBJECT_MEMORY)
176 HLTError("Insufficiant Tracker Object Memory");
181 if (CudaFailedMsg(cudaSetDevice(cudaDevice)))
183 HLTError("Could not set CUDA Device!");
188 if (fgkNSlices * AliHLTTPCCATracker::CommonMemorySize() > HLTCA_GPU_COMMON_MEMORY)
190 HLTError("Insufficiant Common Memory");
196 if (fgkNSlices * (HLTCA_ROW_COUNT + 1) * sizeof(AliHLTTPCCARow) > HLTCA_GPU_ROWS_MEMORY)
198 HLTError("Insufficiant Row Memory");
204 if (fGPUMemSize > fCudaDeviceProp.totalGlobalMem || CudaFailedMsg(cudaMalloc(&fGPUMemory, (size_t) fGPUMemSize)))
206 HLTError("CUDA Memory Allocation Error");
211 if (fDebugLevel >= 1) HLTInfo("GPU Memory used: %d", (int) fGPUMemSize);
212 int hostMemSize = HLTCA_GPU_ROWS_MEMORY + HLTCA_GPU_COMMON_MEMORY + sliceCount * (HLTCA_GPU_SLICE_DATA_MEMORY + HLTCA_GPU_TRACKS_MEMORY) + HLTCA_GPU_TRACKER_OBJECT_MEMORY;
213 if (CudaFailedMsg(cudaMallocHost(&fHostLockedMemory, hostMemSize)))
215 cudaFree(fGPUMemory);
217 HLTError("Error allocating Page Locked Host Memory");
221 if (fDebugLevel >= 1) HLTInfo("Host Memory used: %d", hostMemSize);
223 if (fDebugLevel >= 1)
225 CudaFailedMsg(cudaMemset(fGPUMemory, 143, (size_t) fGPUMemSize));
228 fSliceCount = sliceCount;
229 //Don't run constructor / destructor here, this will be just local memcopy of Tracker in GPU Memory
230 fGpuTracker = (AliHLTTPCCATracker*) TrackerMemory(fHostLockedMemory, 0);
232 for (int i = 0;i < fgkNSlices;i++)
234 fSlaveTrackers[i].SetGPUTracker();
235 fSlaveTrackers[i].SetGPUTrackerCommonMemory((char*) CommonMemory(fHostLockedMemory, i));
236 fSlaveTrackers[i].SetGPUSliceDataMemory(SliceDataMemory(fHostLockedMemory, i), RowMemory(fHostLockedMemory, i));
239 fpCudaStreams = malloc(CAMath::Max(3, fSliceCount) * sizeof(cudaStream_t));
240 cudaStream_t* const cudaStreams = (cudaStream_t*) fpCudaStreams;
241 for (int i = 0;i < CAMath::Max(3, fSliceCount);i++)
243 if (CudaFailedMsg(cudaStreamCreate(&cudaStreams[i])))
245 cudaFree(fGPUMemory);
246 cudaFreeHost(fHostLockedMemory);
248 HLTError("Error creating CUDA Stream");
254 fCudaInitialized = 1;
255 HLTImportant("CUDA Initialisation successfull (Device %d: %s, Thread %dd)", cudaDevice, fCudaDeviceProp.name, fThreadId);
257 #if defined(HLTCA_STANDALONE) & !defined(CUDA_DEVICE_EMULATION)
260 //Do one initial run for Benchmark reasons
261 const int useDebugLevel = fDebugLevel;
263 AliHLTTPCCAClusterData* tmpCluster = new AliHLTTPCCAClusterData[sliceCount];
267 AliHLTTPCCAParam tmpParam;
268 AliHLTTPCCASliceOutput::outputControlStruct tmpOutputControl;
270 fin.open("events/settings.dump");
273 for (int i = 0;i < sliceCount;i++)
275 fSlaveTrackers[i].SetOutputControl(&tmpOutputControl);
276 tmpParam.ReadSettings(fin);
277 InitializeSliceParam(i, tmpParam);
281 fin.open("eventspbpbc/event.0.dump", std::ifstream::binary);
282 for (int i = 0;i < sliceCount;i++)
284 tmpCluster[i].StartReading(i, 0);
285 tmpCluster[i].ReadEvent(fin);
286 tmpCluster[i].FinishReading();
290 AliHLTTPCCASliceOutput **tmpOutput = new AliHLTTPCCASliceOutput*[sliceCount];
291 memset(tmpOutput, 0, sliceCount * sizeof(AliHLTTPCCASliceOutput*));
293 Reconstruct(tmpOutput, tmpCluster, 0, sliceCount);
294 for (int i = 0;i < sliceCount;i++)
298 fSlaveTrackers[i].SetOutputControl(NULL);
302 fDebugLevel = useDebugLevel;
308 template <class T> inline T* AliHLTTPCCAGPUTracker::alignPointer(T* ptr, int alignment)
310 //Macro to align Pointers.
311 //Will align to start at 1 MB segments, this should be consistent with every alignment in the tracker
312 //(As long as every single data structure is <= 1 MB)
314 size_t adr = (size_t) ptr;
317 adr += alignment - (adr % alignment);
322 bool AliHLTTPCCAGPUTracker::CudaFailedMsg(cudaError_t error)
324 //Check for CUDA Error and in the case of an error display the corresponding error string
325 if (error == cudaSuccess) return(false);
326 HLTWarning("CUDA Error: %d / %s", error, cudaGetErrorString(error));
330 int AliHLTTPCCAGPUTracker::CUDASync(char* state)
332 //Wait for CUDA-Kernel to finish and check for CUDA errors afterwards
334 if (fDebugLevel == 0) return(0);
336 cuErr = cudaGetLastError();
337 if (cuErr != cudaSuccess)
339 HLTError("Cuda Error %s while invoking kernel (%s)", cudaGetErrorString(cuErr), state);
342 if (CudaFailedMsg(cudaThreadSynchronize()))
344 HLTError("CUDA Error while synchronizing (%s)", state);
347 if (fDebugLevel >= 3) HLTInfo("CUDA Sync Done");
351 void AliHLTTPCCAGPUTracker::SetDebugLevel(const int dwLevel, std::ostream* const NewOutFile)
353 //Set Debug Level and Debug output File if applicable
354 fDebugLevel = dwLevel;
355 if (NewOutFile) fOutFile = NewOutFile;
358 int AliHLTTPCCAGPUTracker::SetGPUTrackerOption(char* OptionName, int /*OptionValue*/)
360 //Set a specific GPU Tracker Option
362 HLTError("Unknown Option: %s", OptionName);
365 //No Options used at the moment
369 #ifdef HLTCA_STANDALONE
370 void AliHLTTPCCAGPUTracker::StandalonePerfTime(int iSlice, int i)
372 //Run Performance Query for timer i of slice iSlice
373 if (fDebugLevel >= 1)
375 AliHLTTPCCAStandaloneFramework::StandaloneQueryTime( fSlaveTrackers[iSlice].PerfTimer(i));
379 void AliHLTTPCCAGPUTracker::StandalonePerfTime(int /*iSlice*/, int /*i*/) {}
382 void AliHLTTPCCAGPUTracker::DumpRowBlocks(AliHLTTPCCATracker* tracker, int iSlice, bool check)
384 //Dump Rowblocks to File
385 if (fDebugLevel >= 4)
387 *fOutFile << "RowBlock Tracklets" << std::endl;
389 int4* rowBlockPos = (int4*) malloc(sizeof(int4) * (tracker[iSlice].Param().NRows() / HLTCA_GPU_SCHED_ROW_STEP + 1) * 2);
390 int* rowBlockTracklets = (int*) malloc(sizeof(int) * (tracker[iSlice].Param().NRows() / HLTCA_GPU_SCHED_ROW_STEP + 1) * HLTCA_GPU_MAX_TRACKLETS * 2);
391 uint2* blockStartingTracklet = (uint2*) malloc(sizeof(uint2) * HLTCA_GPU_BLOCK_COUNT);
392 CudaFailedMsg(cudaMemcpy(rowBlockPos, fGpuTracker[iSlice].RowBlockPos(), sizeof(int4) * (tracker[iSlice].Param().NRows() / HLTCA_GPU_SCHED_ROW_STEP + 1) * 2, cudaMemcpyDeviceToHost));
393 CudaFailedMsg(cudaMemcpy(rowBlockTracklets, fGpuTracker[iSlice].RowBlockTracklets(), sizeof(int) * (tracker[iSlice].Param().NRows() / HLTCA_GPU_SCHED_ROW_STEP + 1) * HLTCA_GPU_MAX_TRACKLETS * 2, cudaMemcpyDeviceToHost));
394 CudaFailedMsg(cudaMemcpy(blockStartingTracklet, fGpuTracker[iSlice].BlockStartingTracklet(), sizeof(uint2) * HLTCA_GPU_BLOCK_COUNT, cudaMemcpyDeviceToHost));
395 CudaFailedMsg(cudaMemcpy(tracker[iSlice].CommonMemory(), fGpuTracker[iSlice].CommonMemory(), fGpuTracker[iSlice].CommonMemorySize(), cudaMemcpyDeviceToHost));
397 int k = tracker[iSlice].GPUParameters()->fScheduleFirstDynamicTracklet;
398 for (int i = 0; i < tracker[iSlice].Param().NRows() / HLTCA_GPU_SCHED_ROW_STEP + 1;i++)
400 *fOutFile << "Rowblock: " << i << ", up " << rowBlockPos[i].y << "/" << rowBlockPos[i].x << ", down " <<
401 rowBlockPos[tracker[iSlice].Param().NRows() / HLTCA_GPU_SCHED_ROW_STEP + 1 + i].y << "/" << rowBlockPos[tracker[iSlice].Param().NRows() / HLTCA_GPU_SCHED_ROW_STEP + 1 + i].x << std::endl << "Phase 1: ";
402 for (int j = 0;j < rowBlockPos[i].x;j++)
404 //Use Tracker Object to calculate Offset instead of fGpuTracker, since *fNTracklets of fGpuTracker points to GPU Mem!
405 *fOutFile << rowBlockTracklets[(tracker[iSlice].RowBlockTracklets(0, i) - tracker[iSlice].RowBlockTracklets(0, 0)) + j] << ", ";
406 if (check && rowBlockTracklets[(tracker[iSlice].RowBlockTracklets(0, i) - tracker[iSlice].RowBlockTracklets(0, 0)) + j] != k)
408 HLTError("Wrong starting Row Block %d, entry %d, is %d, should be %d", i, j, rowBlockTracklets[(tracker[iSlice].RowBlockTracklets(0, i) - tracker[iSlice].RowBlockTracklets(0, 0)) + j], k);
411 if (rowBlockTracklets[(tracker[iSlice].RowBlockTracklets(0, i) - tracker[iSlice].RowBlockTracklets(0, 0)) + j] == -1)
413 HLTError("Error, -1 Tracklet found");
416 *fOutFile << std::endl << "Phase 2: ";
417 for (int j = 0;j < rowBlockPos[tracker[iSlice].Param().NRows() / HLTCA_GPU_SCHED_ROW_STEP + 1 + i].x;j++)
419 *fOutFile << rowBlockTracklets[(tracker[iSlice].RowBlockTracklets(1, i) - tracker[iSlice].RowBlockTracklets(0, 0)) + j] << ", ";
421 *fOutFile << std::endl;
426 *fOutFile << "Starting Threads: (First Dynamic: " << tracker[iSlice].GPUParameters()->fScheduleFirstDynamicTracklet << ")" << std::endl;
427 for (int i = 0;i < HLTCA_GPU_BLOCK_COUNT;i++)
429 *fOutFile << i << ": " << blockStartingTracklet[i].x << " - " << blockStartingTracklet[i].y << std::endl;
434 free(rowBlockTracklets);
435 free(blockStartingTracklet);
439 __global__ void PreInitRowBlocks(int4* const RowBlockPos, int* const RowBlockTracklets, int* const SliceDataHitWeights, int nSliceDataHits)
441 //Initialize GPU RowBlocks and HitWeights
442 int4* const rowBlockTracklets4 = (int4*) RowBlockTracklets;
443 int4* const sliceDataHitWeights4 = (int4*) SliceDataHitWeights;
444 const int stride = blockDim.x * gridDim.x;
446 i0.x = i0.y = i0.z = i0.w = 0;
447 i1.x = i1.y = i1.z = i1.w = -1;
448 for (int i = blockIdx.x * blockDim.x + threadIdx.x;i < sizeof(int4) * 2 * (HLTCA_ROW_COUNT / HLTCA_GPU_SCHED_ROW_STEP + 1) / sizeof(int4);i += stride)
450 for (int i = blockIdx.x * blockDim.x + threadIdx.x;i < sizeof(int) * (HLTCA_ROW_COUNT / HLTCA_GPU_SCHED_ROW_STEP + 1) * HLTCA_GPU_MAX_TRACKLETS * 2 / sizeof(int4);i += stride)
451 rowBlockTracklets4[i] = i1;
452 for (int i = blockIdx.x * blockDim.x + threadIdx.x;i < nSliceDataHits * sizeof(int) / sizeof(int4);i += stride)
453 sliceDataHitWeights4[i] = i0;
456 int AliHLTTPCCAGPUTracker::Reconstruct(AliHLTTPCCASliceOutput** pOutput, AliHLTTPCCAClusterData* pClusterData, int firstSlice, int sliceCountLocal)
458 //Primary reconstruction function
459 cudaStream_t* const cudaStreams = (cudaStream_t*) fpCudaStreams;
461 if (sliceCountLocal == -1) sliceCountLocal = fSliceCount;
463 if (!fCudaInitialized)
465 HLTError("GPUTracker not initialized");
468 if (sliceCountLocal > fSliceCount)
470 HLTError("GPU Tracker was initialized to run with %d slices but was called to process %d slices", fSliceCount, sliceCountLocal);
473 if (fThreadId != GetThread())
475 HLTError("GPUTracker context was initialized by different thread, Initializing Thread: %d, Processing Thread: %d", fThreadId, GetThread());
479 if (fDebugLevel >= 2) HLTInfo("Running GPU Tracker (Slices %d to %d)", fSlaveTrackers[firstSlice].Param().ISlice(), fSlaveTrackers[firstSlice].Param().ISlice() + sliceCountLocal);
481 if (sliceCountLocal * sizeof(AliHLTTPCCATracker) > HLTCA_GPU_TRACKER_CONSTANT_MEM)
483 HLTError("Insuffissant constant memory (Required %d, Available %d, Tracker %d, Param %d, SliceData %d)", sliceCountLocal * (int) sizeof(AliHLTTPCCATracker), (int) HLTCA_GPU_TRACKER_CONSTANT_MEM, (int) sizeof(AliHLTTPCCATracker), (int) sizeof(AliHLTTPCCAParam), (int) sizeof(AliHLTTPCCASliceData));
487 if (fDebugLevel >= 4)
489 for (int iSlice = 0;iSlice < sliceCountLocal;iSlice++)
491 *fOutFile << std::endl << std::endl << "Slice: " << fSlaveTrackers[firstSlice + iSlice].Param().ISlice() << std::endl;
495 memcpy(fGpuTracker, &fSlaveTrackers[firstSlice], sizeof(AliHLTTPCCATracker) * sliceCountLocal);
497 if (fDebugLevel >= 3) HLTInfo("Allocating GPU Tracker memory and initializing constants");
499 #ifdef HLTCA_GPU_TIME_PROFILE
500 unsigned __int64 a, b, c, d;
501 AliHLTTPCCAStandaloneFramework::StandaloneQueryFreq(&c);
502 AliHLTTPCCAStandaloneFramework::StandaloneQueryTime(&d);
505 for (int iSlice = 0;iSlice < sliceCountLocal;iSlice++)
507 //Make this a GPU Tracker
508 fGpuTracker[iSlice].SetGPUTracker();
509 fGpuTracker[iSlice].SetGPUTrackerCommonMemory((char*) CommonMemory(fGPUMemory, iSlice));
510 fGpuTracker[iSlice].SetGPUSliceDataMemory(SliceDataMemory(fGPUMemory, iSlice), RowMemory(fGPUMemory, iSlice));
511 fGpuTracker[iSlice].SetPointersSliceData(&pClusterData[iSlice], false);
513 //Set Pointers to GPU Memory
514 char* tmpMem = (char*) GlobalMemory(fGPUMemory, iSlice);
516 if (fDebugLevel >= 3) HLTInfo("Initialising GPU Hits Memory");
517 tmpMem = fGpuTracker[iSlice].SetGPUTrackerHitsMemory(tmpMem, pClusterData[iSlice].NumberOfClusters());
518 tmpMem = alignPointer(tmpMem, 1024 * 1024);
520 if (fDebugLevel >= 3) HLTInfo("Initialising GPU Tracklet Memory");
521 tmpMem = fGpuTracker[iSlice].SetGPUTrackerTrackletsMemory(tmpMem, HLTCA_GPU_MAX_TRACKLETS /* *fSlaveTrackers[firstSlice + iSlice].NTracklets()*/);
522 tmpMem = alignPointer(tmpMem, 1024 * 1024);
524 if (fDebugLevel >= 3) HLTInfo("Initialising GPU Track Memory");
525 tmpMem = fGpuTracker[iSlice].SetGPUTrackerTracksMemory(tmpMem, HLTCA_GPU_MAX_TRACKS /* *fSlaveTrackers[firstSlice + iSlice].NTracklets()*/, pClusterData[iSlice].NumberOfClusters());
526 tmpMem = alignPointer(tmpMem, 1024 * 1024);
528 if (fGpuTracker[iSlice].TrackMemorySize() >= HLTCA_GPU_TRACKS_MEMORY)
530 HLTError("Insufficiant Track Memory");
534 if (tmpMem - (char*) GlobalMemory(fGPUMemory, iSlice) > HLTCA_GPU_GLOBAL_MEMORY)
536 HLTError("Insufficiant Global Memory");
540 //Initialize Startup Constants
541 *fSlaveTrackers[firstSlice + iSlice].NTracklets() = 0;
542 *fSlaveTrackers[firstSlice + iSlice].NTracks() = 0;
543 *fSlaveTrackers[firstSlice + iSlice].NTrackHits() = 0;
544 fGpuTracker[iSlice].GPUParametersConst()->fGPUFixedBlockCount = HLTCA_GPU_BLOCK_COUNT * (iSlice + 1) / sliceCountLocal - HLTCA_GPU_BLOCK_COUNT * (iSlice) / sliceCountLocal;
545 if (fDebugLevel >= 3) HLTInfo("Blocks for Slice %d: %d", iSlice, fGpuTracker[iSlice].GPUParametersConst()->fGPUFixedBlockCount);
546 fGpuTracker[iSlice].GPUParametersConst()->fGPUiSlice = iSlice;
547 fGpuTracker[iSlice].GPUParametersConst()->fGPUnSlices = sliceCountLocal;
548 fSlaveTrackers[firstSlice + iSlice].GPUParameters()->fGPUError = 0;
549 fGpuTracker[iSlice].SetGPUTextureBase(fGpuTracker[0].Data().Memory());
552 #ifdef HLTCA_GPU_TEXTURE_FETCH
553 cudaChannelFormatDesc channelDescu2 = cudaCreateChannelDesc<ushort2>();
555 if (CudaFailedMsg(cudaBindTexture(&offset, &gAliTexRefu2, fGpuTracker[0].Data().Memory(), &channelDescu2, sliceCountLocal * HLTCA_GPU_SLICE_DATA_MEMORY)) || offset)
557 HLTError("Error binding CUDA Texture ushort2 (Offset %d)", (int) offset);
560 cudaChannelFormatDesc channelDescu = cudaCreateChannelDesc<unsigned short>();
561 if (CudaFailedMsg(cudaBindTexture(&offset, &gAliTexRefu, fGpuTracker[0].Data().Memory(), &channelDescu, sliceCountLocal * HLTCA_GPU_SLICE_DATA_MEMORY)) || offset)
563 HLTError("Error binding CUDA Texture ushort (Offset %d)", (int) offset);
566 cudaChannelFormatDesc channelDescs = cudaCreateChannelDesc<signed short>();
567 if (CudaFailedMsg(cudaBindTexture(&offset, &gAliTexRefs, fGpuTracker[0].Data().Memory(), &channelDescs, sliceCountLocal * HLTCA_GPU_SLICE_DATA_MEMORY)) || offset)
569 HLTError("Error binding CUDA Texture short (Offset %d)", (int) offset);
574 //Copy Tracker Object to GPU Memory
575 if (fDebugLevel >= 3) HLTInfo("Copying Tracker objects to GPU");
576 #ifdef HLTCA_GPU_TRACKLET_CONSTRUCTOR_DO_PROFILE
577 if (CudaFailedMsg(cudaMalloc(&fGpuTracker[0].fStageAtSync, 100000000))) return(1);
578 CudaFailedMsg(cudaMemset(fGpuTracker[0].fStageAtSync, 0, 100000000));
580 CudaFailedMsg(cudaMemcpyToSymbolAsync(gAliHLTTPCCATracker, fGpuTracker, sizeof(AliHLTTPCCATracker) * sliceCountLocal, 0, cudaMemcpyHostToDevice, cudaStreams[0]));
582 for (int iSlice = 0;iSlice < sliceCountLocal;iSlice++)
584 StandalonePerfTime(firstSlice + iSlice, 0);
586 //Initialize GPU Slave Tracker
587 if (fDebugLevel >= 3) HLTInfo("Creating Slice Data");
588 fSlaveTrackers[firstSlice + iSlice].SetGPUSliceDataMemory(SliceDataMemory(fHostLockedMemory, iSlice), RowMemory(fHostLockedMemory, firstSlice + iSlice));
589 #ifdef HLTCA_GPU_TIME_PROFILE
590 AliHLTTPCCAStandaloneFramework::StandaloneQueryTime(&a);
592 fSlaveTrackers[firstSlice + iSlice].ReadEvent(&pClusterData[iSlice]);
593 #ifdef HLTCA_GPU_TIME_PROFILE
594 AliHLTTPCCAStandaloneFramework::StandaloneQueryTime(&b);
595 printf("Read %f %f\n", ((double) b - (double) a) / (double) c, ((double) a - (double) d) / (double) c);
597 if (fSlaveTrackers[firstSlice + iSlice].Data().MemorySize() > HLTCA_GPU_SLICE_DATA_MEMORY)
599 HLTError("Insufficiant Slice Data Memory");
603 /*if (fSlaveTrackers[firstSlice + iSlice].CheckEmptySlice())
605 if (fDebugLevel >= 3) HLTInfo("Slice Empty, not running GPU Tracker");
606 if (sliceCountLocal == 1)
610 //Initialize temporary memory where needed
611 if (fDebugLevel >= 3) HLTInfo("Copying Slice Data to GPU and initializing temporary memory");
612 PreInitRowBlocks<<<30, 256, 0, cudaStreams[2]>>>(fGpuTracker[iSlice].RowBlockPos(), fGpuTracker[iSlice].RowBlockTracklets(), fGpuTracker[iSlice].Data().HitWeights(), fSlaveTrackers[firstSlice + iSlice].Data().NumberOfHitsPlusAlign());
614 //Copy Data to GPU Global Memory
615 CudaFailedMsg(cudaMemcpyAsync(fGpuTracker[iSlice].CommonMemory(), fSlaveTrackers[firstSlice + iSlice].CommonMemory(), fSlaveTrackers[firstSlice + iSlice].CommonMemorySize(), cudaMemcpyHostToDevice, cudaStreams[iSlice & 1]));
616 CudaFailedMsg(cudaMemcpyAsync(fGpuTracker[iSlice].Data().Memory(), fSlaveTrackers[firstSlice + iSlice].Data().Memory(), fSlaveTrackers[firstSlice + iSlice].Data().GpuMemorySize(), cudaMemcpyHostToDevice, cudaStreams[iSlice & 1]));
617 CudaFailedMsg(cudaMemcpyAsync(fGpuTracker[iSlice].SliceDataRows(), fSlaveTrackers[firstSlice + iSlice].SliceDataRows(), (HLTCA_ROW_COUNT + 1) * sizeof(AliHLTTPCCARow), cudaMemcpyHostToDevice, cudaStreams[iSlice & 1]));
619 if (fDebugLevel >= 4)
621 if (fDebugLevel >= 5) HLTInfo("Allocating Debug Output Memory");
622 fSlaveTrackers[firstSlice + iSlice].TrackletMemory() = reinterpret_cast<char*> ( new uint4 [ fGpuTracker[iSlice].TrackletMemorySize()/sizeof( uint4 ) + 100] );
623 fSlaveTrackers[firstSlice + iSlice].SetPointersTracklets( HLTCA_GPU_MAX_TRACKLETS );
624 fSlaveTrackers[firstSlice + iSlice].HitMemory() = reinterpret_cast<char*> ( new uint4 [ fGpuTracker[iSlice].HitMemorySize()/sizeof( uint4 ) + 100] );
625 fSlaveTrackers[firstSlice + iSlice].SetPointersHits( pClusterData[iSlice].NumberOfClusters() );
628 if (CUDASync("Initialization")) return(1);
629 StandalonePerfTime(firstSlice + iSlice, 1);
631 if (fDebugLevel >= 3) HLTInfo("Running GPU Neighbours Finder");
632 AliHLTTPCCAProcess<AliHLTTPCCANeighboursFinder> <<<fSlaveTrackers[firstSlice + iSlice].Param().NRows(), 256, 0, cudaStreams[iSlice & 1]>>>(iSlice);
634 if (CUDASync("Neighbours finder")) return 1;
636 StandalonePerfTime(firstSlice + iSlice, 2);
638 if (fDebugLevel >= 4)
640 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].Data().Memory(), fGpuTracker[iSlice].Data().Memory(), fSlaveTrackers[firstSlice + iSlice].Data().GpuMemorySize(), cudaMemcpyDeviceToHost));
641 fSlaveTrackers[firstSlice + iSlice].DumpLinks(*fOutFile);
644 if (fDebugLevel >= 3) HLTInfo("Running GPU Neighbours Cleaner");
645 AliHLTTPCCAProcess<AliHLTTPCCANeighboursCleaner> <<<fSlaveTrackers[firstSlice + iSlice].Param().NRows()-2, 256, 0, cudaStreams[iSlice & 1]>>>(iSlice);
646 if (CUDASync("Neighbours Cleaner")) return 1;
648 StandalonePerfTime(firstSlice + iSlice, 3);
650 if (fDebugLevel >= 4)
652 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].Data().Memory(), fGpuTracker[iSlice].Data().Memory(), fSlaveTrackers[firstSlice + iSlice].Data().GpuMemorySize(), cudaMemcpyDeviceToHost));
653 fSlaveTrackers[firstSlice + iSlice].DumpLinks(*fOutFile);
656 if (fDebugLevel >= 3) HLTInfo("Running GPU Start Hits Finder");
657 AliHLTTPCCAProcess<AliHLTTPCCAStartHitsFinder> <<<fSlaveTrackers[firstSlice + iSlice].Param().NRows()-6, 256, 0, cudaStreams[iSlice & 1]>>>(iSlice);
658 if (CUDASync("Start Hits Finder")) return 1;
660 StandalonePerfTime(firstSlice + iSlice, 4);
662 if (fDebugLevel >= 3) HLTInfo("Running GPU Start Hits Sorter");
663 AliHLTTPCCAProcess<AliHLTTPCCAStartHitsSorter> <<<30, 256, 0, cudaStreams[iSlice & 1]>>>(iSlice);
664 if (CUDASync("Start Hits Sorter")) return 1;
666 StandalonePerfTime(firstSlice + iSlice, 5);
668 if (fDebugLevel >= 2)
670 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].CommonMemory(), fGpuTracker[iSlice].CommonMemory(), fGpuTracker[iSlice].CommonMemorySize(), cudaMemcpyDeviceToHost));
671 if (fDebugLevel >= 3) HLTInfo("Obtaining Number of Start Hits from GPU: %d", *fSlaveTrackers[firstSlice + iSlice].NTracklets());
672 if (*fSlaveTrackers[firstSlice + iSlice].NTracklets() > HLTCA_GPU_MAX_TRACKLETS)
674 HLTError("HLTCA_GPU_MAX_TRACKLETS constant insuffisant");
679 if (fDebugLevel >= 4)
681 *fOutFile << "Temporary ";
682 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].TrackletStartHits(), fGpuTracker[iSlice].TrackletTmpStartHits(), pClusterData[iSlice].NumberOfClusters() * sizeof(AliHLTTPCCAHitId), cudaMemcpyDeviceToHost));
683 fSlaveTrackers[firstSlice + iSlice].DumpStartHits(*fOutFile);
684 uint3* tmpMemory = (uint3*) malloc(sizeof(uint3) * fSlaveTrackers[firstSlice + iSlice].Param().NRows());
685 CudaFailedMsg(cudaMemcpy(tmpMemory, fGpuTracker[iSlice].RowStartHitCountOffset(), fSlaveTrackers[firstSlice + iSlice].Param().NRows() * sizeof(uint3), cudaMemcpyDeviceToHost));
686 *fOutFile << "Start Hits Sort Vector:" << std::endl;
687 for (int i = 0;i < fSlaveTrackers[firstSlice + iSlice].Param().NRows();i++)
689 *fOutFile << "Row: " << i << ", Len: " << tmpMemory[i].x << ", Offset: " << tmpMemory[i].y << ", New Offset: " << tmpMemory[i].z << std::endl;
693 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].HitMemory(), fGpuTracker[iSlice].HitMemory(), fSlaveTrackers[firstSlice + iSlice].HitMemorySize(), cudaMemcpyDeviceToHost));
694 fSlaveTrackers[firstSlice + iSlice].DumpStartHits(*fOutFile);
697 StandalonePerfTime(firstSlice + iSlice, 6);
699 fSlaveTrackers[firstSlice + iSlice].SetGPUTrackerTracksMemory((char*) TracksMemory(fHostLockedMemory, iSlice), HLTCA_GPU_MAX_TRACKS, pClusterData[iSlice].NumberOfClusters());
702 StandalonePerfTime(firstSlice, 7);
703 #ifdef HLTCA_GPU_PREFETCHDATA
704 for (int iSlice = 0;iSlice < sliceCountLocal;iSlice++)
706 if (fSlaveTrackers[firstSlice + iSlice].Data().GPUSharedDataReq() * sizeof(ushort_v) > ALIHLTTPCCATRACKLET_CONSTRUCTOR_TEMP_MEM / 4 * sizeof(uint4))
708 HLTError("Insufficiant GPU shared Memory, required: %d, available %d", fSlaveTrackers[firstSlice + iSlice].Data().GPUSharedDataReq() * sizeof(ushort_v), ALIHLTTPCCATRACKLET_CONSTRUCTOR_TEMP_MEM / 4 * sizeof(uint4));
711 if (fDebugLevel >= 1)
713 static int infoShown = 0;
716 HLTInfo("GPU Shared Memory Cache Size: %d", 2 * fSlaveTrackers[firstSlice + iSlice].Data().GPUSharedDataReq() * sizeof(ushort_v));
723 int nHardCollisions = 0;
725 RestartTrackletConstructor:
726 if (fDebugLevel >= 3) HLTInfo("Initialising Tracklet Constructor Scheduler");
727 for (int iSlice = 0;iSlice < sliceCountLocal;iSlice++)
729 AliHLTTPCCATrackletConstructorInit<<<HLTCA_GPU_MAX_TRACKLETS /* *fSlaveTrackers[firstSlice + iSlice].NTracklets() */ / HLTCA_GPU_THREAD_COUNT + 1, HLTCA_GPU_THREAD_COUNT>>>(iSlice);
730 if (CUDASync("Tracklet Initializer")) return 1;
731 DumpRowBlocks(fSlaveTrackers, iSlice);
734 if (fDebugLevel >= 3) HLTInfo("Running GPU Tracklet Constructor");
735 AliHLTTPCCATrackletConstructorNewGPU<<<HLTCA_GPU_BLOCK_COUNT, HLTCA_GPU_THREAD_COUNT>>>();
736 if (CUDASync("Tracklet Constructor (new)")) return 1;
738 StandalonePerfTime(firstSlice, 8);
740 if (fDebugLevel >= 4)
742 for (int iSlice = 0;iSlice < sliceCountLocal;iSlice++)
744 DumpRowBlocks(&fSlaveTrackers[firstSlice], iSlice, false);
745 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].CommonMemory(), fGpuTracker[iSlice].CommonMemory(), fGpuTracker[iSlice].CommonMemorySize(), cudaMemcpyDeviceToHost));
746 if (fDebugLevel >= 5)
748 HLTInfo("Obtained %d tracklets", *fSlaveTrackers[firstSlice + iSlice].NTracklets());
750 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].TrackletMemory(), fGpuTracker[iSlice].TrackletMemory(), fGpuTracker[iSlice].TrackletMemorySize(), cudaMemcpyDeviceToHost));
751 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].HitMemory(), fGpuTracker[iSlice].HitMemory(), fGpuTracker[iSlice].HitMemorySize(), cudaMemcpyDeviceToHost));
752 fSlaveTrackers[firstSlice + iSlice].DumpTrackletHits(*fOutFile);
757 for (int iSlice = 0;iSlice < sliceCountLocal;iSlice += runSlices)
759 if (runSlices < HLTCA_GPU_TRACKLET_SELECTOR_SLICE_COUNT) runSlices++;
760 if (fDebugLevel >= 3) HLTInfo("Running HLT Tracklet selector (Slice %d to %d)", iSlice, iSlice + runSlices);
761 AliHLTTPCCAProcessMulti<AliHLTTPCCATrackletSelector><<<HLTCA_GPU_BLOCK_COUNT, HLTCA_GPU_THREAD_COUNT, 0, cudaStreams[iSlice]>>>(iSlice, CAMath::Min(runSlices, sliceCountLocal - iSlice));
763 if (CUDASync("Tracklet Selector")) return 1;
764 StandalonePerfTime(firstSlice, 9);
766 CudaFailedMsg(cudaMemcpyAsync(fSlaveTrackers[firstSlice + 0].CommonMemory(), fGpuTracker[0].CommonMemory(), fGpuTracker[0].CommonMemorySize(), cudaMemcpyDeviceToHost, cudaStreams[0]));
767 for (int iSliceTmp = 0;iSliceTmp <= sliceCountLocal;iSliceTmp++)
769 if (iSliceTmp < sliceCountLocal)
771 int iSlice = iSliceTmp;
772 if (fDebugLevel >= 3) HLTInfo("Transfering Tracks from GPU to Host");
773 cudaStreamSynchronize(cudaStreams[iSlice]);
774 CudaFailedMsg(cudaMemcpyAsync(fSlaveTrackers[firstSlice + iSlice].Tracks(), fGpuTracker[iSlice].Tracks(), sizeof(AliHLTTPCCATrack) * *fSlaveTrackers[firstSlice + iSlice].NTracks(), cudaMemcpyDeviceToHost, cudaStreams[iSlice]));
775 CudaFailedMsg(cudaMemcpyAsync(fSlaveTrackers[firstSlice + iSlice].TrackHits(), fGpuTracker[iSlice].TrackHits(), sizeof(AliHLTTPCCAHitId) * *fSlaveTrackers[firstSlice + iSlice].NTrackHits(), cudaMemcpyDeviceToHost, cudaStreams[iSlice]));
776 if (iSlice + 1 < sliceCountLocal)
777 CudaFailedMsg(cudaMemcpyAsync(fSlaveTrackers[firstSlice + iSlice + 1].CommonMemory(), fGpuTracker[iSlice + 1].CommonMemory(), fGpuTracker[iSlice + 1].CommonMemorySize(), cudaMemcpyDeviceToHost, cudaStreams[iSlice + 1]));
782 int iSlice = iSliceTmp - 1;
783 cudaStreamSynchronize(cudaStreams[iSlice]);
785 if (fDebugLevel >= 4)
787 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + iSlice].Data().HitWeights(), fGpuTracker[iSlice].Data().HitWeights(), fSlaveTrackers[firstSlice + iSlice].Data().NumberOfHitsPlusAlign() * sizeof(int), cudaMemcpyDeviceToHost));
788 fSlaveTrackers[firstSlice + iSlice].DumpHitWeights(*fOutFile);
789 fSlaveTrackers[firstSlice + iSlice].DumpTrackHits(*fOutFile);
792 if (fSlaveTrackers[firstSlice + iSlice].GPUParameters()->fGPUError)
794 if (fSlaveTrackers[firstSlice + iSlice].GPUParameters()->fGPUError == HLTCA_GPU_ERROR_SCHEDULE_COLLISION && nHardCollisions++ < 10)
796 HLTWarning("Hard scheduling collision occured, rerunning Tracklet Constructor");
797 for (int i = 0;i < sliceCountLocal;i++)
799 cudaThreadSynchronize();
800 CudaFailedMsg(cudaMemcpy(fSlaveTrackers[firstSlice + i].CommonMemory(), fGpuTracker[i].CommonMemory(), fGpuTracker[i].CommonMemorySize(), cudaMemcpyDeviceToHost));
801 *fSlaveTrackers[firstSlice + i].NTracks() = 0;
802 *fSlaveTrackers[firstSlice + i].NTrackHits() = 0;
803 fSlaveTrackers[firstSlice + i].GPUParameters()->fGPUError = HLTCA_GPU_ERROR_NONE;
804 CudaFailedMsg(cudaMemcpy(fGpuTracker[i].CommonMemory(), fSlaveTrackers[firstSlice + i].CommonMemory(), fGpuTracker[i].CommonMemorySize(), cudaMemcpyHostToDevice));
805 PreInitRowBlocks<<<30, 256>>>(fGpuTracker[i].RowBlockPos(), fGpuTracker[i].RowBlockTracklets(), fGpuTracker[i].Data().HitWeights(), fSlaveTrackers[firstSlice + i].Data().NumberOfHitsPlusAlign());
807 goto RestartTrackletConstructor;
809 HLTError("GPU Tracker returned Error Code %d", fSlaveTrackers[firstSlice + iSlice].GPUParameters()->fGPUError);
812 if (fDebugLevel >= 3) HLTInfo("Tracks Transfered: %d / %d", *fSlaveTrackers[firstSlice + iSlice].NTracks(), *fSlaveTrackers[firstSlice + iSlice].NTrackHits());
814 fSlaveTrackers[firstSlice + iSlice].SetOutput(&pOutput[iSlice]);
815 #ifdef HLTCA_GPU_TIME_PROFILE
816 AliHLTTPCCAStandaloneFramework::StandaloneQueryTime(&a);
818 fSlaveTrackers[firstSlice + iSlice].WriteOutput();
819 #ifdef HLTCA_GPU_TIME_PROFILE
820 AliHLTTPCCAStandaloneFramework::StandaloneQueryTime(&b);
821 printf("Write %f %f\n", ((double) b - (double) a) / (double) c, ((double) a - (double) d) / (double) c);
824 if (fDebugLevel >= 4)
826 delete[] fSlaveTrackers[firstSlice + iSlice].HitMemory();
827 delete[] fSlaveTrackers[firstSlice + iSlice].TrackletMemory();
832 StandalonePerfTime(firstSlice, 10);
834 if (fDebugLevel >= 3) HLTInfo("GPU Reconstruction finished");
836 #ifdef HLTCA_GPU_TRACKLET_CONSTRUCTOR_DO_PROFILE
837 char* stageAtSync = (char*) malloc(100000000);
838 CudaFailedMsg(cudaMemcpy(stageAtSync, fGpuTracker[0].fStageAtSync, 100 * 1000 * 1000, cudaMemcpyDeviceToHost));
839 cudaFree(fGpuTracker[0].fStageAtSync);
841 FILE* fp = fopen("profile.txt", "w+");
842 FILE* fp2 = fopen("profile.bmp", "w+b");
843 int nEmptySync = 0, fEmpty;
845 const int bmpheight = 1000;
846 BITMAPFILEHEADER bmpFH;
847 BITMAPINFOHEADER bmpIH;
848 ZeroMemory(&bmpFH, sizeof(bmpFH));
849 ZeroMemory(&bmpIH, sizeof(bmpIH));
851 bmpFH.bfType = 19778; //"BM"
852 bmpFH.bfSize = sizeof(bmpFH) + sizeof(bmpIH) + (HLTCA_GPU_BLOCK_COUNT * HLTCA_GPU_THREAD_COUNT / 32 * 33 - 1) * bmpheight ;
853 bmpFH.bfOffBits = sizeof(bmpFH) + sizeof(bmpIH);
855 bmpIH.biSize = sizeof(bmpIH);
856 bmpIH.biWidth = HLTCA_GPU_BLOCK_COUNT * HLTCA_GPU_THREAD_COUNT / 32 * 33 - 1;
857 bmpIH.biHeight = bmpheight;
859 bmpIH.biBitCount = 32;
861 fwrite(&bmpFH, 1, sizeof(bmpFH), fp2);
862 fwrite(&bmpIH, 1, sizeof(bmpIH), fp2);
864 for (int i = 0;i < bmpheight * HLTCA_GPU_BLOCK_COUNT * HLTCA_GPU_THREAD_COUNT;i += HLTCA_GPU_BLOCK_COUNT * HLTCA_GPU_THREAD_COUNT)
867 for (int j = 0;j < HLTCA_GPU_BLOCK_COUNT * HLTCA_GPU_THREAD_COUNT;j++)
869 fprintf(fp, "%d\t", stageAtSync[i + j]);
871 if (stageAtSync[i + j] == 1) color = RGB(255, 0, 0);
872 if (stageAtSync[i + j] == 2) color = RGB(0, 255, 0);
873 if (stageAtSync[i + j] == 3) color = RGB(0, 0, 255);
874 if (stageAtSync[i + j] == 4) color = RGB(255, 255, 0);
875 fwrite(&color, 1, sizeof(int), fp2);
876 if (j > 0 && j % 32 == 0)
878 color = RGB(255, 255, 255);
879 fwrite(&color, 1, 4, fp2);
881 if (stageAtSync[i + j]) fEmpty = 0;
884 if (fEmpty) nEmptySync++;
886 //if (nEmptySync == HLTCA_GPU_SCHED_ROW_STEP + 2) break;
897 int AliHLTTPCCAGPUTracker::InitializeSliceParam(int iSlice, AliHLTTPCCAParam ¶m)
899 //Initialize Slice Tracker Parameter for a slave tracker
900 fSlaveTrackers[iSlice].Initialize(param);
901 if (fSlaveTrackers[iSlice].Param().NRows() != HLTCA_ROW_COUNT)
903 HLTError("Error, Slice Tracker %d Row Count of %d exceeds Constant of %d", iSlice, fSlaveTrackers[iSlice].Param().NRows(), HLTCA_ROW_COUNT);
909 int AliHLTTPCCAGPUTracker::ExitGPU()
912 cudaThreadSynchronize();
915 cudaFree(fGPUMemory);
918 if (fHostLockedMemory)
920 for (int i = 0;i < CAMath::Max(3, fSliceCount);i++)
922 cudaStreamDestroy(((cudaStream_t*) fpCudaStreams)[i]);
926 cudaFreeHost(fHostLockedMemory);
929 if (CudaFailedMsg(cudaThreadExit()))
931 HLTError("Could not uninitialize GPU");
934 HLTInfo("CUDA Uninitialized");
936 fCudaInitialized = 0;
940 void AliHLTTPCCAGPUTracker::SetOutputControl( AliHLTTPCCASliceOutput::outputControlStruct* val)
942 fOutputControl = val;
943 for (int i = 0;i < fgkNSlices;i++)
945 fSlaveTrackers[i].SetOutputControl(val);
949 int AliHLTTPCCAGPUTracker::GetThread()
952 return((int) (size_t) GetCurrentThread());
954 return((int) syscall (SYS_gettid));