Line data Source code
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 : //* Primary Authors: Felix Rettig, Stefan Kirsch *
7 : //* for The ALICE HLT Project. *
8 : //* *
9 : //* Permission to use, copy, modify and distribute this software and its *
10 : //* documentation strictly for non-commercial purposes is hereby granted *
11 : //* without fee, provided that the above copyright notice appears in all *
12 : //* copies and that both the copyright notice and this permission notice *
13 : //* appear in the supporting documentation. The authors make no claims *
14 : //* about the suitability of this software for any purpose. It is *
15 : //* provided "as is" without express or implied warranty. *
16 : //**************************************************************************
17 :
18 : /// @file AliHLTTRDMonitorComponent.cxx
19 : /// @author Felix Rettig, Stefan Kirsch
20 : /// @date 2012-08-16
21 : /// @brief The TRD monitoring component
22 : ///
23 :
24 : #include <cstdlib>
25 : #include "TH1I.h"
26 : #include "TH2I.h"
27 : #include "TH2F.h"
28 : #include "TObjString.h"
29 : #include "TObjArray.h"
30 : #include "TClonesArray.h"
31 : #include "TFile.h"
32 : #include "AliESDtrack.h"
33 : #include "AliHLTDataTypes.h"
34 : #include "AliHLTTRDDefinitions.h"
35 : #include "AliHLTLogging.h"
36 : #include "AliHLTTRDMonitorComponent.h"
37 : #include "AliTRDonlineTrackingDataContainer.h"
38 :
39 : #define TRDMODULES 18
40 :
41 6 : ClassImp(AliHLTTRDMonitorComponent)
42 :
43 : #define LogError( ... ) { HLTError(__VA_ARGS__); if (fDebugLevel >= 1) { DbgLog("ERROR", __VA_ARGS__); } }
44 : #define LogInfo( ... ) { HLTInfo(__VA_ARGS__); if (fDebugLevel >= 1) { DbgLog("INFO", __VA_ARGS__); } }
45 : #define LogInspect( ... ) { HLTDebug(__VA_ARGS__); if (fDebugLevel >= 1) { DbgLog("INSPECT", __VA_ARGS__); } }
46 : #define LogDebug( ... ) { if (fDebugLevel >= 1) { HLTInfo(__VA_ARGS__); DbgLog("DEBUG", __VA_ARGS__); } }
47 :
48 3 : AliHLTTRDMonitorComponent::AliHLTTRDMonitorComponent() : AliHLTProcessor(),
49 3 : fTrackHighPtThreshold(2.3),
50 3 : fHistoMode(1),
51 3 : fTrackingDataDebugOutput(kFALSE),
52 3 : fDebugLevel(0),
53 3 : fWriteHistos(kFALSE),
54 3 : fEventId(fgkInvalidEventId),
55 3 : fTrackingData(NULL),
56 3 : fHistArray(NULL),
57 3 : fHistEventTypes(NULL),
58 3 : fHistTrackletY(NULL),
59 3 : fHistTrackletDy(NULL),
60 3 : fHistTrackletYDy(NULL),
61 3 : fHistTrackletZ(NULL),
62 3 : fHistTrackletPID(NULL),
63 3 : fHistTrackletsHCId(NULL),
64 3 : fHistTrackPt(NULL),
65 3 : fHistTrackPID(NULL),
66 3 : fHistTrackLayers(NULL),
67 3 : fHistTrackLayersHighPt(NULL),
68 3 : fHistTracksStack(NULL),
69 3 : fHistTrackletTimingStack(NULL),
70 3 : fHistTrackingTiming(NULL),
71 3 : fHistTriggerContribs(NULL)
72 15 : {
73 6 : }
74 :
75 12 : AliHLTTRDMonitorComponent::~AliHLTTRDMonitorComponent() {
76 12 : }
77 :
78 : const char* AliHLTTRDMonitorComponent::GetComponentID() {
79 90 : return "TRDMonitorComponent";
80 : }
81 :
82 : void AliHLTTRDMonitorComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
83 0 : list.push_back(kAliHLTDataTypeTObject | kAliHLTDataOriginTRD);
84 0 : }
85 :
86 : AliHLTComponentDataType AliHLTTRDMonitorComponent::GetOutputDataType() {
87 0 : return (kAliHLTDataTypeTObjArray | kAliHLTDataOriginTRD);
88 : }
89 :
90 : void AliHLTTRDMonitorComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
91 0 : constBase = 10000000;
92 0 : inputMultiplier = 0;
93 0 : }
94 :
95 : AliHLTComponent* AliHLTTRDMonitorComponent::Spawn() {
96 0 : return new AliHLTTRDMonitorComponent;
97 0 : }
98 :
99 : int AliHLTTRDMonitorComponent::Configure(const char* /*arguments*/) {
100 0 : return 0;
101 : }
102 :
103 : int AliHLTTRDMonitorComponent::Reconfigure(const char* /*cdbEntry*/, const char* /*chainId*/) {
104 0 : return 0;
105 : }
106 :
107 : int AliHLTTRDMonitorComponent::DoInit(int argc, const char** argv) {
108 :
109 : int iResult = 0;
110 :
111 : do {
112 :
113 0 : fTrackingData = new AliTRDonlineTrackingDataContainer();
114 0 : if (!fTrackingData) {
115 : iResult = -ENOMEM;
116 0 : break;
117 : }
118 0 : fTrackingData->SetGtuPtMultiplier(-1.); // this component does not know about the B-field direction
119 :
120 0 : fHistArray = new TObjArray(25);
121 0 : if(!fHistArray)
122 0 : return -ENOMEM;
123 0 : fHistArray->SetOwner(kTRUE);
124 :
125 0 : } while (0);
126 :
127 0 : if (iResult < 0){
128 :
129 0 : if (fHistArray) delete fHistArray;
130 0 : fHistArray = NULL;
131 :
132 0 : if (fTrackingData) delete fTrackingData;
133 0 : fTrackingData = NULL;
134 :
135 0 : }
136 :
137 0 : fHistEventTypes = new TH1I("trdmon_event_types", "event types analysed;event type;abundance", 10, 0, 10);
138 0 : fHistEventTypes->GetXaxis()->SetBinLabel(1, "any");
139 0 : fHistEventTypes->GetXaxis()->SetBinLabel(2, "data");
140 0 : fHistEventTypes->GetXaxis()->SetBinLabel(3, "track(lets) present");
141 0 : fHistArray->AddLast(fHistEventTypes);
142 :
143 0 : fHistTrackletY = new TH1I("trdmon_tracklet_y", "tracklet y-position;y (cm);abundance", 125, -75, 75);
144 0 : fHistArray->AddLast(fHistTrackletY);
145 :
146 0 : fHistTrackletDy = new TH1I("trdmon_tracklet_dy", "tracklet deflection;#Delta y (140 #mu m);abundance", 140, -70, 70);
147 0 : fHistArray->AddLast(fHistTrackletDy);
148 :
149 0 : fHistTrackletYDy = new TH2I("trdmon_tracklet_y_dy", "tracklet y-deflection vs. y-position;y (160 #mu m);#Delta y (140 #mu m);", 256, -4096, 4096, 140, -70, 70);
150 0 : fHistArray->AddLast(fHistTrackletYDy);
151 :
152 0 : fHistTrackletZ = new TH1I("trdmon_tracklet_z", "tracklet z-position;z (padrow);abundance", 16, 0, 16);
153 0 : fHistArray->AddLast(fHistTrackletZ);
154 :
155 0 : fHistTrackletPID = new TH1I("trdmon_tracklet_pid", "tracklet PID;PID (a.u.);abundance", 256, 0, 256);
156 0 : fHistArray->AddLast(fHistTrackletPID);
157 :
158 0 : fHistTrackletsHCId = new TH2F("trdmon_tracklets_hc", "tracklet number by HC;TRD sector;TRD half-chamber",
159 : 18, 0, 18, 60, 0, 60);
160 0 : fHistArray->AddLast(fHistTrackletsHCId);
161 :
162 0 : fHistTrackletTimingStack = new TH2I("trdmon_tracklet_timing_stack", "tracklet timing;TRD stack; time after L0 (us)",
163 : 90, 0., 90, 270, 0., 9.);
164 0 : fHistArray->AddLast(fHistTrackletTimingStack);
165 :
166 0 : fHistTrackPt = new TH1I("trdmon_track_pt", "p_{T} of TRD online tracks;p_{T}^{TRD, online} (GeV/c);abundance", 200, -20, 20.);
167 0 : fHistArray->AddLast(fHistTrackPt);
168 :
169 0 : fHistTrackPID = new TH1I("trdmon_track_pid", "PID of TRD online tracks;PID^{TRD, online} (a.u.);abundance", 256, 0, 256);
170 0 : fHistArray->AddLast(fHistTrackPID);
171 :
172 0 : fHistTrackLayers = new TH1I("trdmon_track_layers", "contributing layers to TRD online tracks;contributing layers;abundance", 7, 0, 7);
173 0 : fHistArray->AddLast(fHistTrackLayers);
174 :
175 0 : fHistTrackLayersHighPt = new TH1I("trdmon_track_layers_hpt", "contributing layers to TRD online tracks;contributing layers;abundance", 7, 0, 7);
176 0 : fHistArray->AddLast(fHistTrackLayersHighPt);
177 :
178 0 : fHistTracksStack = new TH2F("trdmon_tracks_stack", "tracks by stack;TRD sector;TRD stack",
179 : 18, 0, 18, 5, 0, 5);
180 0 : fHistArray->AddLast(fHistTracksStack);
181 :
182 : const Double_t trackingTimesTimeBin = 0.025; // us
183 : const Double_t trackingTimesMaxTime = 12.; // us
184 0 : fHistTrackingTiming = new TH2I("trdmon_tracking_timing", "tracking timing;;time after L0 (#mu s);",
185 : 4, 0, 4, (Int_t)(trackingTimesMaxTime/trackingTimesTimeBin), 0., trackingTimesMaxTime);
186 0 : fHistTrackingTiming->GetXaxis()->SetBinLabel(1, "tracklet start");
187 0 : fHistTrackingTiming->GetXaxis()->SetBinLabel(2, "tracklet end");
188 0 : fHistTrackingTiming->GetXaxis()->SetBinLabel(3, "stack done");
189 0 : fHistTrackingTiming->GetXaxis()->SetBinLabel(4, "sector done");
190 0 : fHistArray->AddLast(fHistTrackingTiming);
191 :
192 0 : fHistTriggerContribs = new TH2I("trdmon_trigger_contribs", "TRD internal contributions by sector;TRD sector;trigger contribution;",
193 : 18, 0, 18, 12, 0, 12);
194 0 : fHistTrackingTiming->GetYaxis()->SetBinLabel(1, "trg0");
195 0 : fHistTrackingTiming->GetYaxis()->SetBinLabel(2, "trg1");
196 0 : fHistTrackingTiming->GetYaxis()->SetBinLabel(3, "trg2");
197 0 : fHistTrackingTiming->GetYaxis()->SetBinLabel(4, "trg3");
198 0 : fHistTrackingTiming->GetYaxis()->SetBinLabel(5, "trg4");
199 0 : fHistTrackingTiming->GetYaxis()->SetBinLabel(6, "trg5");
200 0 : fHistTrackingTiming->GetYaxis()->SetBinLabel(7, "trg5");
201 0 : fHistTrackingTiming->GetYaxis()->SetBinLabel(8, "T");
202 0 : fHistArray->AddLast(fHistTriggerContribs);
203 :
204 0 : vector<const char*> remainingArgs;
205 0 : for (int i = 0; i < argc; ++i)
206 0 : remainingArgs.push_back(argv[i]);
207 :
208 0 : if (argc > 0)
209 0 : ConfigureFromArgumentString(remainingArgs.size(), &(remainingArgs[0]));
210 :
211 : return 0;
212 0 : }
213 :
214 : int AliHLTTRDMonitorComponent::DoDeinit() {
215 :
216 0 : if ((fHistoMode == 1) && (fWriteHistos)){
217 0 : TFile out("mon_out/mon_hists.root", "RECREATE");
218 0 : if (!out.IsZombie()) {
219 0 : out.cd();
220 0 : UInt_t numHists = fHistArray->GetEntries();
221 0 : for (UInt_t iHist = 0; iHist < numHists; ++iHist)
222 0 : if (fHistArray->At(iHist))
223 0 : fHistArray->At(iHist)->Write();
224 0 : out.Close();
225 0 : }
226 0 : }
227 :
228 0 : if (fHistArray) delete fHistArray;
229 0 : fHistArray = NULL;
230 :
231 0 : if (fTrackingData) delete fTrackingData;
232 0 : fTrackingData = NULL;
233 :
234 0 : return 0;
235 0 : }
236 :
237 : int AliHLTTRDMonitorComponent::ScanConfigurationArgument(int argc, const char** argv)
238 : {
239 :
240 0 : if (argc <= 0)
241 0 : return 0;
242 :
243 : UShort_t iArg = 0;
244 0 : TString argument(argv[iArg]);
245 :
246 0 : if (!argument.CompareTo("-write-histograms")){
247 0 : LogInfo("writing of histograms enabled.");
248 0 : fWriteHistos = kTRUE; // enable histogram writing, for debugging/tuning only!
249 0 : return 1;
250 : }
251 :
252 0 : if (!argument.CompareTo("-debug")){
253 0 : if (++iArg >= argc) return -EPROTO;
254 0 : argument = argv[iArg];
255 0 : fDebugLevel = argument.Atoi();
256 0 : LogInfo("debug level set to %d.", fDebugLevel);
257 0 : return 2;
258 : }
259 :
260 0 : return 0;
261 :
262 0 : }
263 :
264 : void AliHLTTRDMonitorComponent::DbgLog(const char* prefix, ...){
265 : #ifdef __TRDHLTDEBUG
266 : AliHLTEventID_t eventNumber = fEventId;
267 : Int_t fRunNumber = -1;
268 : printf("TRDHLTGM %s-X-%s: [MON] %s",
269 : (fRunNumber >= 0) ? Form("%06d", fRunNumber) : "XXXXXX",
270 : (eventNumber != fgkInvalidEventId) ? Form("%05llu", eventNumber) : "XXXXX",
271 : (strlen(prefix) > 0) ? Form("<%s> ", prefix) : "");
272 : #endif
273 0 : va_list args;
274 0 : va_start(args, prefix);
275 0 : char* fmt = va_arg(args, char*);
276 0 : vprintf(fmt, args);
277 0 : printf("\n");
278 0 : va_end(args);
279 0 : }
280 :
281 : int AliHLTTRDMonitorComponent::PrepareTRDData() {
282 :
283 : int result = 1;
284 :
285 0 : fTrackingData->Clear();
286 0 : for (const AliHLTComponentBlockData* datablock = GetFirstInputBlock(AliHLTTRDDefinitions::fgkOnlineDataType);
287 0 : datablock != NULL;
288 0 : datablock = GetNextInputBlock())
289 : {
290 0 : fTrackingData->Decompress(datablock->fPtr, datablock->fSize, kTRUE);
291 : }
292 :
293 0 : fTrackingData->PrintSummary("monitor component");
294 :
295 0 : return result;
296 :
297 : }
298 :
299 : void AliHLTTRDMonitorComponent::DumpTrackingData(){
300 0 : TString trklStr("");
301 0 : TString matchStr("");
302 : UShort_t layerMask;
303 :
304 0 : if (fTrackingData->GetNumTracklets() + fTrackingData->GetNumTracks() == 0)
305 0 : return;
306 :
307 0 : for (UShort_t iStack = 0; iStack < fkTRDStacks; ++iStack){
308 0 : for (Int_t iTrk = 0; iTrk < fTrackingData->GetNumTracks(iStack); ++iTrk){
309 :
310 0 : layerMask = fTrackingData->GetTrackLayerMask(iStack, iTrk);
311 0 : trklStr = Form("trkl: ");
312 0 : for (Short_t iLayer = 5; iLayer >= 0; --iLayer){
313 0 : if ((layerMask >> iLayer) & 1)
314 0 : trklStr += Form("0x%08x (%+8.3f) ",
315 0 : fTrackingData->GetTrackTrackletWord(iStack, iTrk, iLayer),
316 0 : fTrackingData->GetTrackTrackletLocalY(iStack, iTrk, iLayer));
317 : else
318 0 : trklStr += "--------------------- ";
319 : } // loop over layers
320 0 : trklStr.Remove(trklStr.Length() - 2, 2);
321 :
322 0 : if (fTrackingDataDebugOutput){
323 :
324 0 : printf("###DOTDB EV%04llu GTU TRACK - S%02d-%d pt: %+7.2f pid: %3d lm: 0x%02x %s\n",
325 0 : fEventId,
326 0 : iStack/5, iStack%5, fTrackingData->GetTrackPt(iStack, iTrk), fTrackingData->GetTrackPID(iStack, iTrk),
327 0 : layerMask, trklStr.Data());
328 : }
329 :
330 : // paranoia checks
331 0 : for (Short_t iLayer = 5; iLayer >= 0; --iLayer){
332 0 : if (((layerMask >> iLayer) & 1) && (fTrackingData->GetTrackTrackletWord(iStack, iTrk, iLayer) == 0x10001000))
333 0 : LogError("invalid layer mask / tracklet value combination A", "");
334 :
335 0 : if ((((layerMask >> iLayer) & 1) == 0) && (fTrackingData->GetTrackTrackletWord(iStack, iTrk, iLayer) != 0x10001000))
336 0 : LogError("invalid layer mask / tracklet value combination B", "");
337 : }
338 :
339 : } // loop over tracks in stack
340 : } // loop over stacks
341 :
342 0 : }
343 :
344 :
345 : int AliHLTTRDMonitorComponent::ProcessTRDData(){
346 :
347 : UInt_t numTracklets;
348 : UInt_t numTracks;
349 :
350 0 : if (fHistoMode == 0){
351 0 : UInt_t numHists = fHistArray->GetEntries();
352 0 : for (UInt_t iHist = 0; iHist < numHists; ++iHist)
353 0 : if (fHistArray->At(iHist))
354 0 : ((TH1*) (fHistArray->At(iHist)))->Reset();
355 0 : }
356 :
357 : // tracklets
358 0 : for (UInt_t iDet = 0; iDet < fkTRDChambers; ++iDet){
359 0 : numTracklets = fTrackingData->GetNumTracklets(iDet);
360 0 : for (UInt_t iTrkl = 0; iTrkl < numTracklets; ++iTrkl){
361 0 : fHistTrackletY->Fill(fTrackingData->GetTrackletLocalY(iDet, iTrkl));
362 0 : fHistTrackletDy->Fill(fTrackingData->GetTrackletBinDy(iDet, iTrkl));
363 0 : fHistTrackletYDy->Fill(fTrackingData->GetTrackletBinY(iDet, iTrkl), fTrackingData->GetTrackletBinDy(iDet, iTrkl));
364 0 : fHistTrackletZ->Fill(fTrackingData->GetTrackletBinZ(iDet, iTrkl));
365 0 : fHistTrackletPID->Fill(fTrackingData->GetTrackletPID(iDet, iTrkl));
366 0 : Int_t hc = fTrackingData->GetTrackletHCId(iDet, iTrkl);
367 0 : fHistTrackletsHCId->Fill(hc/60, hc%60);
368 : }
369 : }
370 :
371 0 : for (UShort_t iStack = 0; iStack < fkTRDStacks; ++iStack){
372 : // timing
373 0 : fHistTrackletTimingStack->Fill(iStack, fTrackingData->GetTrackletStartTime(iStack/5, iStack%5));
374 0 : fHistTrackletTimingStack->Fill(iStack, fTrackingData->GetTrackletEndTime(iStack/5, iStack%5));
375 0 : fHistTrackingTiming->Fill(0., fTrackingData->GetTrackletStartTime(iStack/5, iStack%5));
376 0 : fHistTrackingTiming->Fill(1., fTrackingData->GetTrackletEndTime(iStack/5, iStack%5));
377 0 : fHistTrackingTiming->Fill(2., fTrackingData->GetTMUTrackingDoneTime(iStack/5, iStack%5));
378 :
379 : // GTU tracks
380 0 : numTracks = fTrackingData->GetNumTracks(iStack);
381 0 : fHistTracksStack->Fill(iStack/5, iStack%5, numTracks);
382 0 : for (UInt_t iTrk = 0; iTrk < numTracks; ++iTrk){
383 0 : Double_t gpt = fTrackingData->GetTrackPt(iStack, iTrk);
384 0 : fHistTrackPt->Fill(gpt);
385 0 : fHistTrackPID->Fill(fTrackingData->GetTrackPID(iStack, iTrk));
386 0 : fHistTrackLayers->Fill(fTrackingData->GetTrackLayerNum(iStack, iTrk));
387 0 : if (gpt >= fTrackHighPtThreshold)
388 0 : fHistTrackLayersHighPt->Fill(fTrackingData->GetTrackLayerNum(iStack, iTrk));
389 : } // loop over tracks in stack
390 : } // loop over stacks
391 :
392 0 : for (UShort_t iSector = 0; iSector < fkTRDSectors; ++iSector){
393 0 : fHistTrackingTiming->Fill(3., fTrackingData->GetSMUTrackingDoneTime(iSector), fkTRDStacksPerSector);
394 :
395 0 : UInt_t sectorTrgFlags = fTrackingData->GetSectorTrgContribs(iSector);
396 0 : for (UShort_t iTrgCtb = 0; iTrgCtb < 12; ++iTrgCtb)
397 0 : if ((sectorTrgFlags >> iTrgCtb) & 1)
398 0 : fHistTriggerContribs->Fill(iSector, iTrgCtb);
399 : }
400 :
401 0 : return kTRUE;
402 : }
403 :
404 : int AliHLTTRDMonitorComponent::DoEvent(const AliHLTComponentEventData& evtData,
405 : AliHLTComponentTriggerData& /*trigData*/) {
406 :
407 : int iResult = 0;
408 0 : fEventId = evtData.fEventID;
409 0 : fHistEventTypes->Fill(0.);
410 :
411 0 : LogDebug("### START DoEvent [event id: %llu, %d blocks, size: %d]",
412 : fEventId, evtData.fBlockCnt, evtData.fStructSize);
413 :
414 0 : if (!IsDataEvent()) { // process data events only
415 0 : LogDebug("### END DoEvent [event id: %llu, %d blocks, size: %d] (skipped: no data event)",
416 : fEventId, evtData.fBlockCnt, evtData.fStructSize);
417 0 : return iResult;
418 : }
419 :
420 0 : fHistEventTypes->Fill(1.);
421 :
422 0 : fTrackingData->SetLogPrefix(Form("TRDHLTGM XXXXXX-%05llu: [MON] {TrkDat} ", fEventId));
423 :
424 : do {
425 :
426 : // access to TRD specific data from AliHLTTRDPreprocessorComponent
427 0 : if (!PrepareTRDData()){
428 0 : LogError("access to TRD data failed. Skipping event...", "");
429 : break;
430 : }
431 :
432 0 : if (fTrackingData->GetNumTracks() + fTrackingData->GetNumTracklets() == 0) {
433 0 : LogDebug("no TRD-relevant information, skipping further event processing");
434 : break;
435 : }
436 :
437 0 : fHistEventTypes->Fill(2.);
438 :
439 : // DumpTrackingData();
440 :
441 0 : if (!ProcessTRDData()) {
442 0 : LogError("processing of TRD data failed, skipping further event processing");
443 : break;
444 : }
445 :
446 : break;
447 :
448 : } while (1);
449 :
450 0 : PushBack(fHistArray, (kAliHLTDataTypeTObjArray | kAliHLTDataOriginTRD), 0x3fffff);
451 :
452 0 : LogDebug("### END DoEvent [event id: %llu, %d blocks, size: %d]",
453 : fEventId, evtData.fBlockCnt, evtData.fStructSize);
454 :
455 0 : return iResult;
456 0 : }
|