]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/TRD/AliHLTTRDPreprocessorComponent.cxx
cleaning up HLT defines
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDPreprocessorComponent.cxx
CommitLineData
c7b7f445 1// $Id$
2//**************************************************************************
3//* This file is property of and copyright by the ALICE HLT Project *
4//* ALICE Experiment at CERN, All rights reserved. *
5//* *
6//* Permission to use, copy, modify and distribute this software and its *
7//* documentation strictly for non-commercial purposes is hereby granted *
8//* without fee, provided that the above copyright notice appears in all *
9//* copies and that both the copyright notice and this permission notice *
10//* appear in the supporting documentation. The authors make no claims *
11//* about the suitability of this software for any purpose. It is *
12//* provided "as is" without express or implied warranty. *
13//**************************************************************************
14
15/// @file AliHLTTRDPreprocessorComponent.cxx
16/// @author Felix Rettig, Stefan Kirsch
17/// @date 2012-08-16
18/// @brief A pre-processing component for TRD tracking/trigger data on FEP-level
19/// @ingroup alihlt_trd_components
20
21#include <cstdlib>
22#include "TList.h"
23#include "AliLog.h"
24#include "AliHLTDataTypes.h"
25#include "AliHLTTRDDefinitions.h"
26#include "AliHLTTRDPreprocessorComponent.h"
27#include "AliTRDdigitsManager.h"
28#include "AliRawReaderMemory.h"
29#include "AliTRDrawStream.h"
30#include "AliTRDtrackletWord.h"
31#include "AliESDTrdTracklet.h"
32#include "AliESDTrdTrack.h"
33#include "AliTRDonlineTrackingDataContainer.h"
34
35ClassImp(AliHLTTRDPreprocessorComponent)
36
37#define LogError( ... ) { HLTError(__VA_ARGS__); if (fDebugLevel >= 1) { DbgLog("ERROR", __VA_ARGS__); } }
38#define LogInfo( ... ) { HLTInfo(__VA_ARGS__); if (fDebugLevel >= 1) { DbgLog("INFO", __VA_ARGS__); } }
39#define LogInspect( ... ) { HLTDebug(__VA_ARGS__); if (fDebugLevel >= 1) { DbgLog("INSPECT", __VA_ARGS__); } }
40#define LogDebug( ... ) { if (fDebugLevel >= 1) { HLTInfo(__VA_ARGS__); DbgLog("DEBUG", __VA_ARGS__); } }
41
42AliHLTTRDPreprocessorComponent::AliHLTTRDPreprocessorComponent() :
43 AliHLTProcessor(),
44 fDebugLevel(0),
45 fEventId(fgkInvalidEventId),
46 fTrackletArray(NULL),
47 fGtuTrackArray(NULL),
48 fRawReaderMem(NULL),
49 fDigitsManagerTrd(NULL),
50 fRawReaderTrd(NULL),
51 fTrackingData(NULL)
52{
53 // constructor
54}
55
56AliHLTTRDPreprocessorComponent::~AliHLTTRDPreprocessorComponent() {
57 // destructor
58}
59
60const char* AliHLTTRDPreprocessorComponent::GetComponentID() {
61 return "TRDPreprocessorComponent";
62}
63
64void AliHLTTRDPreprocessorComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
65 list.push_back(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
66}
67
68AliHLTComponentDataType AliHLTTRDPreprocessorComponent::GetOutputDataType() {
69 return kAliHLTMultipleDataType;
70}
71
72int AliHLTTRDPreprocessorComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList) {
73 tgtList.clear();
74 tgtList.push_back(kAliHLTDataTypeTObject | kAliHLTDataOriginTRD);
75 return tgtList.size();
76}
77
78void AliHLTTRDPreprocessorComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier ) {
79 constBase = 5000000;
80 inputMultiplier = 0;
81}
82
83void AliHLTTRDPreprocessorComponent::GetOCDBObjectDescription( TMap* const /*targetMap*/) {
84}
85
86AliHLTComponent* AliHLTTRDPreprocessorComponent::Spawn(){
87 return new AliHLTTRDPreprocessorComponent;
88}
89
90int AliHLTTRDPreprocessorComponent::Reconfigure(const char* /*cdbEntry*/, const char* /*chainId*/) {
91 return 0;
92}
93
94int AliHLTTRDPreprocessorComponent::ReadPreprocessorValues(const char* /*modules*/){
95 return 0;
96}
97
98int AliHLTTRDPreprocessorComponent::ScanConfigurationArgument(int argc, const char** argv){
99
100 if (argc <= 0)
101 return 0;
102
103 UShort_t iArg = 0;
104 TString argument(argv[iArg]);
105
106 if (!argument.CompareTo("-debug")){
107 if (++iArg >= argc) return -EPROTO;
108 argument = argv[iArg];
109 fDebugLevel = argument.Atoi();
110 LogInfo("debug level set to %d.", fDebugLevel);
111 return 2;
112 }
113
114 return 0;
115}
116
117int AliHLTTRDPreprocessorComponent::DoInit(int argc, const char** argv){
118
119 int iResult = 0;
120
121 do {
122
123 fRawReaderMem = new AliRawReaderMemory;
124 if (!fRawReaderMem) {
125 iResult=-ENOMEM;
126 break;
127 }
128
129 fTrackletArray = new TClonesArray("AliTRDtrackletWord", 500);
130 if (!fTrackletArray) {
131 iResult=-ENOMEM;
132 break;
133 }
134 fTrackletArray->BypassStreamer(kTRUE); //## needed for performance improvement?
135
136 fGtuTrackArray = new TClonesArray("AliESDTrdTrack", 50);
137 if (!fGtuTrackArray){
138 iResult=-ENOMEM;
139 break;
140 }
141
142 fRawReaderTrd = new AliTRDrawStream(fRawReaderMem);
143 if (!fRawReaderTrd) {
144 iResult=-ENOMEM;
145 break;
146 }
147
148 if (kFALSE){ // do not create digits manager for reader -> do not process ADC raw data
149 fDigitsManagerTrd = new AliTRDdigitsManager();
150 if (!fDigitsManagerTrd) {
151 iResult=-ENOMEM;
152 break;
153 }
154 fDigitsManagerTrd->CreateArrays();
155 fRawReaderTrd->SetDigitsManager(fDigitsManagerTrd);
156 }
157
158 fRawReaderTrd->SetDigitsManager(fDigitsManagerTrd);
159 fRawReaderTrd->SetTrackletArray(fTrackletArray);
160 fRawReaderTrd->SetTrackArray(fGtuTrackArray);
161
162 // Disable raw reader error messages that could flood HLT logbook
163 AliLog::SetClassDebugLevel("AliTRDrawStream", 0);
164 fRawReaderTrd->SetErrorDebugLevel(AliTRDrawStream::kLinkMonitor, 1);
165
166 fTrackingData = new AliTRDonlineTrackingDataContainer();
167 if (!fTrackingData){
168 iResult=-ENOMEM;
169 break;
170 }
171
172 } while (0);
173
174 if (iResult < 0) {
175
176 if (fTrackingData) delete fTrackingData;
177 fTrackingData = NULL;
178
179 if (fRawReaderTrd) delete fRawReaderTrd;
180 fRawReaderTrd = NULL;
181
182 if (fRawReaderMem) delete fRawReaderMem;
183 fRawReaderMem = NULL;
184
185 if (fDigitsManagerTrd) delete fDigitsManagerTrd;
186 fDigitsManagerTrd = NULL;
187
188 if (fGtuTrackArray) delete fGtuTrackArray;
189 fGtuTrackArray = NULL;
190
191 if (fTrackletArray) delete fTrackletArray;
192 fTrackletArray = NULL;
193
194 }
195
196 vector<const char*> remainingArgs;
197 for (int i = 0; i < argc; ++i)
198 remainingArgs.push_back(argv[i]);
199
200 if (argc > 0)
201 ConfigureFromArgumentString(remainingArgs.size(), &(remainingArgs[0]));
202
203 return iResult;
204}
205
206int AliHLTTRDPreprocessorComponent::DoDeinit() {
207
208 if (fTrackingData) delete fTrackingData;
209 fTrackingData = NULL;
210
211 if (fRawReaderTrd) delete fRawReaderTrd;
212 fRawReaderTrd = NULL;
213
214 if (fRawReaderMem) delete fRawReaderMem;
215 fRawReaderMem = NULL;
216
217 if (fGtuTrackArray) delete fGtuTrackArray;
218 fGtuTrackArray = NULL;
219
220 if (fTrackletArray) delete fTrackletArray;
221 fTrackletArray = NULL;
222
223 return 0;
224}
225
226//void AliHLTTRDPreprocessorComponent::DbgLog(const char* prefix, const char* msg){
227// AliHLTEventID_t eventNumber = fEventId;
228// int runNumber = -1;
229// printf("TRDGM %s-%s: [PRE] %s%s\n",
230// (runNumber >= 0) ? Form("%06d", runNumber) : "XXXXXX",
231// (eventNumber != fgkInvalidEventId) ? Form("%05llu", eventNumber) : "XXXXX",
232// (strlen(prefix) > 0) ? Form("<%s> ", prefix) : "", msg);
233//}
234
235
236void AliHLTTRDPreprocessorComponent::DbgLog(const char* prefix, ...){
237#ifdef __TRDHLTDEBUG
238 AliHLTEventID_t eventNumber = fEventId;
239 int runNumber = -1;
240 printf("TRDHLTGM %s-X-%s: [PRE] %s",
241 (runNumber >= 0) ? Form("%06d", runNumber) : "XXXXXX",
242 (eventNumber != fgkInvalidEventId) ? Form("%05llu", eventNumber) : "XXXXX",
243 (strlen(prefix) > 0) ? Form("<%s> ", prefix) : "");
244#endif
245 va_list args;
246 va_start(args, prefix);
247 char* fmt = va_arg(args, char*);
248 vprintf(fmt, args);
249 printf("\n");
250 va_end(args);
251}
252
253
254int AliHLTTRDPreprocessorComponent::DoEvent(const AliHLTComponentEventData& hltEventData,
255 AliHLTComponentTriggerData& /*trigData*/) {
256
257 fEventId = hltEventData.fEventID;
258 fTrackingData->SetLogPrefix(Form("TRDHLTGM XXXXXX-%05llu: [PRE] {TrkDat} ", hltEventData.fEventID));
259
260 LogDebug("### START DoEvent [event id: %llu, %d blocks, size: %d]",
261 hltEventData.fEventID, hltEventData.fBlockCnt, hltEventData.fStructSize);
262
263 // event processing function
264 int iResult = 0;
265 UInt_t sourceSectors = 0;
266
267 fTrackingData->Clear();
268 fTrackletArray->Clear();
269 fGtuTrackArray->Clear();
270
271 if (!IsDataEvent()) { // process data events only
272 LogDebug("### END DoEvent [event id: %llu, %d blocks, size: %d] (skipped: no data event)",
273 hltEventData.fEventID, hltEventData.fBlockCnt, hltEventData.fStructSize);
274 return iResult;
275 }
276
277 TString infoStr("");
278
279 // #FIXME: Also take care of SOR, EOR, etc...
280
281 // loop over all incoming TRD raw data blocks
282 for (const AliHLTComponentBlockData* pBlock = GetFirstInputBlock(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTRD);
283 pBlock != NULL && iResult >= 0;
284 pBlock = GetNextInputBlock()) {
285
286 int trdSector = -1;
287
288 // determine sector from block specification
289 for (unsigned pos = 0; pos < 8*sizeof(AliHLTUInt32_t); pos++) {
290 if (pBlock->fSpecification & (0x1 << pos)) {
291 if (trdSector >= 0) {
292 HLTWarning("Cannot uniquely identify DDL number from specification, skipping data block %s 0x%08x",
293 DataType2Text(pBlock->fDataType).c_str(),
294 pBlock->fSpecification);
295 trdSector = -1;
296 break;
297 }
298 trdSector = pos;
299 }
300 }
301 if (trdSector < 0) continue;
302
303 // add data block to rawreader
304 infoStr += Form("%02d, ", trdSector);
305 sourceSectors |= pBlock->fSpecification;
306 if(!fRawReaderMem->AddBuffer((UChar_t*) pBlock->fPtr, pBlock->fSize, trdSector + 1024)){
307 LogError("Could not add buffer of data block %s, 0x%08x to rawreader",
308 DataType2Text(pBlock->fDataType).c_str(),
309 pBlock->fSpecification);
310 continue;
311 }
312
313 } // loop over all incoming TRD raw data blocks
314
315 if (sourceSectors){
316 infoStr.Remove(infoStr.Length() - 2, 2);
317 LogDebug("preprocessing data from sectors: %s...", infoStr.Data());
318 } else {
319 LogDebug("### END DoEvent [event id: %llu, %d blocks, size: %d] (skipping: no TRD data)",
320 hltEventData.fEventID, hltEventData.fBlockCnt, hltEventData.fStructSize);
321 return iResult;
322 }
323
324 // extract header info, TRD tracklets and tracks from raw data
325 fRawReaderTrd->ReadEvent();
326
327 // read and store header info
328 for (UShort_t iSector = 0; iSector < 18; ++iSector){
329 if ((sourceSectors >> iSector) & 1){
330 fTrackingData->SetSectorTrgWord(iSector, fRawReaderTrd->GetTriggerFlags(iSector));
331 for (UShort_t iStack = 0; iStack < 5; ++iStack)
332 fTrackingData->SetStackTrgWord(iSector, iStack, fRawReaderTrd->GetTrkFlags(iSector, iStack));
333 }
334 }
335
336// //## dbg only!!!!
337// for (UShort_t iSector = 0; iSector < 18; ++iSector){
338// fTrackingData->SetSectorTrgWord(iSector, 0x2345352 ^ (iSector + 34));
339// for (UShort_t iStack = 0; iStack < 5; ++iStack){
340// ULong64_t dwl = ((ULong64_t)0xaffe << 16) | (ULong64_t)iSector << 8 | iStack;
341// ULong64_t dw = (((ULong64_t)0xbead << 16) | (ULong64_t)iSector << 8 | iStack) | ((ULong64_t)dwl << 32);
342// fTrackingData->SetStackTrgWord(iSector, iStack, dw);
343// }
344// }
345
346 // read and process TRD tracklets
347 Int_t trackletIndex[1080] = { 0 };
348 TList trklList;
349 Int_t iTrkl = 0;
350 trklList.SetOwner(kFALSE);
351 AliTRDrawStream::SortTracklets(fTrackletArray, trklList, trackletIndex);
352 TIter trackletIter(&trklList);
353 while (AliTRDtrackletBase* tracklet = (AliTRDtrackletBase*) trackletIter()) {
354
355 // tracklet data to tracking data for transport to other HLT components
356 fTrackingData->AddTracklet(tracklet->GetHCId(), tracklet->GetTrackletWord());
357
358// // conversion to AliESDTrdTracklet only for testing
359// // label -2, we assume it's always raw tracklets here
360// AliESDTrdTracklet esdTracklet(tracklet->GetTrackletWord(), tracklet->GetHCId(), -2);
361// Int_t det = esdTracklet.GetDetector();
362// if (kFALSE){
363// printf("TRDPreproc: TRD tracklet %3d - S%02d-%d-%d (det %3d): 0x%08x - y=%+5d dy=%+3d pid=%3d\n",
364// iTrkl, det/30, (det%30)/6, det%6, det,
365// esdTracklet.GetTrackletWord(), esdTracklet.GetBinY(), esdTracklet.GetBinDy(), esdTracklet.GetPID());
366// }
367
368 ++iTrkl;
369 }
370
371 // read and process GTU tracks
372 UInt_t numGtuTracks = fGtuTrackArray->GetEntriesFast();
373 AliESDTrdTrack *trdTrack = NULL;
374 UInt_t stack;
375 Int_t refIndex[6];
376
377 for (UInt_t iTrack = 0; iTrack < numGtuTracks; iTrack++) {
378 trdTrack = (AliESDTrdTrack*) ((*fGtuTrackArray)[iTrack]);
379 if (trdTrack){
380
381 AliTRDrawStream::AssignTracklets(trdTrack, trackletIndex, refIndex);
382 stack = trdTrack->GetStack();
383
384 for (Int_t iLayer = 0; iLayer < 6; ++iLayer) {
385 AliESDTrdTracklet *trkl = (refIndex[iLayer] > -1) ? (AliESDTrdTracklet*) trklList.At(refIndex[iLayer]) : 0x0;
386 if (trkl)
387 trdTrack->AddTrackletReference(trkl, iLayer);
388 }
389
390 fTrackingData->AddTrack(trdTrack);
391
392// for (Int_t iLayer = 0; iLayer < 6; ++iLayer) {
393// Int_t det = trdTrack->GetSector()*30 + stack*6 + iLayer;
394//
395// AliESDTrdTracklet *trkl = refIndex[iLayer] > -1 ? (AliESDTrdTracklet*) trklList.At(refIndex[iLayer]) : 0x0;
396// if (trkl) {
397// AliDebugClass(5, Form("adding tracklet with index %i: 0x%08x",
398// refIndex[iLayer], trkl->GetTrackletWord()));
399// if (trkl->GetDetector() != det)
400// AliErrorClass(Form("inconsistent assignment of tracklet 0x%08x in det %i to track in %i",
401// trkl->GetTrackletWord(), trkl->GetDetector(), det));
402// trdTrack->AddTrackletReference(trkl, iLayer);
403// }
404// }
405//
406// UInt_t pid = 0;
407// UShort_t lyrNum = 0;
408// TString trackletInfo("");
409// for (UShort_t iLayer = 0; iLayer < 6; ++iLayer){
410// AliESDTrdTracklet* trkl = trdTrack->GetTracklet(5 - iLayer);
411// if (trkl){
412// trackletInfo += Form("0x%08x ", trkl->GetTrackletWord());
413// pid += trkl->GetPID();
414// lyrNum++;
415// } else
416// trackletInfo += "---------- ";
417// }
418// trackletInfo.Remove(trackletInfo.Length() - 2, 2);
419// pid /= lyrNum;
420//
421// Int_t pidDiff = trdTrack->GetPID() - pid;
422//
423// UInt_t lm = trdTrack->GetLayerMask();
424// DbgLog("", Form("GTU track %d - S%02d-%d pt=%+.3f lm=(L5)%d%d%d%d%d%d(L0) [%s] [pid %s: %d|%d]",
425// iTrack, trdTrack->GetSector(), trdTrack->GetStack(), trdTrack->Pt(),
426// (lm >> 5) & 0x1, (lm >> 4) & 0x1, (lm >> 3) & 0x1,
427// (lm >> 2) & 0x1, (lm >> 1) & 0x1, (lm >> 0) & 0x1,
428// trackletInfo.Data(),
429// (TMath::Abs(pidDiff) <= 1) ? "ok" : "mismatch",
430// trdTrack->GetPID(), pid));
431//
432
433 } else
434 LogError("GTU track %d is NULL!\n", iTrack);
435 } // loop over GTU tracks
436
437 fRawReaderMem->ClearBuffers();
438
439 if (sourceSectors){
440 LogDebug("pushing data for sectors: 0x%05x", sourceSectors);
441 // transport of TRD tracklets and GTU tracks via tracking data container
442 void* buffer;
443 UInt_t bufferSize;
444 fTrackingData->PrintSummary("preproc component");
445 fTrackingData->Compress(buffer, bufferSize);
446 iResult += PushBack(buffer, bufferSize, AliHLTTRDDefinitions::fgkOnlineDataType, sourceSectors);
447 free(buffer);
448 }
449
450 LogDebug("### END DoEvent [event id: %llu, %d blocks, size: %d]",
451 hltEventData.fEventID, hltEventData.fBlockCnt, hltEventData.fStructSize);
452
453 return iResult;
454}