]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDrawStream.cxx
- propagate track in-time flag to AliESDTrdTrack
[u/mrichter/AliRoot.git] / TRD / AliTRDrawStream.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 ////////////////////////////////////////////////////////////////////////////
17 //                                                                        //
18 //  Decoding data from the TRD raw stream                                 //
19 //  and translation into ADC values, on-line tracklets and tracks         //
20 //                                                                        //
21 //  Author: J. Klein (jochen.klein@cern.ch)                               //
22 //                                                                        //
23 ////////////////////////////////////////////////////////////////////////////
24
25 #include <cstdio>
26 #include <cstdarg>
27
28 #include "TClonesArray.h"
29 #include "TTree.h"
30
31 #include "AliLog.h"
32 #include "AliRawReader.h"
33 #include "AliTRDcalibDB.h"
34 #include "AliTRDdigitsManager.h"
35 #include "AliTRDdigitsParam.h"
36 #include "AliTRDcalibDB.h"
37 #include "AliTRDmcmSim.h"
38 #include "AliTRDtrapConfig.h"
39 #include "AliTRDarrayADC.h"
40 #include "AliTRDarrayDictionary.h"
41 #include "AliTRDSignalIndex.h"
42 #include "AliTRDtrackletWord.h"
43 #include "AliTRDtrackletMCM.h"
44 #include "AliESDTrdTrack.h"
45 #include "AliTreeLoader.h"
46 #include "AliLoader.h"
47
48 #include "AliTRDrawStream.h"
49
50 // temporary
51 #include "AliRunLoader.h"
52
53 ClassImp(AliTRDrawStream)
54
55 // some static information
56 Int_t AliTRDrawStream::fgMcmOrder[] = {12, 13, 14, 15,
57                                        8, 9, 10, 11,
58                                        4, 5, 6, 7,
59                                        0, 1, 2, 3};
60 Int_t  AliTRDrawStream::fgRobOrder [] = {0, 1, 2, 3};
61 const Int_t  AliTRDrawStream::fgkNlinks = 12;
62 const Int_t  AliTRDrawStream::fgkNstacks = 5;
63 const Int_t  AliTRDrawStream::fgkNsectors = 18;
64 const Int_t  AliTRDrawStream::fgkNtriggers = 12;
65 const UInt_t AliTRDrawStream::fgkDataEndmarker     = 0x00000000;
66 const UInt_t AliTRDrawStream::fgkTrackletEndmarker = 0x10001000;
67 const UInt_t AliTRDrawStream::fgkStackEndmarker[] = { 0xe0d01000, 0xe0d10000 };
68
69 const char* AliTRDrawStream::fgkErrorMessages[] = {
70   "Unknown error",
71   "Link monitor active",
72   "Event counter mismatch",
73   "not a TRD equipment (1024-1041)",
74   "Invalid Stack header",
75   "Invalid detector number",
76   "No digits could be retrieved from the digitsmanager",
77   "HC header mismatch",
78   "HC check bits wrong",
79   "Unexpected position in readout stream",
80   "Invalid testpattern mode",
81   "Testpattern mismatch",
82   "Number of timebins changed",
83   "ADC mask inconsistent",
84   "ADC check bits invalid",
85   "Missing ADC data",
86   "Missing expected ADC channels",
87   "Missing MCM headers",
88   "Missing TP data",
89   "CRC mismatch"
90 };
91
92 Int_t AliTRDrawStream::fgErrorDebugLevel[] = {
93   0,
94   0,
95   2,
96   1,
97   0,
98   1,
99   1,
100   0,
101   1,
102   2,
103   1,
104   0,
105   1,
106   1,
107   2,
108   1,
109   1,
110   1,
111   0,
112   0
113 };
114
115 AliTRDrawStream::ErrorBehav_t AliTRDrawStream::fgErrorBehav[] = {
116   AliTRDrawStream::kTolerate,
117   AliTRDrawStream::kDiscardHC,
118   AliTRDrawStream::kTolerate,
119   AliTRDrawStream::kAbort,
120   AliTRDrawStream::kAbort,
121   AliTRDrawStream::kAbort,
122   AliTRDrawStream::kAbort,
123   AliTRDrawStream::kDiscardHC,
124   AliTRDrawStream::kDiscardHC,
125   AliTRDrawStream::kTolerate,
126   AliTRDrawStream::kTolerate,
127   AliTRDrawStream::kTolerate,
128   AliTRDrawStream::kTolerate,
129   AliTRDrawStream::kTolerate,
130   AliTRDrawStream::kTolerate,
131   AliTRDrawStream::kTolerate,
132   AliTRDrawStream::kTolerate,
133   AliTRDrawStream::kTolerate,
134   AliTRDrawStream::kTolerate,
135   AliTRDrawStream::kTolerate
136 };
137
138 AliTRDrawStream::AliTRDrawStream(AliRawReader *rawReader) :
139   fStoreError(&AliTRDrawStream::ForgetError),
140   fRawReader(rawReader),
141   fDigitsManager(0x0),
142   fDigitsParam(0x0),
143   fErrors(0x0),
144   fLastError(),
145   fErrorFlags(0),
146   fStats(),
147   fPayloadStart(0x0),
148   fPayloadCurr(0x0),
149   fPayloadSize(0),
150   fNtimebins(-1),
151   fLastEvId(-1),
152   fCurrSlot(-1),
153   fCurrLink(-1),
154   fCurrRobPos(-1),
155   fCurrMcmPos(-1),
156   fCurrEquipmentId(0),
157   fCurrSmHeaderSize(0),
158   fCurrSmHeaderVersion(0),
159   fCurrTrailerReadout(0),
160   fCurrTrgHeaderAvail(0),
161   fCurrTrgHeaderReadout(0),
162   fCurrTrkHeaderAvail(0),
163   fCurrStackEndmarkerAvail(0),
164   fCurrEvType(0),
165   fCurrTriggerEnable(0),
166   fCurrTriggerFired(0),
167   fCurrTrackEnable(0),
168   fCurrTrackletEnable(0),
169   fCurrStackMask(0),
170   fCurrTrkHeaderIndexWord(0x0),
171   fCurrTrkHeaderSize(0x0),
172   fCurrTrkFlags(0x0),
173   fCurrTrgHeaderIndexWord(0x0),
174   fCurrTrgHeaderSize(0x0),
175   fCurrTrgFlags(0x0),
176   fCurrStackIndexWord(0x0),
177   fCurrStackHeaderSize(0x0),
178   fCurrStackHeaderVersion(0x0),
179   fCurrLinkMask(0x0),
180   fCurrCleanCheckout(0x0),
181   fCurrBoardId(0x0),
182   fCurrHwRev(-1),
183   fCurrHwRevTMU(0x0),
184   fCurrLinkMonitorFlags(0x0),
185   fCurrLinkDataTypeFlags(0x0),
186   fCurrLinkDebugFlags(0x0),
187   fCurrMatchFlagsSRAM(0),
188   fCurrMatchFlagsPostBP(0),
189   fCurrChecksumStack(),
190   fCurrChecksumSIU(0),
191   fCurrSpecial(-1),
192   fCurrMajor(-1),
193   fCurrMinor(-1),
194   fCurrAddHcWords(-1),
195   fCurrSm(-1),
196   fCurrStack(-1),
197   fCurrLayer(-1),
198   fCurrSide(-1),
199   fCurrHC(-1),
200   fCurrCheck(-1),
201   fCurrNtimebins(-1),
202   fCurrBC(-1),
203   fCurrPtrgCnt(-1),
204   fCurrPtrgPhase(-1),
205   fNDumpMCMs(0),
206   fTrackletArray(0x0),
207   fAdcArray(0x0),
208   fSignalIndex(0x0),
209   fTrackletTree(0x0),
210   fTracklets(0x0),
211   fTracks(0x0),
212   fMarkers(0x0)
213 {
214   // default constructor
215
216   fCurrTrkHeaderIndexWord = new UInt_t[fgkNstacks];
217   fCurrTrkHeaderSize      = new UInt_t[fgkNstacks];
218   fCurrTrkFlags           = new ULong64_t[fgkNsectors*fgkNstacks];
219   fCurrTrgHeaderIndexWord = new UInt_t[fgkNtriggers];
220   fCurrTrgHeaderSize      = new UInt_t[fgkNtriggers];
221   fCurrTrgFlags           = new UInt_t[fgkNsectors];
222   fCurrStackIndexWord     = new UInt_t[fgkNstacks];
223   fCurrStackHeaderSize    = new UInt_t[fgkNstacks];
224   fCurrStackHeaderVersion = new UInt_t[fgkNstacks];
225   fCurrLinkMask           = new UInt_t[fgkNstacks];
226   fCurrCleanCheckout      = new UInt_t[fgkNstacks];
227   fCurrBoardId            = new UInt_t[fgkNstacks];
228   fCurrHwRevTMU           = new UInt_t[fgkNstacks];
229   fCurrLinkMonitorFlags   = new UInt_t[fgkNstacks * fgkNlinks];
230   fCurrLinkDataTypeFlags  = new UInt_t[fgkNstacks * fgkNlinks];
231   fCurrLinkDebugFlags     = new UInt_t[fgkNstacks * fgkNlinks];
232   for (Int_t iSector = 0; iSector < fgkNsectors; iSector++)
233     fCurrTrgFlags[iSector] = 0;
234   for (Int_t i = 0; i < 100; i++)
235     fDumpMCM[i] = 0;
236
237   // preparing TClonesArray
238   fTrackletArray = new TClonesArray("AliTRDtrackletWord", 256);
239
240   // setting up the error tree
241   fErrors = new TTree("errorStats", "Error statistics");
242   fErrors->SetDirectory(0x0);
243   fErrors->Branch("error", &fLastError);
244   fErrors->SetCircular(1000);
245   for (Int_t i = 0; i < 100; i++) {
246     fErrorBuffer[i] = 0;
247   }
248
249 }
250
251 AliTRDrawStream::~AliTRDrawStream()
252 {
253   // destructor
254
255   delete fErrors;
256
257   delete [] fCurrTrkHeaderIndexWord;
258   delete [] fCurrTrkHeaderSize;
259   delete [] fCurrTrkFlags;
260   delete [] fCurrTrgHeaderIndexWord;
261   delete [] fCurrTrgHeaderSize;
262   delete [] fCurrTrgFlags;
263   delete [] fCurrStackIndexWord;
264   delete [] fCurrStackHeaderSize;
265   delete [] fCurrStackHeaderVersion;
266   delete [] fCurrLinkMask;
267   delete [] fCurrCleanCheckout;
268   delete [] fCurrBoardId;
269   delete [] fCurrHwRevTMU;
270   delete [] fCurrLinkMonitorFlags;
271   delete [] fCurrLinkDataTypeFlags;
272   delete [] fCurrLinkDebugFlags;
273 }
274
275 Bool_t AliTRDrawStream::ReadEvent(TTree *trackletTree)
276 {
277   // read the current event from the raw reader and fill it to the digits manager
278
279   if (!fRawReader) {
280     AliError("No raw reader available");
281     return kFALSE;
282   }
283
284   // tracklet output
285   ConnectTracklets(trackletTree);
286
287   // some preparations
288   fDigitsParam = 0x0;
289
290   // loop over all DDLs
291   // data starts with GTU payload, i.e. SM index word
292   UChar_t *buffer = 0x0;
293
294   while (fRawReader->ReadNextData(buffer)) {
295
296     fCurrEquipmentId = fRawReader->GetEquipmentId();
297     AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
298
299     if (fCurrEquipmentId < kDDLOffset || fCurrEquipmentId > kDDLMax) {
300       EquipmentError(kNonTrdEq, "Skipping");
301       continue;
302     }
303
304     if (fMarkers)
305       new ((*fMarkers)[fMarkers->GetEntriesFast()])
306         AliTRDrawStreamError(-kSecactive, fCurrEquipmentId - kDDLOffset);
307
308     ReadGTUHeaders((UInt_t*) buffer);
309
310     if (fCurrTrailerReadout)
311       ReadGTUTrailer();
312
313     // loop over all active links
314     AliDebug(2, Form("Stack mask 0x%02x", fCurrStackMask));
315     for (Int_t iStack = 0; iStack < fgkNstacks; iStack++) {
316       fCurrSlot = iStack;
317       if ((fCurrStackMask & (1 << fCurrSlot)) == 0)
318         continue;
319
320       AliDebug(2, Form("Stack %i, Link mask: 0x%02x", fCurrSlot, fCurrLinkMask[fCurrSlot]));
321       for (Int_t iLink = 0; iLink < fgkNlinks; iLink++) {
322         fCurrLink = iLink;
323         fCurrHC   = (fCurrEquipmentId - kDDLOffset) * fgkNstacks * fgkNlinks +
324           fCurrSlot * fgkNlinks + iLink;
325         if ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) == 0)
326           continue;
327
328         Int_t   size  = 0;
329
330         fErrorFlags = 0;
331         // check for link monitor error flag
332         if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0) {
333           LinkError(kLinkMonitor);
334           if (fgErrorBehav[kLinkMonitor] == kTolerate)
335             size += ReadLinkData();
336         }
337         else
338           // read the data from one HC
339           size += ReadLinkData();
340
341         // read all data endmarkers
342         size += SeekNextLink();
343       }
344
345       // continue with next stack
346       SeekNextStack();
347     }
348   }
349
350   return kTRUE;
351 }
352
353
354 Bool_t AliTRDrawStream::NextDDL()
355 {
356   // continue reading with the next equipment
357
358   if (!fRawReader)
359     return kFALSE;
360
361   fCurrEquipmentId = 0;
362   fCurrSlot = 0;
363   fCurrLink = 0;
364
365   UChar_t *buffer = 0x0;
366
367   while (fRawReader->ReadNextData(buffer)) {
368
369     fCurrEquipmentId = fRawReader->GetEquipmentId();
370     AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
371
372     if (fCurrEquipmentId < kDDLOffset || fCurrEquipmentId > kDDLMax) {
373       EquipmentError(kNonTrdEq, "Skipping");
374       continue;
375     }
376
377     if (fMarkers)
378       new ((*fMarkers)[fMarkers->GetEntriesFast()])
379         AliTRDrawStreamError(-kSecactive, fCurrEquipmentId - kDDLOffset);
380
381     ReadGTUHeaders((UInt_t*) buffer);
382
383     if (fCurrTrailerReadout)
384       ReadGTUTrailer();
385
386     return kTRUE;
387   }
388
389   return kFALSE;
390 }
391
392
393 Int_t AliTRDrawStream::NextChamber(AliTRDdigitsManager *digMgr)
394 {
395   // read the data for the next chamber
396   // in case you only want to read the data of a single chamber
397   // to read all data ReadEvent(...) is recommended
398
399   fDigitsManager = digMgr;
400   fDigitsParam   = 0x0;
401
402   fErrorFlags = 0;
403
404   // tracklet output preparation
405   TTree *trklTree = 0x0;
406   AliRunLoader *rl = AliRunLoader::Instance();
407   AliLoader* trdLoader = rl ? rl->GetLoader("TRDLoader") : NULL;
408   AliDataLoader *trklLoader = trdLoader ? trdLoader->GetDataLoader("tracklets") : NULL;
409   if (trklLoader) {
410     AliTreeLoader *trklTreeLoader = (AliTreeLoader*) trklLoader->GetBaseLoader("tracklets-raw");
411     if (trklTreeLoader)
412       trklTree = trklTreeLoader->Tree();
413     else
414       trklTree = trklLoader->Tree();
415   }
416
417   if (fTrackletTree != trklTree)
418     ConnectTracklets(trklTree);
419
420   if (!fRawReader) {
421     AliError("No raw reader available");
422     return -1;
423   }
424
425   while (fCurrSlot < 0 || fCurrSlot >= fgkNstacks) {
426     if (!NextDDL()) {
427       fCurrSlot = -1;
428       return -1;
429     }
430     while ((fCurrSlot < fgkNstacks) &&
431            (((fCurrStackMask & (1 << fCurrSlot)) == 0) ||
432             ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink))) == 0)) {
433       if ((fCurrStackMask & (1 << fCurrSlot)) == 0) {
434         ++fCurrSlot;
435         fCurrSlot = 0;
436         continue;
437       }
438       fCurrLink++;
439       if (fCurrLink >= fgkNlinks) {
440         SeekNextStack();
441         fCurrLink = 0;
442         fCurrSlot++;
443       }
444     }
445   }
446
447   AliDebug(2, Form("Stack %i, Link %i, mask: 0x%02x", fCurrSlot, fCurrLink, fCurrLinkMask[fCurrSlot]));
448   fCurrHC   = (fCurrEquipmentId - kDDLOffset) * fgkNlinks * fgkNstacks +
449     fCurrSlot * fgkNlinks + fCurrLink;
450
451   if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0) {
452     LinkError(kLinkMonitor);
453     if (fgErrorBehav[kLinkMonitor] == kTolerate)
454       ReadLinkData();
455   }
456   else
457     // read the data from one HC
458     ReadLinkData();
459
460   // read all data endmarkers
461   SeekNextLink();
462
463   if (fCurrLink % 2 == 0) {
464     // if we just read the A-side HC then also check the B-side
465     fCurrLink++;
466     fCurrHC++;
467     if (fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) {
468       if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0) {
469         LinkError(kLinkMonitor);
470         if (fgErrorBehav[kLinkMonitor] == kTolerate)
471           ReadLinkData();
472       }
473       else {
474         ReadLinkData();
475       }
476       SeekNextLink();
477     }
478   }
479
480   do {
481     if ((fCurrStackMask & (1 << fCurrSlot)) == 0) {
482       fCurrLink = 0;
483       fCurrSlot++;
484     }
485     else {
486       fCurrLink++;
487       if (fCurrLink >= fgkNlinks) {
488         SeekNextStack();
489         fCurrLink = 0;
490         fCurrSlot++;
491       }
492     }
493   } while ((fCurrSlot < fgkNstacks) &&
494            (((fCurrStackMask & (1 << fCurrSlot)) == 0) ||
495             ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink))) == 0));
496
497   // return chamber information from HC if it is valid
498   // otherwise return information from link position
499   if (fCurrSm < 0 || fCurrSm >= fgkNsectors || fCurrStack < 0 || fCurrStack >= fgkNstacks || fCurrLayer < 0 || fCurrLayer >= fgkNlinks/2)
500     return ((fCurrEquipmentId-kDDLOffset) + fCurrSlot * fgkNlinks/2 + fCurrLink/2);
501   else
502     return (fCurrSm * fgkNstacks*fgkNlinks/2 + fCurrStack * fgkNlinks/2 + fCurrLayer);
503 }
504
505
506 Int_t AliTRDrawStream::ReadGTUHeaders(UInt_t *buffer)
507 {
508   // check the data source and read the headers
509
510   if (fCurrEquipmentId >= kDDLOffset && fCurrEquipmentId <= kDDLMax) {
511     // this is ROC data
512
513     // setting the pointer to data and current reading position
514     fPayloadCurr = fPayloadStart = buffer;
515     fPayloadSize = fRawReader->GetDataSize() / sizeof(UInt_t);
516     fStats.fStatsSector[fCurrEquipmentId - kDDLOffset].fBytes = fRawReader->GetDataSize();
517     AliDebug(2, Form("Read buffer of size: %i", fRawReader->GetDataSize()));
518
519     AliDebug(1, DumpRaw("raw data", fPayloadCurr, TMath::Min(fPayloadSize, 1000)));
520
521     // read SM header
522     if (ReadSmHeader() < 0) {
523       AliError(Form("Reading SM header failed, skipping this DDL %i", fCurrEquipmentId));
524       return -1;
525     }
526
527     // read tracking headers (if available)
528     if (fCurrTrkHeaderAvail) {
529       for (Int_t iStack = 0; iStack < fgkNstacks; iStack++) {
530         if ((fCurrStackMask & (1 << iStack)) != 0)
531           ReadTrackingHeader(iStack);
532       }
533     }
534
535     // read trigger header(s) (if available)
536     if (fCurrTrgHeaderAvail)
537       ReadTriggerHeaders();
538
539     // read stack header
540     for (Int_t iStack = 0; iStack < fgkNstacks; iStack++) {
541       if ((fCurrStackMask & (1 << iStack)) != 0)
542         ReadStackHeader(iStack);
543     }
544
545     return 0;
546   }
547   else
548     return -1;
549 }
550
551 Int_t AliTRDrawStream::ReadSmHeader()
552 {
553   // read the SMU index header at the current reading position
554   // and store the information in the corresponding variables
555
556   if (fPayloadCurr - fPayloadStart >= fPayloadSize - 1) {
557     EquipmentError(kUnknown, "SM Header incomplete");
558     return -1;
559   }
560
561   fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] = 0;
562
563   fCurrSmHeaderSize           = ((*fPayloadCurr) >> 16) & 0xffff;
564   fCurrSmHeaderVersion        = ((*fPayloadCurr) >> 12) &    0xf;
565   fCurrTrackEnable            = ((*fPayloadCurr) >>  6) &    0x1;
566   fCurrTrackletEnable         = ((*fPayloadCurr) >>  5) &    0x1;
567   fCurrStackMask              = ((*fPayloadCurr)      ) &   0x1f;
568   fCurrHwRev                  = (fPayloadCurr[1] >> 12) & 0xffff;
569   fCurrStackEndmarkerAvail    = 0;
570
571   switch (fCurrSmHeaderVersion) {
572   case 0xb:
573     fCurrTrailerReadout = 0;
574     fCurrTrgHeaderAvail = 0;
575     fCurrEvType = 0;
576     fCurrTrkHeaderAvail = 0;
577
578     DecodeGTUtracks();
579     break;
580
581   case 0xc:
582     fCurrTrailerReadout = ((*fPayloadCurr) >> 10) &    0x1;
583     fCurrTrgHeaderAvail = 1;
584     fCurrTrgHeaderReadout = ((*fPayloadCurr) >>  9) &    0x1;
585     fCurrEvType         = ((*fPayloadCurr) >>  7) &    0x3;
586     fCurrTrkHeaderAvail = fCurrTrackEnable;
587     fCurrTriggerEnable  = (fPayloadCurr[2] >>  8) &  0xfff;
588     fCurrTriggerFired   = (fPayloadCurr[2] >>  20) &  0xfff;
589     fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] = fCurrTriggerFired;
590     break;
591
592   case 0xd:
593     fCurrTrailerReadout = ((*fPayloadCurr) >> 10) &    0x1;
594     fCurrTrgHeaderAvail = 1;
595     fCurrTrgHeaderReadout = ((*fPayloadCurr) >>  9) &    0x1;
596     fCurrEvType         = ((*fPayloadCurr) >>  7) &    0x3;
597     fCurrTrkHeaderAvail = fCurrTrackEnable;
598     fCurrTriggerEnable  = (fPayloadCurr[2] >>  8) &  0xfff;
599     fCurrTriggerFired   = (fPayloadCurr[2] >>  20) &  0xfff;
600     fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] = fCurrTriggerFired;
601     fCurrStackEndmarkerAvail    = 1;
602     break;
603
604   default:
605     AliError(Form("unknown SM header version: 0x%x", fCurrSmHeaderVersion));
606   }
607
608   AliDebug(5, Form("SM header: size: %i, version: %i, track enable: %i, tracklet enable: %i, stack mask: %2x, trailer: %i, trgheader: %i, trkheader: %i",
609                    fCurrSmHeaderSize,
610                    fCurrSmHeaderVersion,
611                    fCurrTrackEnable,
612                    fCurrTrackletEnable,
613                    fCurrStackMask,
614                    fCurrTrailerReadout,
615                    fCurrTrgHeaderAvail,
616                    fCurrTrkHeaderAvail ));
617
618   // jump to the first word after the SM header
619   fPayloadCurr += fCurrSmHeaderSize + 1;
620
621   return fCurrSmHeaderSize + 1;
622 }
623
624 Int_t AliTRDrawStream::DecodeGTUtracks()
625 {
626   // decode GTU track words
627   // this depends on the hardware revision of the SMU
628
629   Int_t sector = fCurrEquipmentId-kDDLOffset;
630
631   if ((sector < 0) || (sector > 17)) {
632     AliError(Form("Invalid sector %i for GTU tracks", sector));
633     return -1;
634   }
635
636   AliDebug(1, DumpRaw(Form("GTU tracks in sector %2i (hw rev %i)", sector, fCurrHwRev),
637                       fPayloadCurr + 4, 10, 0xffe0ffff));
638
639   fCurrTrgFlags[sector] = 0;
640
641   if (fCurrHwRev < 1772) {
642     UInt_t    fastWord;         // fast trigger word
643     ULong64_t trackWord = 0;    // extended track word
644     Int_t stack = 0;
645     Int_t idx = 0;
646     for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
647       if (fPayloadCurr[iWord] == 0x10000000) { // stack boundary marker
648         stack++;
649         idx = 0;
650       }
651       else {
652         if ((idx == 0) &&
653             ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
654           fastWord = fPayloadCurr[iWord];
655           fCurrTrgFlags[sector] |= 1 << (stack+11); // assume tracking done if fast word sent
656           AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fastWord));
657           continue;
658         }
659         else if ((idx & 0x1) == 0x1) {
660           trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
661           AliDebug(1,Form("track debug word: 0x%016llx", trackWord));
662           if (fTracks) {
663             AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
664               AliESDTrdTrack();
665
666             trk->SetSector(sector);
667             trk->SetStack((trackWord >> 60) & 0x7);
668             trk->SetA(0);
669             trk->SetB(0);
670             trk->SetPID(0);
671             trk->SetLayerMask((trackWord >> 16) & 0x3f);
672             trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
673             trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
674             trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
675             trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
676             trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
677             trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
678
679             trk->SetFlags(0);
680             trk->SetReserved(0);
681             trk->SetLabel(-3);
682
683             Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
684             if (TMath::Abs(pt) > 0.1) {
685               trk->SetA((Int_t) (0.15*51625./100./pt / 160e-4 * 2));
686             }
687           }
688         }
689         else {
690           trackWord = fPayloadCurr[iWord];
691         }
692         idx++;
693       }
694     }
695   }
696   else if (fCurrHwRev < 1804) {
697     UInt_t    fastWord;         // fast trigger word
698     ULong64_t trackWord = 0;    // extended track word
699     Int_t stack = 0;
700     Int_t idx = 0;
701     for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
702       if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
703         stack++;
704         idx = 0;
705       }
706       else {
707         if ((idx == 0) &&
708             ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
709           fastWord = fPayloadCurr[iWord];
710           fCurrTrgFlags[sector] |= 1 << (stack+11); // assume tracking done if fast word sent
711           AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fastWord));
712           continue;
713         }
714         else if ((idx & 0x1) == 0x1) {
715           trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
716           AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
717           if (fTracks) {
718             AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
719               AliESDTrdTrack();
720
721             trk->SetSector(fCurrEquipmentId-kDDLOffset);
722             trk->SetStack((trackWord >> 60) & 0x7);
723             trk->SetA(0);
724             trk->SetB(0);
725             trk->SetPID(0);
726             trk->SetLayerMask((trackWord >> 16) & 0x3f);
727             trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
728             trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
729             trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
730             trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
731             trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
732             trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
733
734             trk->SetFlags(0);
735             trk->SetReserved(0);
736             trk->SetLabel(-3);
737
738             Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
739             if (TMath::Abs(pt) > 0.1) {
740               trk->SetA((Int_t) (-0.15*51625./100./pt / 160e-4 * 2));
741             }
742           }
743         }
744         else {
745           trackWord = fPayloadCurr[iWord];
746         }
747         idx++;
748       }
749     }
750   }
751   else if (fCurrHwRev < 1819) {
752     UInt_t    fastWord;         // fast trigger word
753     ULong64_t trackWord = 0;    // extended track word
754     Int_t stack = 0;
755     Int_t idx = 0;
756     for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
757       if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
758         stack++;
759         idx = 0;
760       }
761       else {
762         if ((idx == 0) &&
763             ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
764           fastWord = fPayloadCurr[iWord];
765           if (fastWord & (1 << 13))
766             fCurrTrgFlags[sector] |= 1 << (stack+11);
767           AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fastWord));
768           continue;
769         }
770         else if ((idx & 0x1) == 0x1) {
771           trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
772           AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
773
774           if (fTracks) {
775             AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
776               AliESDTrdTrack();
777
778             trk->SetSector(fCurrEquipmentId-kDDLOffset);
779             trk->SetStack((trackWord >> 60) & 0x7);
780             trk->SetA(0);
781             trk->SetB(0);
782             // trk->SetPt(((trackWord & 0xffff) ^ 0x8000) - 0x8000);
783             trk->SetPID(0);
784             trk->SetLayerMask((trackWord >> 16) & 0x3f);
785             trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
786             trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
787             trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
788             trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
789             trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
790             trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
791
792             trk->SetFlags(0);
793             trk->SetReserved(0);
794             trk->SetLabel(-3);
795
796             Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
797             if (TMath::Abs(pt) > 0.1) {
798               trk->SetA((Int_t) (0.15*51625./100./trk->Pt() / 160e-4 * 2));
799             }
800           }
801         }
802         else {
803           trackWord = fPayloadCurr[iWord];
804         }
805         idx++;
806       }
807     }
808   }
809   else if (fCurrHwRev < 1860) {
810     UInt_t    fastWord;         // fast trigger word
811     ULong64_t trackWord = 0;    // extended track word
812     Int_t stack = 0;
813     Int_t idx = 0;
814     Bool_t upperWord = kFALSE;
815     Int_t word = 0;
816     for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
817       if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
818         stack++;
819         idx = 0;
820         upperWord = kFALSE;
821       }
822       else {
823         // assemble the 32-bit words out of 16-bit blocks
824         if (upperWord) {
825           word |= (fPayloadCurr[iWord] & 0xffff0000);
826           upperWord = kFALSE;
827         }
828         else {
829           // lower word is read first
830           word = (fPayloadCurr[iWord] & 0xffff0000) >> 16;
831           upperWord = kTRUE;
832           continue;
833         }
834
835         if ((word & 0xffff0008) == 0x13370008) {
836           fastWord = word;
837           AliDebug(1, Form("stack %i: fast track word: 0x%08x", stack, fastWord));
838           if (fastWord & (1 << 13))
839             fCurrTrgFlags[sector] |= 1 << (stack+11);
840           continue;
841         }
842         else if ((idx & 0x1) == 0x1) {
843           trackWord |= ((ULong64_t) word) << 32;
844           AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
845           if (fTracks) {
846             AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
847               AliESDTrdTrack();
848
849             trk->SetSector(fCurrEquipmentId-kDDLOffset);
850             trk->SetStack((trackWord >> 60) & 0x7);
851             trk->SetA(0);
852             trk->SetB(0);
853             trk->SetPID(0);
854             trk->SetLayerMask((trackWord >> 16) & 0x3f);
855             trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
856             trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
857             trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
858             trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
859             trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
860             trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
861
862             trk->SetFlags(0);
863             trk->SetReserved(0);
864             trk->SetLabel(-3);
865
866             Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
867             if (TMath::Abs(pt) > 0.1) {
868               trk->SetA((Int_t) (0.15*51625./100./pt / 160e-4 * 2));
869             }
870           }
871         }
872         else {
873           trackWord = word;
874         }
875         idx++;
876       }
877     }
878
879   }
880   else {
881     ULong64_t trackWord = 0; // this is the debug word
882     Int_t stack = 0;
883     Int_t idx = 0;
884     Bool_t upperWord = kFALSE;
885     Int_t word = 0;
886     for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
887       if (fPayloadCurr[iWord] == 0xffe0ffff) {
888         stack++;
889         idx = 0;
890         upperWord = kFALSE;
891       }
892       else {
893         // assemble the 32-bit words out of 16-bit blocks
894         if (upperWord) {
895           word |= (fPayloadCurr[iWord] & 0xffff0000);
896           upperWord = kFALSE;
897         }
898         else {
899           // lower word is read first
900           word = (fPayloadCurr[iWord] & 0xffff0000) >> 16;
901           upperWord = kTRUE;
902           continue;
903         }
904
905         if ((word & 0xffff0008) == 0x13370008) {
906           AliDebug(1, Form("stack %i: fast track word: 0x%08x", stack, word));
907           continue;
908         }
909         else if ((word & 0xffff0010) == 0x13370010) {
910           AliDebug(1, Form("stack %i: tracking done word: 0x%08x", stack, word));
911           fCurrTrgFlags[sector] |= 1 << (stack+11);
912           continue;
913         }
914         else if ((idx & 0x1) == 0x1) {
915           trackWord |= ((ULong64_t) word) << 32;
916           AliDebug(1, Form("track debug word: 0x%16llx", trackWord));
917           if (fTracks) {
918             AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
919               AliESDTrdTrack();
920             trk->SetSector(fCurrEquipmentId-kDDLOffset);
921             trk->SetStack((trackWord >> 60) & 0x7);
922             trk->SetA(0);
923             trk->SetB(0);
924             trk->SetPID(0);
925             trk->SetLayerMask((trackWord >> 16) & 0x3f);
926             trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
927             trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
928             trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
929             trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
930             trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
931             trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
932
933             trk->SetFlags(0);
934             trk->SetReserved(0);
935             trk->SetLabel(-3);
936
937             Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
938             if (TMath::Abs(pt) > 0.1) {
939               trk->SetA(-(Int_t) (0.15*51625./100./pt / 160e-4 * 2));
940             }
941           }
942         }
943         else {
944           trackWord = word;
945         }
946         idx++;
947       }
948     }
949   }
950   return 0;
951 }
952
953 Int_t AliTRDrawStream::ReadTrackingHeader(Int_t stack)
954 {
955   // read the tracking information and store it for the given stack
956
957   // index word
958
959   fCurrTrkHeaderIndexWord[stack] = *fPayloadCurr;
960   fCurrTrkHeaderSize[stack]      = ((*fPayloadCurr) >> 16) & 0x3ff;
961
962   AliDebug(1, Form("tracking header index word: 0x%08x, size: %i (hw rev: %i)",
963                    fCurrTrkHeaderIndexWord[stack], fCurrTrkHeaderSize[stack], fCurrHwRev));
964   Int_t trackingTime = *fPayloadCurr & 0x3ff;
965
966   fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] |= ((fCurrTrkHeaderIndexWord[stack] >> 10) & 0x1) << (22 + stack);
967   fPayloadCurr++;
968
969   // data words
970   ULong64_t trackWord = 0;
971   Int_t idx = 0;
972   Int_t trackIndex = fTracks ? fTracks->GetEntriesFast() : -1;
973
974   for (UInt_t iWord = 0; iWord < fCurrTrkHeaderSize[stack]; iWord++) {
975
976     if (!(idx & 0x1)) {
977       // first part of 64-bit word
978       trackWord = fPayloadCurr[iWord];
979     }
980     else {
981       trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
982
983       if (trackWord & (1ul << 63)) {
984         if ((trackWord & (0x3ful << 56)) != 0) {
985           // track word
986           AliDebug(2, Form("track word: 0x%016llx", trackWord));
987
988           if (fTracks) {
989             AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
990               AliESDTrdTrack();
991
992             trk->SetSector(fCurrEquipmentId-kDDLOffset);
993             trk->SetLayerMask((trackWord >> 56) & 0x3f);
994             trk->SetA( (((trackWord >> 38) & 0x3ffff) ^ 0x20000) - 0x20000);
995             trk->SetB( (((trackWord >> 20) & 0x3ffff) ^ 0x20000) - 0x20000);
996             trk->SetC( (((trackWord >> 8)  &  0xffff) ^  0x8000) -  0x8000);
997             trk->SetPID((trackWord >>  0) & 0xff);
998             trk->SetStack(stack);
999             trk->SetLabel(-3);
1000
1001             // now compare the track word with the one generated from the ESD information
1002             if (trackWord != trk->GetTrackWord(0)) {
1003               AliError(Form("track word 0x%016llx does not match the read one 0x%016llx",
1004                             trk->GetTrackWord(0), trackWord));
1005             }
1006           }
1007         }
1008         else {
1009           // done marker (so far only used to set trigger flag)
1010           fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] |= 1 << (27 + stack);
1011           fCurrTrkFlags[(fCurrEquipmentId-kDDLOffset)*fgkNstacks + stack] = trackWord;
1012
1013           AliDebug(2, Form("tracking done marker: 0x%016llx, trigger flags: 0x%08x",
1014                            trackWord, fCurrTrgFlags[fCurrEquipmentId-kDDLOffset]));
1015           AliDebug(2, Form("seg / stack / first / last / done / index : %i %i %lli %lli %lli %i",
1016                            fCurrEquipmentId - kDDLOffset, stack,
1017                            (trackWord >> 20) & 0x3ff,
1018                            (trackWord >> 10) & 0x3ff,
1019                            (trackWord >>  0) & 0x3ff,
1020                            trackingTime));
1021         }
1022       }
1023       else {
1024         // extended track word
1025         AliDebug(2, Form("extended track word: 0x%016llx", trackWord));
1026
1027         if (fTracks) {
1028           AliESDTrdTrack *trk = (AliESDTrdTrack*) (*fTracks)[trackIndex];
1029
1030           trk->SetFlags((trackWord >> 52) & 0x7ff);
1031           trk->SetFlagsTiming((trackWord >> 51) & 0x1);
1032           trk->SetReserved((trackWord >> 49) & 0x3);
1033           trk->SetY((trackWord >> 36) & 0x1fff);
1034           trk->SetTrackletIndex((trackWord >>  0) & 0x3f, 0);
1035           trk->SetTrackletIndex((trackWord >>  6) & 0x3f, 1);
1036           trk->SetTrackletIndex((trackWord >> 12) & 0x3f, 2);
1037           trk->SetTrackletIndex((trackWord >> 18) & 0x3f, 3);
1038           trk->SetTrackletIndex((trackWord >> 24) & 0x3f, 4);
1039           trk->SetTrackletIndex((trackWord >> 30) & 0x3f, 5);
1040
1041           if (trackWord != trk->GetExtendedTrackWord(0)) {
1042             AliError(Form("extended track word 0x%016llx does not match the read one 0x%016llx",
1043                           trk->GetExtendedTrackWord(0), trackWord));
1044             }
1045
1046           trackIndex++;
1047         }
1048       }
1049     }
1050     idx++;
1051   }
1052
1053   fPayloadCurr += fCurrTrkHeaderSize[stack];
1054
1055   return fCurrTrkHeaderSize[stack];
1056 }
1057
1058 Int_t AliTRDrawStream::ReadTriggerHeaders()
1059 {
1060   // read all trigger headers present
1061
1062   AliDebug(1, Form("trigger mask: 0x%03x, fired: 0x%03x\n",
1063                    fCurrTriggerEnable, fCurrTriggerFired));
1064   // loop over potential trigger blocks
1065   for (Int_t iTrigger = 0; iTrigger < fgkNtriggers; iTrigger++) {
1066     // check for trigger enable
1067     if (fCurrTriggerEnable & (1 << iTrigger)) {
1068       // check for readout mode and trigger fired
1069       if ((fCurrTrgHeaderReadout == 0) || (fCurrTriggerFired & (1 << iTrigger))) {
1070         // index word
1071         AliDebug(1, Form("trigger index word %i: 0x%08x\n", iTrigger, *fPayloadCurr));
1072         fCurrTrgHeaderIndexWord[iTrigger] = *fPayloadCurr;
1073         fCurrTrgHeaderSize[iTrigger]      = ((*fPayloadCurr) >> 16) & 0x3ff;
1074         if (iTrigger == 7) {
1075           // timeout trigger, use to extract tracking time
1076           fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] |= (*fPayloadCurr & 0x3ff) << 12;
1077         }
1078
1079         fPayloadCurr++;
1080         // data words
1081         fPayloadCurr += fCurrTrgHeaderSize[iTrigger];
1082       }
1083     }
1084   }
1085
1086   return 0;
1087 }
1088
1089 Int_t AliTRDrawStream::ReadStackHeader(Int_t stack)
1090 {
1091   // read the stack header
1092   // and store the information in the corresponding variables
1093
1094   fCurrStackIndexWord[stack]     = *fPayloadCurr;
1095   fCurrStackHeaderSize[stack]    = (((*fPayloadCurr) >> 16) & 0xffff) + 1;
1096   fCurrStackHeaderVersion[stack] = ((*fPayloadCurr) >> 12) & 0xf;
1097   fCurrLinkMask[stack]           = (*fPayloadCurr) & 0xfff;
1098
1099   // dumping stack header
1100   AliDebug(1, DumpRaw(Form("stack %i header", stack), fPayloadCurr, fCurrStackHeaderSize[stack]));
1101
1102   if (fPayloadCurr - fPayloadStart >= fPayloadSize - (Int_t) fCurrStackHeaderSize[stack]) {
1103     EquipmentError(kStackHeaderInvalid, "Stack index header %i incomplete", stack);
1104     // dumping stack header
1105     AliError(DumpRaw(Form("stack %i header", stack), fPayloadCurr, fCurrStackHeaderSize[stack]));
1106
1107     return -1;
1108   }
1109
1110   switch (fCurrStackHeaderVersion[stack]) {
1111   case 0xa:
1112     if (fCurrStackHeaderSize[stack] < 8) {
1113       LinkError(kStackHeaderInvalid, "Stack header smaller than expected!");
1114       return -1;
1115     }
1116
1117     fCurrCleanCheckout[stack] = fPayloadCurr[1] & 0x1;
1118     fCurrBoardId[stack]       = (fPayloadCurr[1] >> 8) & 0xff;
1119     fCurrHwRevTMU[stack]      = (fPayloadCurr[1] >> 16) & 0xffff;
1120
1121     for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1122       // A side
1123       fCurrLinkMonitorFlags  [stack*fgkNlinks + iLayer*2]      = fPayloadCurr[iLayer+2] & 0xf;
1124       fCurrLinkDataTypeFlags [stack*fgkNlinks + iLayer*2]      = (fPayloadCurr[iLayer+2] >> 4) & 0x3;
1125       fCurrLinkDebugFlags    [stack*fgkNlinks + iLayer*2]      = (fPayloadCurr[iLayer+2] >> 12) & 0xf;
1126       // B side
1127       fCurrLinkMonitorFlags  [stack*fgkNlinks + iLayer*2 + 1]  = (fPayloadCurr[iLayer+2] >> 16) & 0xf;
1128       fCurrLinkDataTypeFlags [stack*fgkNlinks + iLayer*2 + 1]  = (fPayloadCurr[iLayer+2] >> 20) & 0x3;
1129       fCurrLinkDebugFlags    [stack*fgkNlinks + iLayer*2 + 1]  = (fPayloadCurr[iLayer+2] >> 28) & 0xf;
1130     }
1131     break;
1132
1133   default:
1134     LinkError(kStackHeaderInvalid, "Invalid Stack Header version %x", fCurrStackHeaderVersion[stack]);
1135   }
1136
1137   fPayloadCurr += fCurrStackHeaderSize[stack];
1138
1139   return fCurrStackHeaderSize[stack];
1140 }
1141
1142 Int_t AliTRDrawStream::ReadGTUTrailer()
1143 {
1144   // read the SM trailer containing CRCs from various stages
1145
1146   UInt_t* trailer = fPayloadStart + fPayloadSize -1;
1147
1148   // look for the trailer index word from the end
1149   for (Int_t iWord = 0; iWord < fPayloadSize; iWord++) {
1150     if ((fPayloadStart[fPayloadSize-1-iWord] & 0xffff) == 0x1f51) {
1151       trailer = fPayloadStart + fPayloadSize - 1 - iWord;
1152       break;
1153     }
1154   }
1155
1156   if (((*trailer) & 0xffff) == 0x1f51) {
1157     UInt_t trailerIndexWord = (*trailer);
1158     Int_t trailerSize = (trailerIndexWord >> 16) & 0xffff;
1159     AliDebug(2, DumpRaw("GTU trailer", trailer, trailerSize+1));
1160     // parse the trailer
1161     if (trailerSize >= 4) {
1162       // match flags from GTU
1163       fCurrMatchFlagsSRAM     = (trailer[1] >>  0) &    0x1f;
1164       fCurrMatchFlagsPostBP   = (trailer[1] >>  5) &    0x1f;
1165       // individual checksums
1166       fCurrChecksumStack[0]   = (trailer[1] >> 16) &  0xffff;
1167       fCurrChecksumStack[1]   = (trailer[2] >>  0) &  0xffff;
1168       fCurrChecksumStack[2]   = (trailer[2] >> 16) &  0xffff;
1169       fCurrChecksumStack[3]   = (trailer[3] >>  0) &  0xffff;
1170       fCurrChecksumStack[4]   = (trailer[3] >> 16) &  0xffff;
1171       fCurrChecksumSIU        = trailer[4];
1172
1173       if ((fCurrMatchFlagsSRAM & fCurrStackMask) != fCurrStackMask)
1174         EquipmentError(kCRCmismatch, "CRC mismatch SRAM: 0x%02x", fCurrMatchFlagsSRAM);
1175       if ((fCurrMatchFlagsPostBP & fCurrStackMask) != fCurrStackMask)
1176         EquipmentError(kCRCmismatch, "CRC mismatch BP: 0x%02x", fCurrMatchFlagsPostBP);
1177
1178     }
1179     else {
1180       LinkError(kUnknown, "Invalid GTU trailer");
1181     }
1182   }
1183   else
1184     EquipmentError(kUnknown, "trailer index marker mismatch");
1185
1186   return 0;
1187 }
1188
1189 Int_t AliTRDrawStream::ReadLinkData()
1190 {
1191   // read the data in one link (one HC) until the data endmarker is reached
1192   // returns the number of words read!
1193
1194   Int_t count = 0;
1195   UInt_t* startPosLink = fPayloadCurr;
1196
1197   AliDebug(1, DumpRaw(Form("link data from seg %2i slot %i link %2i", fCurrEquipmentId-kDDLOffset, fCurrSlot, fCurrLink),
1198                       fPayloadCurr, TMath::Min((Int_t) (fPayloadSize - (fPayloadCurr-fPayloadStart)), 100), 0x00000000));
1199
1200   if (fMarkers)
1201     new ((*fMarkers)[fMarkers->GetEntriesFast()])
1202       AliTRDrawStreamError(-kHCactive, fCurrEquipmentId-kDDLOffset, fCurrSlot, fCurrLink);
1203
1204   if (fErrorFlags & kDiscardHC)
1205     return count;
1206
1207   if (fCurrTrackletEnable) {
1208     count += ReadTracklets();
1209     if (fErrorFlags & kDiscardHC)
1210       return count;
1211   }
1212
1213   AliDebug(1, DumpRaw("HC header", fPayloadCurr, 4, 0x00000000));
1214   count += ReadHcHeader();
1215   if (fErrorFlags & kDiscardHC)
1216     return count;
1217
1218   Int_t det = fCurrSm * 30 + fCurrStack * 6 + fCurrLayer;
1219
1220   if (det > -1 && det < 540) {
1221
1222     // ----- check which kind of data -----
1223     if (fCurrMajor & 0x40) {
1224       if ((fCurrMajor & 0x7) == 0x7) {
1225         AliDebug(1, "This is a config event");
1226         UInt_t *startPos = fPayloadCurr;
1227         while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1228                *fPayloadCurr != fgkDataEndmarker)
1229           fPayloadCurr++;
1230         count += fPayloadCurr - startPos;
1231
1232         // feeding TRAP config
1233         AliTRDtrapConfig *trapcfg = AliTRDcalibDB::Instance()->GetTrapConfig();
1234         AliTRDmcmSim::ReadPackedConfig(trapcfg, fCurrHC, startPos, fPayloadCurr - startPos);
1235       }
1236       else {
1237         Int_t tpmode = fCurrMajor & 0x7;
1238         AliDebug(1, Form("Checking testpattern (mode %i) data", tpmode));
1239         count += ReadTPData(tpmode);
1240       }
1241     }
1242     else {
1243       // reading real data
1244       if (fDigitsManager) {
1245         if ((fAdcArray = fDigitsManager->GetDigits(det))) {
1246           //fAdcArray->Expand();
1247           if (fAdcArray->GetNtime() != fCurrNtimebins)
1248             fAdcArray->Allocate(16, 144, fCurrNtimebins);
1249         }
1250         else {
1251           LinkError(kNoDigits);
1252         }
1253
1254         if (!fDigitsParam) {
1255           fDigitsParam = fDigitsManager->GetDigitsParam();
1256         }
1257         if (fDigitsParam) {
1258           fDigitsParam->SetPretriggerPhase(det, fCurrPtrgPhase);
1259           fDigitsParam->SetNTimeBins(det, fCurrNtimebins);
1260           fDigitsParam->SetADCbaseline(det, 10);
1261         }
1262
1263         if (fDigitsManager->UsesDictionaries()) {
1264           fDigitsManager->GetDictionary(det, 0)->Reset();
1265           fDigitsManager->GetDictionary(det, 1)->Reset();
1266           fDigitsManager->GetDictionary(det, 2)->Reset();
1267         }
1268
1269         if ((fSignalIndex = fDigitsManager->GetIndexes(det))) {
1270           fSignalIndex->SetSM(fCurrSm);
1271           fSignalIndex->SetStack(fCurrStack);
1272           fSignalIndex->SetLayer(fCurrLayer);
1273           fSignalIndex->SetDetNumber(det);
1274           if (!fSignalIndex->IsAllocated())
1275             fSignalIndex->Allocate(16, 144, fCurrNtimebins);
1276         }
1277
1278         if (fCurrMajor & 0x20) {
1279           AliDebug(1, "This is a zs event");
1280           count += ReadZSData();
1281         }
1282         else {
1283           AliDebug(1, "This is a nozs event");
1284           count += ReadNonZSData();
1285         }
1286       }
1287       else {
1288         // just read until data endmarkers
1289         while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1290                *fPayloadCurr != fgkDataEndmarker)
1291           fPayloadCurr++;
1292       }
1293     }
1294   }
1295   else {
1296     LinkError(kInvalidDetector, "%i", det);
1297     while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1298            *fPayloadCurr != fgkDataEndmarker)
1299       fPayloadCurr++;
1300   }
1301
1302   if (fCurrSm > -1 && fCurrSm < 18) {
1303     fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fBytes     += (fPayloadCurr - startPosLink) * sizeof(UInt_t);
1304     fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fBytesRead += count * sizeof(UInt_t);
1305     fStats.fStatsSector[fCurrSm].fBytesRead                      += count * sizeof(UInt_t);
1306     fStats.fBytesRead                                            += count * sizeof(UInt_t);
1307   }
1308
1309   return count;
1310 }
1311
1312 Int_t AliTRDrawStream::ReadTracklets()
1313 {
1314   // read the tracklets from one HC
1315
1316   fTrackletArray->Clear();
1317
1318   UInt_t *start = fPayloadCurr;
1319   while (*(fPayloadCurr) != fgkTrackletEndmarker &&
1320          *(fPayloadCurr) != fgkStackEndmarker[0] &&
1321          *(fPayloadCurr) != fgkStackEndmarker[1] &&
1322          fPayloadCurr - fPayloadStart < (fPayloadSize - 1)) {
1323     new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(*(fPayloadCurr), fCurrHC);
1324
1325     fPayloadCurr++;
1326   }
1327
1328   if (fTrackletArray->GetEntriesFast() > 0) {
1329     AliDebug(1, Form("Found %i tracklets in %i %i %i (ev. %i)", fTrackletArray->GetEntriesFast(),
1330                      (fCurrEquipmentId-kDDLOffset), fCurrSlot, fCurrLink, fRawReader->GetEventIndex()));
1331     if (fCurrSm > -1 && fCurrSm < 18) {
1332       fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNTracklets += fTrackletArray->GetEntriesFast();
1333       fStats.fStatsSector[fCurrSm].fNTracklets                      += fTrackletArray->GetEntriesFast();
1334     }
1335     if (fTrackletTree)
1336       fTrackletTree->Fill();
1337     if (fTracklets)
1338       for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
1339         new ((*fTracklets)[fTracklets->GetEntriesFast()]) AliTRDtrackletWord(*((AliTRDtrackletWord*)(*fTrackletArray)[iTracklet]));
1340       }
1341   }
1342
1343   // loop over remaining tracklet endmarkers
1344   while ((*(fPayloadCurr) == fgkTrackletEndmarker &&
1345           fPayloadCurr - fPayloadStart < fPayloadSize))
1346     fPayloadCurr++;
1347
1348   return fPayloadCurr - start;
1349 }
1350
1351 Int_t AliTRDrawStream::ReadHcHeader()
1352 {
1353   // read and parse the HC header of one HC
1354   // and store the information in the corresponding variables
1355
1356   AliDebug(1, Form("HC header: 0x%08x", *fPayloadCurr));
1357   UInt_t *start = fPayloadCurr;
1358   // check not to be at the data endmarker
1359   if (*fPayloadCurr == fgkDataEndmarker ||
1360       *(fPayloadCurr) == fgkStackEndmarker[0] ||
1361       *(fPayloadCurr) == fgkStackEndmarker[1]) {
1362     LinkError(kHCmismatch, "found endmarker where HC header should be");
1363     return 0;
1364   }
1365
1366   fCurrSpecial    = (*fPayloadCurr >> 31) & 0x1;
1367   fCurrMajor      = (*fPayloadCurr >> 24) & 0x7f;
1368   fCurrMinor      = (*fPayloadCurr >> 17) & 0x7f;
1369   fCurrAddHcWords = (*fPayloadCurr >> 14) & 0x7;
1370   fCurrSm         = (*fPayloadCurr >> 9) & 0x1f;
1371   fCurrLayer      = (*fPayloadCurr >> 6) & 0x7;
1372   fCurrStack      = (*fPayloadCurr >> 3) & 0x7;
1373   fCurrSide       = (*fPayloadCurr >> 2) & 0x1;
1374   fCurrCheck      = (*fPayloadCurr) & 0x3;
1375
1376   if ((fCurrSm != (((Int_t) fCurrEquipmentId) - kDDLOffset)) ||
1377       (fCurrStack != fCurrSlot) ||
1378       (fCurrLayer != fCurrLink / 2) ||
1379       (fCurrSide != fCurrLink % 2)) {
1380     LinkError(kHCmismatch,
1381               "HC: %i, %i, %i, %i\n 0x%08x 0x%08x 0x%08x 0x%08x",
1382               fCurrSm, fCurrStack, fCurrLayer, fCurrSide,
1383               fPayloadCurr[0], fPayloadCurr[1], fPayloadCurr[2], fPayloadCurr[3]);
1384   }
1385   if (fCurrCheck != 0x1) {
1386     LinkError(kHCcheckFailed);
1387   }
1388
1389   if (fCurrAddHcWords > 0) {
1390     fCurrNtimebins = (fPayloadCurr[1] >> 26) & 0x3f;
1391     fCurrBC = (fPayloadCurr[1] >> 10) & 0xffff;
1392     fCurrPtrgCnt = (fPayloadCurr[1] >> 6) & 0xf;
1393     fCurrPtrgPhase = (fPayloadCurr[1] >> 2) & 0xf;
1394   }
1395
1396   fPayloadCurr += 1 + fCurrAddHcWords;
1397
1398   return (fPayloadCurr - start);
1399 }
1400
1401 Int_t AliTRDrawStream::ReadTPData(Int_t mode)
1402 {
1403   // testing of testpattern 1 to 3 (hardcoded), 0 missing
1404   // evcnt checking missing
1405   Int_t cpu = 0;
1406   Int_t cpufromchannel[] = {0, 0, 0, 0, 0,  1, 1, 1, 1, 1,  2, 2, 2, 2, 2,  3, 3, 3, 3, 3, 3};
1407   Int_t evno  = -1;
1408   Int_t evcnt = 0;
1409   Int_t count = 0;
1410   Int_t mcmcount = -1;
1411   Int_t wordcount = 0;
1412   Int_t channelcount = 0;
1413   UInt_t expword = 0;
1414   UInt_t expadcval = 0;
1415   UInt_t diff = 0;
1416   Int_t lastmcmpos = -1;
1417   Int_t lastrobpos = -1;
1418
1419   UInt_t* start = fPayloadCurr;
1420
1421   while (*(fPayloadCurr) != fgkDataEndmarker &&
1422          fPayloadCurr - fPayloadStart < fPayloadSize - 1) {
1423
1424     // ----- Checking MCM Header -----
1425     AliDebug(2, DumpMcmHeader("MCM header: ", *fPayloadCurr));
1426     UInt_t *startPosMCM = fPayloadCurr;
1427     mcmcount++;
1428
1429     // ----- checking for proper readout order - ROB -----
1430     fCurrRobPos = ROB(*fPayloadCurr);
1431     if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1432       if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) > lastrobpos)
1433         lastmcmpos = -1;
1434       lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1435     }
1436     else {
1437       ROBError(kPosUnexp, Form("#%i after #%i in readout order", GetROBReadoutPos(ROB(*fPayloadCurr) / 2), lastrobpos));
1438     }
1439
1440     // ----- checking for proper readout order - MCM -----
1441     fCurrMcmPos = MCM(*fPayloadCurr);
1442     if (GetMCMReadoutPos(MCM(*fPayloadCurr)) > lastmcmpos) {
1443       lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
1444     }
1445     else {
1446       MCMError(kPosUnexp, Form("#%i after #%i in readout order", GetMCMReadoutPos(MCM(*fPayloadCurr)), lastmcmpos));
1447     }
1448
1449     if (EvNo(*fPayloadCurr) != (evno & 0xfffff)) {
1450       if (evno == -1) {
1451         evno = EvNo(*fPayloadCurr);
1452       }
1453       else {
1454         MCMError(kEvCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
1455       }
1456     }
1457
1458     fPayloadCurr++;
1459
1460     evcnt = 0x3f & *fPayloadCurr >> 26;
1461     cpu = -1;
1462     channelcount = 0;
1463     while (channelcount < 21) {
1464       count = 0;
1465       if (cpu != cpufromchannel[channelcount]) {
1466         cpu = cpufromchannel[channelcount];
1467         expadcval = (1 << 9) | (fCurrRobPos << 6) | (fCurrMcmPos << 2) | cpu;
1468         wordcount = 0;
1469       }
1470
1471       while (count < 10) {
1472         if (*fPayloadCurr == fgkDataEndmarker) {
1473           MCMError(kMissTpData);
1474           return (fPayloadCurr - start);
1475         }
1476
1477         if (channelcount % 2 == 0)
1478           expword = 0x3;
1479         else
1480           expword = 0x2;
1481
1482         if (mode == 1) {
1483           // ----- TP 1 -----
1484           expword |= expadcval << 2;
1485           expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
1486           expword |= expadcval << 12;
1487           expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
1488           expword |= expadcval << 22;
1489           expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
1490         }
1491         else if (mode == 2) {
1492           // ----- TP 2 ------
1493           expword = ((0x3f & evcnt) << 26) | ((fCurrSm + 1) << 21) | ((fCurrLayer + 1) << 18) |
1494             ((fCurrStack + 1) << 15) |
1495             (fCurrRobPos << 12) | (fCurrMcmPos << 8) | (cpu << 6) | (wordcount + 1);
1496         }
1497         else if (mode == 3) {
1498           // ----- TP 3 -----
1499           expword = ((0xfff & evcnt) << 20) | (fCurrSm << 15) | (fCurrLink/2 << 12) | (fCurrStack << 9) |
1500             (fCurrRobPos << 6) | (fCurrMcmPos << 2) | (cpu << 0);
1501         }
1502         else {
1503           expword = 0;
1504           LinkError(kTPmodeInvalid, "Just reading");
1505         }
1506
1507         diff = *fPayloadCurr ^ expword;
1508         AliDebug(11, Form("Comparing ch %2i, word %2i (cpu %i): 0x%08x <-> 0x%08x",
1509                           channelcount, wordcount, cpu, *fPayloadCurr, expword));
1510
1511         if (diff != 0) {
1512           MCMError(kTPmismatch,
1513                    "Seen 0x%08x, expected 0x%08x, diff: 0x%08x, 0x%04x, 0x%02x - word %2i (cpu %i, ch %i)",
1514                    *fPayloadCurr, expword, diff,
1515                    0xffff & (diff | diff >> 16),
1516                    0xff & (diff | diff >> 8 | diff >> 16 | diff >> 24),
1517                    wordcount, cpu, channelcount);;
1518         }
1519         fPayloadCurr++;
1520         count++;
1521         wordcount++;
1522         if (*fPayloadCurr == fgkDataEndmarker)
1523           return (fPayloadCurr - start);
1524       }
1525       channelcount++;
1526     }
1527     // continue with next MCM
1528
1529     if (IsDumping() && DumpingMCM(fCurrHC/2, fCurrRobPos, fCurrMcmPos)) {
1530       AliInfo(DumpRaw(Form("Event %i: Det %3i ROB %i MCM %2i", fRawReader->GetEventIndex(), fCurrHC/2, fCurrRobPos, fCurrMcmPos),
1531                       startPosMCM, fPayloadCurr - startPosMCM));
1532     }
1533
1534   }
1535   return fPayloadCurr - start;
1536 }
1537
1538
1539 Int_t AliTRDrawStream::ReadZSData()
1540 {
1541   // read the zs data from one link from the current reading position
1542
1543   UInt_t *start = fPayloadCurr;
1544
1545   Int_t mcmcount = 0;
1546   Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
1547   Int_t channelcount = 0;
1548   Int_t channelcountExp = 0;
1549   Int_t channelcountMax = 0;
1550   Int_t timebins;
1551   Int_t currentTimebin = 0;
1552   Int_t adcwc = 0;
1553   Int_t evno = -1;
1554   Int_t lastmcmpos = -1;
1555   Int_t lastrobpos = -1;
1556
1557   if (fCurrNtimebins != fNtimebins) {
1558     if (fNtimebins > 0)
1559       LinkError(kNtimebinsChanged,
1560                 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
1561     fNtimebins = fCurrNtimebins;
1562   }
1563
1564   timebins = fNtimebins;
1565
1566   while (*(fPayloadCurr) != fgkDataEndmarker &&
1567          *(fPayloadCurr) != fgkStackEndmarker[0] &&
1568          *(fPayloadCurr) != fgkStackEndmarker[1] &&
1569          fPayloadCurr - fPayloadStart < fPayloadSize) {
1570
1571     // ----- Checking MCM Header -----
1572     AliDebug(2, DumpMcmHeader("MCM header: ", *fPayloadCurr));
1573     UInt_t *startPosMCM = fPayloadCurr;
1574
1575     // ----- checking for proper readout order - ROB -----
1576     fCurrRobPos = ROB(*fPayloadCurr);
1577     if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1578       if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) > lastrobpos)
1579         lastmcmpos = -1;
1580       lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1581     }
1582     else {
1583       ROBError(kPosUnexp, Form("#%i after #%i and #%i in readout order",
1584                                GetROBReadoutPos(ROB(*fPayloadCurr) / 2), lastrobpos, GetROBReadoutPos(fCurrRobPos)));
1585     }
1586
1587     // ----- checking for proper readout order - MCM -----
1588     fCurrMcmPos = MCM(*fPayloadCurr);
1589     if (GetMCMReadoutPos(MCM(*fPayloadCurr)) > lastmcmpos) {
1590       lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
1591     }
1592     else {
1593       MCMError(kPosUnexp, Form("#%i after #%i and #%i in readout order",
1594                                GetMCMReadoutPos(MCM(*fPayloadCurr)), lastmcmpos, GetMCMReadoutPos(fCurrMcmPos)));
1595     }
1596
1597     if (EvNo(*fPayloadCurr) != (evno & 0xfffff)) {
1598       if (evno == -1) {
1599         evno = EvNo(*fPayloadCurr);
1600       }
1601       else {
1602         MCMError(kEvCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
1603       }
1604     }
1605     Int_t adccoloff = AdcColOffset(*fPayloadCurr);
1606     Int_t padcoloff = PadColOffset(*fPayloadCurr);
1607     Int_t row = Row(*fPayloadCurr);
1608     fPayloadCurr++;
1609
1610     if ((row > 11) && (fCurrStack == 2)) {
1611       MCMError(kUnknown, "Data in padrow > 11 for stack 2");
1612     }
1613
1614     if (fErrorFlags & (kDiscardMCM | kDiscardHC | kDiscardDDL))
1615       break; //???
1616
1617     // ----- Reading ADC channels -----
1618     AliDebug(2, DumpAdcMask("ADC mask: ", *fPayloadCurr));
1619
1620     // ----- analysing the ADC mask -----
1621     channelcount = 0;
1622     channelcountExp = GetNActiveChannelsFromMask(*fPayloadCurr);
1623     channelcountMax = GetNActiveChannels(*fPayloadCurr);
1624     Int_t channelmask = GetActiveChannels(*fPayloadCurr);
1625     Int_t channelno = -1;
1626     fPayloadCurr++;
1627
1628     if (channelcountExp != channelcountMax) {
1629       if (channelcountExp > channelcountMax) {
1630         Int_t temp = channelcountExp;
1631         channelcountExp = channelcountMax;
1632         channelcountMax = temp;
1633       }
1634       while (channelcountExp < channelcountMax && channelcountExp < 21 &&
1635              fPayloadCurr - fPayloadStart < fPayloadSize - 10 * channelcountExp - 1) {
1636         MCMError(kAdcMaskInconsistent,
1637                  "Possible MCM-H: 0x%08x, possible ADC-mask: 0x%08x",
1638                  *(fPayloadCurr + 10 * channelcountExp),
1639                  *(fPayloadCurr + 10 * channelcountExp + 1) );
1640         if (!CouldBeMCMhdr( *(fPayloadCurr + 10 * channelcountExp)) && !CouldBeADCmask( *(fPayloadCurr + 10 * channelcountExp + 1)))
1641           channelcountExp++;
1642         else {
1643           break;
1644         }
1645       }
1646       MCMError(kAdcMaskInconsistent,
1647                "Inconsistency in no. of active channels: Counter: %i, Mask: %i, chosen: %i!",
1648                GetNActiveChannels(fPayloadCurr[-1]), GetNActiveChannelsFromMask(fPayloadCurr[-1]), channelcountExp);
1649     }
1650     AliDebug(2, Form("expecting %i active channels, %i timebins", channelcountExp, fCurrNtimebins));
1651
1652     // ----- reading marked ADC channels -----
1653     while (channelcount < channelcountExp && *(fPayloadCurr) != fgkDataEndmarker) {
1654       if (channelno < 20)
1655         channelno++;
1656       while (channelno < 20 && (channelmask & 1 << channelno) == 0)
1657         channelno++;
1658
1659       if (fCurrNtimebins > 30) {
1660         currentTimebin = ((*fPayloadCurr >> 2) & 0x3f);
1661         timebins = ((*fPayloadCurr >> 8) & 0xf) * 3;
1662       }
1663       else {
1664         currentTimebin = 0;
1665       }
1666
1667       adcwc = 0;
1668       Int_t nADCwords = (timebins + 2) / 3;
1669       AliDebug(3, Form("Now reading %i words for channel %2i", nADCwords, channelno));
1670       Int_t adccol = adccoloff - channelno;
1671       Int_t padcol = padcoloff - channelno;
1672 //      if (adccol < 3 || adccol > 165)
1673 //      AliInfo(Form("writing channel %i of det %3i %i:%2i to adcrow/-col: %i/%i padcol: %i",
1674 //                   channelno, fCurrHC/2, fCurrRobPos, fCurrMcmPos, row, adccol, padcol));
1675
1676       while ((adcwc < nADCwords) &&
1677              (*(fPayloadCurr) != fgkDataEndmarker) &&
1678              (fPayloadCurr - fPayloadStart < fPayloadSize)) {
1679         int check = 0x3 & *fPayloadCurr;
1680         if (channelno % 2 != 0) { // odd channel
1681           if (check != 0x2 && channelno < 21) {
1682             MCMError(kAdcCheckInvalid,
1683                      "%i for %2i. ADC word in odd channel %i",
1684                      check, adcwc+1, channelno);
1685           }
1686         }
1687         else {                  // even channel
1688           if (check != 0x3 && channelno < 21) {
1689             MCMError(kAdcCheckInvalid,
1690                      "%i for %2i. ADC word in even channel %i",
1691                      check, adcwc+1, channelno);
1692           }
1693         }
1694
1695         // filling the actual timebin data
1696         int tb2 = 0x3ff & *fPayloadCurr >> 22;
1697         int tb1 = 0x3ff & *fPayloadCurr >> 12;
1698         int tb0 = 0x3ff & *fPayloadCurr >> 2;
1699         if (adcwc != 0 || fCurrNtimebins <= 30)
1700           fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
1701         else
1702           tb0 = -1;
1703         if (currentTimebin >= fCurrNtimebins)
1704           break;
1705         fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
1706         if (currentTimebin >= fCurrNtimebins)
1707           break;
1708         fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
1709
1710         adcwc++;
1711         fPayloadCurr++;
1712       }
1713
1714       if (adcwc != nADCwords)
1715         MCMError(kAdcDataAbort);
1716
1717       // adding index
1718       if (padcol > 0 && padcol < 144) {
1719         fSignalIndex->AddIndexRC(row, padcol);
1720       }
1721
1722       channelcount++;
1723     }
1724
1725     if (fCurrSm > -1 && fCurrSm < 18) {
1726       fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNChannels += channelcount;
1727       fStats.fStatsSector[fCurrSm].fNChannels                      += channelcount;
1728     }
1729     if (channelcount != channelcountExp)
1730       MCMError(kAdcChannelsMiss);
1731
1732     mcmcount++;
1733     if (fCurrSm > -1 && fCurrSm < 18) {
1734       fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNMCMs++;
1735       fStats.fStatsSector[fCurrSm].fNMCMs++;
1736     }
1737
1738     if (IsDumping() && DumpingMCM(fCurrHC/2, fCurrRobPos, fCurrMcmPos)) {
1739       AliInfo(DumpRaw(Form("Event %i: Det %3i ROB %i MCM %2i", fRawReader->GetEventIndex(), fCurrHC/2, fCurrRobPos, fCurrMcmPos),
1740                       startPosMCM, fPayloadCurr - startPosMCM));
1741     }
1742
1743     // continue with next MCM
1744   }
1745
1746   // check for missing MCMs (if header suppression is inactive)
1747   if (((fCurrMajor & 0x1) == 0) && (mcmcount != mcmcountExp)) {
1748     LinkError(kMissMcmHeaders,
1749               "No. of MCM headers %i not as expected: %i",
1750               mcmcount, mcmcountExp);
1751   }
1752
1753   return (fPayloadCurr - start);
1754 }
1755
1756 Int_t AliTRDrawStream::ReadNonZSData()
1757 {
1758   // read the non-zs data from one link from the current reading position
1759
1760   UInt_t *start = fPayloadCurr;
1761
1762   Int_t mcmcount = 0;
1763   Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
1764   Int_t channelcount = 0;
1765   Int_t channelcountExp = 0;
1766   Int_t timebins;
1767   Int_t currentTimebin = 0;
1768   Int_t adcwc = 0;
1769   Int_t evno = -1;
1770   Int_t lastmcmpos = -1;
1771   Int_t lastrobpos = -1;
1772
1773   if (fCurrNtimebins != fNtimebins) {
1774     if (fNtimebins > 0)
1775       LinkError(kNtimebinsChanged,
1776                 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
1777     fNtimebins = fCurrNtimebins;
1778   }
1779
1780   timebins = fNtimebins;
1781
1782   while (*(fPayloadCurr) != fgkDataEndmarker &&
1783          fPayloadCurr - fPayloadStart < fPayloadSize - 2) {
1784
1785     // ----- Checking MCM Header -----
1786     AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
1787
1788     // ----- checking for proper readout order - ROB -----
1789     fCurrRobPos = ROB(*fPayloadCurr);
1790     if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1791       if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) > lastrobpos)
1792         lastmcmpos = -1;
1793       lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1794     }
1795     else {
1796       ROBError(kPosUnexp, Form("#%i after #%i in readout order", GetROBReadoutPos(ROB(*fPayloadCurr) / 2), lastrobpos));
1797     }
1798
1799     // ----- checking for proper readout order - MCM -----
1800     fCurrMcmPos = MCM(*fPayloadCurr);
1801     if (GetMCMReadoutPos(MCM(*fPayloadCurr)) > lastmcmpos) {
1802       lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
1803     }
1804     else {
1805       MCMError(kPosUnexp, Form("#%i after #%i in readout order", GetMCMReadoutPos(MCM(*fPayloadCurr)), lastmcmpos));
1806     }
1807
1808     if (EvNo(*fPayloadCurr) != (evno & 0xfffff)) {
1809       if (evno == -1) {
1810         evno = EvNo(*fPayloadCurr);
1811       }
1812       else {
1813         MCMError(kEvCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
1814       }
1815     }
1816
1817     channelcount = 0;
1818     channelcountExp = 21;
1819     int channelno = -1;
1820
1821     Int_t adccoloff = AdcColOffset(*fPayloadCurr);
1822     Int_t padcoloff = PadColOffset(*fPayloadCurr);
1823     Int_t row = Row(*fPayloadCurr);
1824
1825     fPayloadCurr++;
1826
1827     // ----- reading marked ADC channels -----
1828     while (channelcount < channelcountExp &&
1829            *(fPayloadCurr) != fgkDataEndmarker) {
1830       if (channelno < 20)
1831         channelno++;
1832
1833       currentTimebin = 0;
1834
1835       adcwc = 0;
1836       Int_t nADCwords = (timebins + 2) / 3;
1837       AliDebug(2, Form("Now looking %i words", nADCwords));
1838       Int_t adccol = adccoloff - channelno;
1839       Int_t padcol = padcoloff - channelno;
1840       while ((adcwc < nADCwords) &&
1841              (*(fPayloadCurr) != fgkDataEndmarker) &&
1842              (fPayloadCurr - fPayloadStart < fPayloadSize)) {
1843         int check = 0x3 & *fPayloadCurr;
1844         if (channelno % 2 != 0) { // odd channel
1845           if (check != 0x2 && channelno < 21) {
1846             MCMError(kAdcCheckInvalid,
1847                      "%i for %2i. ADC word in odd channel %i",
1848                      check, adcwc+1, channelno);
1849           }
1850         }
1851         else {                  // even channel
1852           if (check != 0x3 && channelno < 21) {
1853             MCMError(kAdcCheckInvalid,
1854                      "%i for %2i. ADC word in even channel %i",
1855                      check, adcwc+1, channelno);
1856           }
1857         }
1858
1859         // filling the actual timebin data
1860         int tb2 = 0x3ff & *fPayloadCurr >> 22;
1861         int tb1 = 0x3ff & *fPayloadCurr >> 12;
1862         int tb0 = 0x3ff & *fPayloadCurr >> 2;
1863         if (adcwc != 0 || fCurrNtimebins <= 30)
1864           fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
1865         else
1866           tb0 = -1;
1867         if (currentTimebin >= fCurrNtimebins)
1868           break;
1869         fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
1870         if (currentTimebin >= fCurrNtimebins)
1871           break;
1872         fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
1873
1874         adcwc++;
1875         fPayloadCurr++;
1876       }
1877
1878       if (adcwc != nADCwords)
1879         MCMError(kAdcDataAbort);
1880
1881       // adding index
1882       if (padcol > 0 && padcol < 144) {
1883         fSignalIndex->AddIndexRC(row, padcol);
1884       }
1885
1886       channelcount++;
1887     }
1888
1889     if (channelcount != channelcountExp)
1890       MCMError(kAdcChannelsMiss);
1891     mcmcount++;
1892     // continue with next MCM
1893   }
1894
1895   // check for missing MCMs (if header suppression is inactive)
1896   if (mcmcount != mcmcountExp) {
1897     LinkError(kMissMcmHeaders,
1898               "%i not as expected: %i", mcmcount, mcmcountExp);
1899   }
1900
1901   return (fPayloadCurr - start);
1902 }
1903
1904 Int_t AliTRDrawStream::SeekNextStack()
1905 {
1906   // proceed in raw data stream till the next stack
1907
1908   if (!fCurrStackEndmarkerAvail)
1909     return 0;
1910
1911   UInt_t *start = fPayloadCurr;
1912
1913   // read until data endmarkers
1914   while ((fPayloadCurr - fPayloadStart < fPayloadSize-1) &&
1915          ((fPayloadCurr[0] != fgkStackEndmarker[0]) ||
1916           (fPayloadCurr[1] != fgkStackEndmarker[1])))
1917     fPayloadCurr++;
1918
1919   if ((fPayloadCurr - start) != 0)
1920     StackError(kUnknown, "skipped %i words to reach stack endmarker", fPayloadCurr - start);
1921
1922   AliDebug(2, Form("stack endmarker: 0x%08x 0x%08x", fPayloadCurr[0], fPayloadCurr[1]));
1923
1924   // goto next stack
1925   fPayloadCurr++;
1926   fPayloadCurr++;
1927
1928   return (fPayloadCurr-start);
1929 }
1930
1931 Int_t AliTRDrawStream::SeekNextLink()
1932 {
1933   // proceed in raw data stream till the next link
1934
1935   UInt_t *start = fPayloadCurr;
1936
1937   // read until data endmarkers
1938   while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1939          ((fPayloadCurr[0] != fgkStackEndmarker[0]) ||
1940           (fPayloadCurr[1] != fgkStackEndmarker[1])) &&
1941          *fPayloadCurr != fgkDataEndmarker)
1942     fPayloadCurr++;
1943
1944   // read all data endmarkers
1945   while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1946          *fPayloadCurr == fgkDataEndmarker)
1947     fPayloadCurr++;
1948
1949   return (fPayloadCurr - start);
1950 }
1951
1952 Bool_t AliTRDrawStream::ConnectTracklets(TTree *trklTree)
1953 {
1954   // connect the tracklet tree used to store the tracklet output
1955
1956   fTrackletTree = trklTree;
1957   if (!fTrackletTree)
1958     return kTRUE;
1959
1960   if (!fTrackletTree->GetBranch("hc"))
1961     fTrackletTree->Branch("hc", &fCurrHC, "hc/I");
1962   else
1963     fTrackletTree->SetBranchAddress("hc", &fCurrHC);
1964
1965   if (!fTrackletTree->GetBranch("trkl"))
1966     fTrackletTree->Branch("trkl", &fTrackletArray);
1967   else
1968     fTrackletTree->SetBranchAddress("trkl", &fTrackletArray);
1969
1970   return kTRUE;
1971 }
1972
1973
1974 void AliTRDrawStream::EquipmentError(ErrorCode_t err, const char *const msg, ...)
1975 {
1976   // register error according to error code on equipment level
1977   // and return the corresponding error message
1978
1979   fLastError.fSector = fCurrEquipmentId - kDDLOffset;
1980   fLastError.fStack  = -1;
1981   fLastError.fLink   = -1;
1982   fLastError.fRob    = -1;
1983   fLastError.fMcm    = -1;
1984   fLastError.fError  = err;
1985   (this->*fStoreError)();
1986
1987   va_list ap;
1988   if (fgErrorDebugLevel[err] > 10)
1989     AliDebug(fgErrorDebugLevel[err],
1990              Form("Event %6i: Eq. %2d - %s : %s",
1991                   fRawReader->GetEventIndex(), fCurrEquipmentId, fgkErrorMessages[err],
1992                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1993   else
1994     AliError(Form("Event %6i: Eq. %2d - %s : %s",
1995                   fRawReader->GetEventIndex(), fCurrEquipmentId, fgkErrorMessages[err],
1996                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1997   fErrorFlags |= fgErrorBehav[err];
1998 }
1999
2000
2001 void AliTRDrawStream::StackError(ErrorCode_t err, const char *const msg, ...)
2002 {
2003   // register error according to error code on stack level
2004   // and return the corresponding error message
2005
2006   fLastError.fSector = fCurrEquipmentId - kDDLOffset;
2007   fLastError.fStack  = fCurrSlot;
2008   fLastError.fLink   = -1;
2009   fLastError.fRob    = -1;
2010   fLastError.fMcm    = -1;
2011   fLastError.fError  = err;
2012   (this->*fStoreError)();
2013
2014   va_list ap;
2015   if (fgErrorDebugLevel[err] > 0)
2016     AliDebug(fgErrorDebugLevel[err],
2017              Form("Event %6i: Eq. %2d S %i - %s : %s",
2018                   fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgkErrorMessages[err],
2019                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
2020   else
2021     AliError(Form("Event %6i: Eq. %2d S %i - %s : %s",
2022                   fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgkErrorMessages[err],
2023                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
2024   fErrorFlags |= fgErrorBehav[err];
2025 }
2026
2027
2028 void AliTRDrawStream::LinkError(ErrorCode_t err, const char *const msg, ...)
2029 {
2030   // register error according to error code on link level
2031   // and return the corresponding error message
2032
2033   fLastError.fSector = fCurrEquipmentId - kDDLOffset;
2034   fLastError.fStack  = fCurrSlot;
2035   fLastError.fLink   = fCurrLink;
2036   fLastError.fRob    = -1;
2037   fLastError.fMcm    = -1;
2038   fLastError.fError  = err;
2039   (this->*fStoreError)();
2040
2041   va_list ap;
2042   if (fgErrorDebugLevel[err] > 0)
2043     AliDebug(fgErrorDebugLevel[err],
2044              Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
2045                   fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgkErrorMessages[err],
2046                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
2047   else
2048     AliError(Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
2049                   fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgkErrorMessages[err],
2050                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
2051   fErrorFlags |= fgErrorBehav[err];
2052 }
2053
2054
2055 void AliTRDrawStream::ROBError(ErrorCode_t err, const char *const msg, ...)
2056 {
2057   // register error according to error code on ROB level
2058   // and return the corresponding error message
2059
2060   fLastError.fSector = fCurrEquipmentId - kDDLOffset;
2061   fLastError.fStack  = fCurrSlot;
2062   fLastError.fLink   = fCurrLink;
2063   fLastError.fRob    = fCurrRobPos;
2064   fLastError.fMcm    = -1;
2065   fLastError.fError  = err;
2066   (this->*fStoreError)();
2067
2068   va_list ap;
2069   if (fgErrorDebugLevel[err] > 0)
2070     AliDebug(fgErrorDebugLevel[err],
2071              Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : %s",
2072                   fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgkErrorMessages[err],
2073                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
2074   else
2075     AliError(Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : %s",
2076                   fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgkErrorMessages[err],
2077                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
2078   fErrorFlags |= fgErrorBehav[err];
2079 }
2080
2081
2082 void AliTRDrawStream::MCMError(ErrorCode_t err, const char *const msg, ...)
2083 {
2084   // register error according to error code on MCM level
2085   // and return the corresponding error message
2086
2087   fLastError.fSector = fCurrEquipmentId - kDDLOffset;
2088   fLastError.fStack  = fCurrSlot;
2089   fLastError.fLink   = fCurrLink;
2090   fLastError.fRob    = fCurrRobPos;
2091   fLastError.fMcm    = fCurrMcmPos;
2092   fLastError.fError  = err;
2093   (this->*fStoreError)();
2094
2095   va_list ap;
2096   if (fgErrorDebugLevel[err] > 0)
2097     AliDebug(fgErrorDebugLevel[err],
2098              Form("Event %6i: Eq. %2d S %i l %2i ROB %i MCM %2i - %s : %s",
2099                   fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fCurrMcmPos, fgkErrorMessages[err],
2100                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
2101   else
2102     AliError(Form("Event %6i: Eq. %2d S %i l %2i ROB %i MCM %2i - %s : %s",
2103                   fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fCurrMcmPos, fgkErrorMessages[err],
2104                   (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
2105   fErrorFlags |= fgErrorBehav[err];
2106 }
2107
2108 const char* AliTRDrawStream::GetErrorMessage(ErrorCode_t errCode)
2109 {
2110   // return the error message for the given error code
2111
2112   if (errCode > 0 && errCode < kLastErrorCode)
2113     return fgkErrorMessages[errCode];
2114   else
2115     return "";
2116 }
2117
2118 void AliTRDrawStream::AliTRDrawStats::ClearStats()
2119 {
2120   // clear statistics (includes clearing sector-wise statistics)
2121
2122   fBytesRead = 0;
2123   for (Int_t iSector = 0; iSector < 18; iSector++) {
2124     fStatsSector[iSector].ClearStats();
2125   }
2126
2127 }
2128
2129 void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::ClearStats()
2130 {
2131   // clear statistics (includes clearing HC-wise statistics)
2132
2133   fBytes = 0;
2134   fBytesRead = 0;
2135   fNTracklets = 0;
2136   fNMCMs = 0;
2137   fNChannels = 0;
2138
2139   for (Int_t iHC = 0; iHC < 60; iHC++) {
2140     fStatsHC[iHC].ClearStats();
2141   }
2142 }
2143
2144 void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::AliTRDrawStatsHC::ClearStats()
2145 {
2146   // clear statistics
2147
2148   fBytes = 0;
2149   fBytesRead = 0;
2150   fNTracklets = 0;
2151   fNMCMs = 0;
2152   fNChannels = 0;
2153 }
2154
2155 void AliTRDrawStream::SetDumpMCM(Int_t det, Int_t rob, Int_t mcm, Bool_t dump)
2156 {
2157   // mark MCM for dumping of raw data
2158
2159   if (dump) {
2160     fDumpMCM[fNDumpMCMs++] = (det << 7) | (rob << 4) | mcm;
2161   }
2162   else {
2163     Int_t iMCM;
2164     for (iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
2165       if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
2166         fNDumpMCMs--;
2167         break;
2168       }
2169     }
2170     for ( ; iMCM < fNDumpMCMs; iMCM++) {
2171       fDumpMCM[iMCM] = fDumpMCM[iMCM+1];
2172     }
2173   }
2174 }
2175
2176 Bool_t AliTRDrawStream::DumpingMCM(Int_t det, Int_t rob, Int_t mcm)  const
2177 {
2178   // check if MCM data should be dumped
2179
2180   for (Int_t iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
2181     if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
2182       return kTRUE;
2183     }
2184   }
2185   return kFALSE;
2186 }
2187
2188 TString AliTRDrawStream::DumpRaw(TString title, const UInt_t *start, Int_t length, UInt_t endmarker)
2189 {
2190   // dump raw data
2191
2192   title += "\n";
2193   for (Int_t pos = 0; pos < length; pos += 4) {
2194     if ((start[pos+0] != endmarker) && pos+0 < length)
2195       if ((start[pos+1] != endmarker && pos+1 < length))
2196         if ((start[pos+2] != endmarker && pos+2 < length))
2197           if ((start[pos+3] != endmarker && pos+3 < length))
2198             title += Form("   0x%08x 0x%08x 0x%08x 0x%08x\n",
2199                           start[pos+0], start[pos+1], start[pos+2], start[pos+3]);
2200           else {
2201             title += Form("   0x%08x 0x%08x 0x%08x 0x%08x\n",
2202                           start[pos+0], start[pos+1], start[pos+2], start[pos+3]);
2203             return title;
2204           }
2205         else {
2206           title += Form("   0x%08x 0x%08x 0x%08x\n",
2207                         start[pos+0], start[pos+1], start[pos+2]);
2208           return title;
2209         }
2210       else {
2211         title += Form("   0x%08x 0x%08x\n",
2212                       start[pos+0], start[pos+1]);
2213         return title;
2214       }
2215     else {
2216       title += Form("   0x%08x\n",
2217                     start[pos+0]);
2218       return title;
2219     }
2220   }
2221   return title;
2222 }
2223
2224 TString AliTRDrawStream::DumpMcmHeader(TString title, UInt_t word)
2225 {
2226   title += Form("0x%08x -> ROB: %i, MCM: %2i",
2227                 word, ROB(word), MCM(word));
2228   return title;
2229 }
2230
2231 TString AliTRDrawStream::DumpAdcMask(TString title, UInt_t word)
2232 {
2233   title += Form("0x%08x -> #ch : %2i, 0x%06x (%2i ch)",
2234                 word, GetNActiveChannels(word), GetActiveChannels(word), GetNActiveChannelsFromMask(word));
2235   return title;
2236 }
2237
2238 AliTRDrawStream::AliTRDrawStreamError::AliTRDrawStreamError(Int_t error, Int_t sector, Int_t stack, Int_t link, Int_t rob, Int_t mcm) :
2239   fError(error),
2240   fSector(sector),
2241   fStack(stack),
2242   fLink(link),
2243   fRob(rob),
2244   fMcm(mcm)
2245 {
2246   // ctor
2247
2248 }
2249
2250 void AliTRDrawStream::SortTracklets(TClonesArray *trklArray, TList &sortedTracklets, Int_t *indices)
2251 {
2252   // sort tracklets for referencing from GTU tracks
2253
2254   if (!trklArray)
2255     return;
2256
2257   Int_t nTracklets = trklArray->GetEntriesFast();
2258
2259   Int_t lastHC = -1;
2260   for (Int_t iTracklet = 0; iTracklet < nTracklets; iTracklet++) {
2261     AliTRDtrackletBase *trkl = (AliTRDtrackletBase*) ((*trklArray)[iTracklet]);
2262     Int_t hc = trkl->GetHCId();
2263     if ((hc < 0) || (hc >= 1080)) {
2264       AliErrorClass(Form("HC for tracklet: 0x%08x out of range: %i", trkl->GetTrackletWord(), trkl->GetHCId()));
2265       continue;
2266     }
2267     AliDebugClass(5, Form("hc: %4i : 0x%08x z: %2i", hc, trkl->GetTrackletWord(), trkl->GetZbin()));
2268     if (hc != lastHC) {
2269       AliDebugClass(2, Form("set tracklet index for HC %i to %i", hc, iTracklet));
2270       indices[hc] = iTracklet + 1;
2271       lastHC = hc;
2272     }
2273   }
2274
2275   for (Int_t iDet = 0; iDet < 540; iDet++) {
2276     Int_t trklIndexA = indices[2*iDet + 0] - 1;
2277     Int_t trklIndexB = indices[2*iDet + 1] - 1;
2278     Int_t trklIndex  = sortedTracklets.GetEntries();
2279     AliTRDtrackletBase *trklA = trklIndexA > -1 ? (AliTRDtrackletBase*) ((*trklArray)[trklIndexA]) : 0x0;
2280     AliTRDtrackletBase *trklB = trklIndexB > -1 ? (AliTRDtrackletBase*) ((*trklArray)[trklIndexB]) : 0x0;
2281     AliTRDtrackletBase *trklNext = 0x0;
2282     while (trklA != 0x0 || trklB != 0x0) {
2283       AliDebugClass(5, Form("det %i - A: %i/%i -> %p, B: %i/%i -> %p",
2284                        iDet, trklIndexA, nTracklets, trklA, trklIndexB, nTracklets, trklB));
2285       if (trklA == 0x0) {
2286         trklNext = trklB;
2287         trklIndexB++;
2288         trklB = trklIndexB < nTracklets ? (AliTRDtrackletBase*) ((*trklArray)[trklIndexB]) : 0x0;
2289         if (trklB && trklB->GetHCId() != 2*iDet + 1)
2290           trklB = 0x0;
2291       }
2292       else if (trklB == 0x0) {
2293         trklNext = trklA;
2294         trklIndexA++;
2295         trklA = trklIndexA < nTracklets ? (AliTRDtrackletBase*) ((*trklArray)[trklIndexA]) : 0x0;
2296         if (trklA && trklA->GetHCId() != 2*iDet)
2297           trklA = 0x0;
2298       }
2299       else {
2300         if (trklA->GetZbin() <= trklB->GetZbin()) {
2301           trklNext = trklA;
2302           trklIndexA++;
2303           trklA = trklIndexA < nTracklets ? (AliTRDtrackletBase*) ((*trklArray)[trklIndexA]) : 0x0;
2304           if (trklA && trklA->GetHCId() != 2*iDet)
2305             trklA = 0x0;
2306         }
2307         else {
2308           trklNext = trklB;
2309           trklIndexB++;
2310           trklB = trklIndexB < nTracklets ? (AliTRDtrackletBase*) ((*trklArray)[trklIndexB]) : 0x0;
2311           if (trklB && trklB->GetHCId() != 2*iDet + 1)
2312             trklB = 0x0;
2313         }
2314       }
2315       if (trklNext) {
2316         sortedTracklets.Add(trklNext);
2317
2318       }
2319     }
2320
2321     // updating tracklet indices as in output
2322     if (sortedTracklets.GetEntries() != trklIndex) {
2323       indices[2*iDet + 0] = trklIndex;
2324       indices[2*iDet + 1] = sortedTracklets.GetEntries();
2325     } else {
2326       indices[2*iDet + 0] = indices[2*iDet + 1] = -1;
2327     }
2328   }
2329 }
2330
2331 void AliTRDrawStream::AssignTracklets(AliESDTrdTrack *trdTrack, Int_t *trackletIndex, Int_t refIndex[6])
2332 {
2333   UInt_t mask  = trdTrack->GetLayerMask();
2334   UInt_t stack = trdTrack->GetStack();
2335
2336   for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2337     refIndex[iLayer] = -1;
2338
2339     if (mask & (1 << iLayer)) {
2340
2341       Int_t det = trdTrack->GetSector()*30 + stack*6 + iLayer;
2342       Int_t idx = trdTrack->GetTrackletIndex(iLayer);
2343
2344       if ((det < 0) || (det > 539)) {
2345         AliErrorClass(Form("Invalid detector no. from track: %i", 2*det));
2346         continue;
2347       }
2348       if (trackletIndex[2*det] >= 0) {
2349         if ((trackletIndex[2*det] + idx > -1) &&
2350             (trackletIndex[2*det] + idx < trackletIndex[2*det+1])) {
2351           refIndex[iLayer] = trackletIndex[2*det] + idx;
2352         } else {
2353           AliErrorClass(Form("Requested tracklet index %i out of range", idx));
2354         }
2355       } else {
2356         AliErrorClass(Form("Non-existing tracklets requested in det %i", det));
2357       }
2358     }
2359   }
2360 }