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