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