1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 ////////////////////////////////////////////////////////////////////////////
18 // Decoding data from the TRD raw stream //
19 // and translation into ADC values, on-line tracklets and tracks //
21 // Author: J. Klein (jochen.klein@cern.ch) //
23 ////////////////////////////////////////////////////////////////////////////
28 #include "TClonesArray.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"
44 #include "AliTRDrawStream.h"
47 #include "AliRunLoader.h"
49 ClassImp(AliTRDrawStream)
51 // some static information
52 Int_t AliTRDrawStream::fgMcmOrder[] = {12, 13, 14, 15,
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;
64 const char* AliTRDrawStream::fgkErrorMessages[] = {
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",
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",
81 "Missing expected ADC channels",
82 "Missing MCM headers",
86 Int_t AliTRDrawStream::fgErrorDebugLevel[] = {
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
130 AliTRDrawStream::AliTRDrawStream(AliRawReader *rawReader) :
131 fStoreError(&AliTRDrawStream::ForgetError),
132 fRawReader(rawReader),
149 fCurrSmHeaderSize(0),
150 fCurrSmHeaderVersion(0),
151 fCurrTrailerReadout(0),
152 fCurrTrgHeaderAvail(0),
153 fCurrTrgHeaderReadout(0),
154 fCurrTrkHeaderAvail(0),
156 fCurrTriggerEnable(0),
157 fCurrTriggerFired(0),
159 fCurrTrackletEnable(0),
161 fCurrTrkHeaderIndexWord(0x0),
162 fCurrTrkHeaderSize(0x0),
164 fCurrTrgHeaderIndexWord(0x0),
165 fCurrTrgHeaderSize(0x0),
167 fCurrStackIndexWord(0x0),
168 fCurrStackHeaderSize(0x0),
169 fCurrStackHeaderVersion(0x0),
171 fCurrCleanCheckout(0x0),
175 fCurrLinkMonitorFlags(0x0),
176 fCurrLinkDataTypeFlags(0x0),
177 fCurrLinkDebugFlags(0x0),
201 // default constructor
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++)
224 // preparing TClonesArray
225 fTrackletArray = new TClonesArray("AliTRDtrackletWord", 256);
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++) {
238 AliTRDrawStream::~AliTRDrawStream()
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;
262 Bool_t AliTRDrawStream::ReadEvent(TTree *trackletTree)
264 // read the current event from the raw reader and fill it to the digits manager
267 AliError("No raw reader available");
272 ConnectTracklets(trackletTree);
277 // loop over all DDLs
278 // data starts with GTU payload, i.e. SM index word
279 UChar_t *buffer = 0x0;
281 while (fRawReader->ReadNextData(buffer)) {
283 fCurrEquipmentId = fRawReader->GetEquipmentId();
284 AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
286 if (fCurrEquipmentId < kDDLOffset || fCurrEquipmentId > kDDLMax) {
287 EquipmentError(kNonTrdEq, "Skipping");
292 new ((*fMarkers)[fMarkers->GetEntriesFast()])
293 AliTRDrawStreamError(-kSecactive, fCurrEquipmentId - kDDLOffset);
295 ReadGTUHeaders((UInt_t*) buffer);
297 if (fCurrTrailerReadout)
300 // loop over all active links
301 AliDebug(2, Form("Stack mask 0x%02x", fCurrStackMask));
302 for (Int_t iStack = 0; iStack < fgkNstacks; iStack++) {
304 if ((fCurrStackMask & (1 << fCurrSlot)) == 0)
307 AliDebug(2, Form("Stack %i, Link mask: 0x%02x", fCurrSlot, fCurrLinkMask[fCurrSlot]));
308 for (Int_t iLink = 0; iLink < fgkNlinks; iLink++) {
310 fCurrHC = (fCurrEquipmentId - kDDLOffset) * fgkNstacks * fgkNlinks +
311 fCurrSlot * fgkNlinks + iLink;
312 if ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) == 0)
316 // check for link monitor error flag
317 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0) {
318 LinkError(kLinkMonitor);
319 if (fgErrorBehav[kLinkMonitor] == kTolerate)
323 // read the data from one HC
326 // read all data endmarkers
336 Bool_t AliTRDrawStream::NextDDL()
338 // continue reading with the next equipment
343 fCurrEquipmentId = 0;
347 UChar_t *buffer = 0x0;
349 while (fRawReader->ReadNextData(buffer)) {
351 fCurrEquipmentId = fRawReader->GetEquipmentId();
352 AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
354 if (fCurrEquipmentId < kDDLOffset || fCurrEquipmentId > kDDLMax) {
355 EquipmentError(kNonTrdEq, "Skipping");
360 new ((*fMarkers)[fMarkers->GetEntriesFast()])
361 AliTRDrawStreamError(-kSecactive, fCurrEquipmentId - kDDLOffset);
363 ReadGTUHeaders((UInt_t*) buffer);
365 if (fCurrTrailerReadout)
375 Int_t AliTRDrawStream::NextChamber(AliTRDdigitsManager *digMgr)
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
381 fDigitsManager = digMgr;
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;
392 AliTreeLoader *trklTreeLoader = (AliTreeLoader*) trklLoader->GetBaseLoader("tracklets-raw");
394 trklTree = trklTreeLoader->Tree();
396 trklTree = trklLoader->Tree();
399 if (fTrackletTree != trklTree)
400 ConnectTracklets(trklTree);
403 AliError("No raw reader available");
407 while (fCurrSlot < 0 || fCurrSlot >= fgkNstacks) {
412 while ((fCurrSlot < fgkNstacks) &&
413 (((fCurrStackMask & (1 << fCurrSlot)) == 0) ||
414 ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink))) == 0)) {
416 if (fCurrLink >= fgkNlinks) {
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;
427 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0) {
428 LinkError(kLinkMonitor);
429 if (fgErrorBehav[kLinkMonitor] == kTolerate)
433 // read the data from one HC
436 // read all data endmarkers
439 if (fCurrLink % 2 == 0) {
440 // if we just read the A-side HC then also check the B-side
443 if (fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) {
444 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0) {
445 LinkError(kLinkMonitor);
446 if (fgErrorBehav[kLinkMonitor] == kTolerate)
459 if (fCurrLink >= fgkNlinks) {
463 } while ((fCurrSlot < fgkNstacks) &&
464 (((fCurrStackMask & (1 << fCurrSlot)) == 0) ||
465 ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink))) == 0));
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);
472 return (fCurrSm * fgkNstacks*fgkNlinks/2 + fCurrStack * fgkNlinks/2 + fCurrLayer);
476 Int_t AliTRDrawStream::ReadGTUHeaders(UInt_t *buffer)
478 // check the data source and read the headers
480 if (fCurrEquipmentId >= kDDLOffset && fCurrEquipmentId <= kDDLMax) {
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()));
489 AliDebug(1, DumpRaw("raw data", fPayloadCurr, TMath::Min(fPayloadSize, 1000)));
492 if (ReadSmHeader() < 0) {
493 AliError(Form("Reading SM header failed, skipping this DDL %i", fCurrEquipmentId));
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);
505 // read trigger header(s) (if available)
506 if (fCurrTrgHeaderAvail)
507 ReadTriggerHeaders();
510 for (Int_t iStack = 0; iStack < fgkNstacks; iStack++) {
511 if ((fCurrStackMask & (1 << iStack)) != 0)
512 ReadStackHeader(iStack);
521 Int_t AliTRDrawStream::ReadSmHeader()
523 // read the SMU index header at the current reading position
524 // and store the information in the corresponding variables
526 if (fPayloadCurr - fPayloadStart >= fPayloadSize - 1) {
527 EquipmentError(kUnknown, "SM Header incomplete");
531 fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] = 0;
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;
540 switch (fCurrSmHeaderVersion) {
542 fCurrTrailerReadout = 0;
543 fCurrTrgHeaderAvail = 0;
545 fCurrTrkHeaderAvail = 0;
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;
562 AliError(Form("unknown SM header version: 0x%x", fCurrSmHeaderVersion));
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",
567 fCurrSmHeaderVersion,
573 fCurrTrkHeaderAvail ));
575 // jump to the first word after the SM header
576 fPayloadCurr += fCurrSmHeaderSize + 1;
578 return fCurrSmHeaderSize + 1;
581 Int_t AliTRDrawStream::DecodeGTUtracks()
583 // decode GTU track words
584 // this depends on the hardware revision of the SMU
586 Int_t sector = fCurrEquipmentId-kDDLOffset;
588 if ((sector < 0) || (sector > 17)) {
589 AliError(Form("Invalid sector %i for GTU tracks", sector));
593 AliDebug(1, DumpRaw(Form("GTU tracks in sector %2i (hw rev %i)", sector, fCurrHwRev),
594 fPayloadCurr + 4, 10, 0xffe0ffff));
596 fCurrTrgFlags[sector] = 0;
598 if (fCurrHwRev < 1772) {
599 UInt_t fastWord; // fast trigger word
600 ULong64_t trackWord = 0; // extended track word
603 for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
604 if (fPayloadCurr[iWord] == 0x10000000) { // stack boundary marker
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));
616 else if ((idx & 0x1) == 0x1) {
617 trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
618 AliDebug(1,Form("track debug word: 0x%016llx", trackWord));
620 AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
623 trk->SetSector(sector);
624 trk->SetStack((trackWord >> 60) & 0x7);
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);
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));
646 trackWord = fPayloadCurr[iWord];
652 else if (fCurrHwRev < 1804) {
653 UInt_t fastWord; // fast trigger word
654 ULong64_t trackWord = 0; // extended track word
657 for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
658 if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
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));
670 else if ((idx & 0x1) == 0x1) {
671 trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
672 AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
674 AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
677 trk->SetSector(fCurrEquipmentId-kDDLOffset);
678 trk->SetStack((trackWord >> 60) & 0x7);
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);
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));
700 trackWord = fPayloadCurr[iWord];
706 else if (fCurrHwRev < 1819) {
707 UInt_t fastWord; // fast trigger word
708 ULong64_t trackWord = 0; // extended track word
711 for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
712 if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
718 ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
719 fastWord = fPayloadCurr[iWord];
720 AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fastWord));
723 else if ((idx & 0x1) == 0x1) {
724 trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
725 AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
728 AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
731 trk->SetSector(fCurrEquipmentId-kDDLOffset);
732 trk->SetStack((trackWord >> 60) & 0x7);
735 // trk->SetPt(((trackWord & 0xffff) ^ 0x8000) - 0x8000);
737 trk->SetLayerMask((trackWord >> 16) & 0x3f);
738 trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
739 trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
740 trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
741 trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
742 trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
743 trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
748 Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
749 if (TMath::Abs(pt) > 0.1) {
750 trk->SetA((Int_t) (0.15*51625./100./trk->Pt() / 160e-4 * 2));
755 trackWord = fPayloadCurr[iWord];
761 else if (fCurrHwRev < 1860) {
762 UInt_t fastWord; // fast trigger word
763 ULong64_t trackWord = 0; // extended track word
766 Bool_t upperWord = kFALSE;
768 for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
769 if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
775 // assemble the 32-bit words out of 16-bit blocks
777 word |= (fPayloadCurr[iWord] & 0xffff0000);
781 // lower word is read first
782 word = (fPayloadCurr[iWord] & 0xffff0000) >> 16;
787 if ((word & 0xffff0008) == 0x13370008) {
789 AliDebug(1, Form("stack %i: fast track word: 0x%08x", stack, fastWord));
790 if (fastWord & (1 << 13))
791 fCurrTrgFlags[sector] |= 1 << (stack+11);
794 else if ((idx & 0x1) == 0x1) {
795 trackWord |= ((ULong64_t) word) << 32;
796 AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
798 AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
801 trk->SetSector(fCurrEquipmentId-kDDLOffset);
802 trk->SetStack((trackWord >> 60) & 0x7);
806 trk->SetLayerMask((trackWord >> 16) & 0x3f);
807 trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
808 trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
809 trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
810 trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
811 trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
812 trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
817 Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
818 if (TMath::Abs(pt) > 0.1) {
819 trk->SetA((Int_t) (0.15*51625./100./pt / 160e-4 * 2));
832 ULong64_t trackWord = 0; // this is the debug word
835 Bool_t upperWord = kFALSE;
837 for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
838 if (fPayloadCurr[iWord] == 0xffe0ffff) {
844 // assemble the 32-bit words out of 16-bit blocks
846 word |= (fPayloadCurr[iWord] & 0xffff0000);
850 // lower word is read first
851 word = (fPayloadCurr[iWord] & 0xffff0000) >> 16;
856 if ((word & 0xffff0008) == 0x13370008) {
857 AliDebug(1, Form("stack %i: fast track word: 0x%08x", stack, word));
860 else if ((word & 0xffff0010) == 0x13370010) {
861 AliDebug(1, Form("stack %i: tracking done word: 0x%08x", stack, word));
862 fCurrTrgFlags[sector] |= 1 << (stack+11);
865 else if ((idx & 0x1) == 0x1) {
866 trackWord |= ((ULong64_t) word) << 32;
867 AliDebug(1, Form("track debug word: 0x%16llx", trackWord));
869 AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
871 trk->SetSector(fCurrEquipmentId-kDDLOffset);
872 trk->SetStack((trackWord >> 60) & 0x7);
876 trk->SetLayerMask((trackWord >> 16) & 0x3f);
877 trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
878 trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
879 trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
880 trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
881 trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
882 trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
887 Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
888 if (TMath::Abs(pt) > 0.1) {
889 trk->SetA(-(Int_t) (0.15*51625./100./pt / 160e-4 * 2));
903 Int_t AliTRDrawStream::ReadTrackingHeader(Int_t stack)
905 // read the tracking information and store it for the given stack
909 fCurrTrkHeaderIndexWord[stack] = *fPayloadCurr;
910 fCurrTrkHeaderSize[stack] = ((*fPayloadCurr) >> 16) & 0x3ff;
912 AliDebug(1, Form("tracking header index word: 0x%08x, size: %i (hw rev: %i)",
913 fCurrTrkHeaderIndexWord[stack], fCurrTrkHeaderSize[stack], fCurrHwRev));
914 Int_t trackingTime = *fPayloadCurr & 0x3ff;
916 fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] |= ((fCurrTrkHeaderIndexWord[stack] >> 10) & 0x1) << (22 + stack);
920 ULong64_t trackWord = 0;
922 Int_t trackIndex = fTracks ? fTracks->GetEntriesFast() : -1;
924 for (UInt_t iWord = 0; iWord < fCurrTrkHeaderSize[stack]; iWord++) {
927 // first part of 64-bit word
928 trackWord = fPayloadCurr[iWord];
931 trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
933 if (trackWord & (1ul << 63)) {
934 if ((trackWord & (0x3ful << 56)) != 0) {
936 AliDebug(2, Form("track word: 0x%016llx", trackWord));
939 AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
942 trk->SetSector(fCurrEquipmentId-kDDLOffset);
943 trk->SetLayerMask((trackWord >> 56) & 0x3f);
944 trk->SetA( (((trackWord >> 38) & 0x3ffff) ^ 0x20000) - 0x20000);
945 trk->SetB( (((trackWord >> 20) & 0x3ffff) ^ 0x20000) - 0x20000);
946 trk->SetC( (((trackWord >> 8) & 0xffff) ^ 0x8000) - 0x8000);
947 trk->SetPID((trackWord >> 0) & 0xff);
948 trk->SetStack(stack);
950 // now compare the track word with the one generated from the ESD information
951 if (trackWord != trk->GetTrackWord(0)) {
952 AliError(Form("track word 0x%016llx does not match the read one 0x%016llx",
953 trk->GetTrackWord(0), trackWord));
958 // done marker (so far only used to set trigger flag)
959 fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] |= 1 << (27 + stack);
960 fCurrTrkFlags[(fCurrEquipmentId-kDDLOffset)*fgkNstacks + fCurrStack] = trackWord;
962 AliDebug(2, Form("tracking done marker: 0x%016llx, trigger flags: 0x%08x",
963 trackWord, fCurrTrgFlags[fCurrEquipmentId-kDDLOffset]));
964 AliDebug(2, Form("seg / stack / first / last / done / index : %i %i %lli %lli %lli %i",
965 fCurrEquipmentId - kDDLOffset, stack,
966 (trackWord >> 20) & 0x3ff,
967 (trackWord >> 10) & 0x3ff,
968 (trackWord >> 0) & 0x3ff,
973 // extended track word
974 AliDebug(2, Form("extended track word: 0x%016llx", trackWord));
977 AliESDTrdTrack *trk = (AliESDTrdTrack*) (*fTracks)[trackIndex];
979 trk->SetFlags((trackWord >> 52) & 0x7ff);
980 trk->SetReserved((trackWord >> 49) & 0x7);
981 trk->SetY((trackWord >> 36) & 0x1fff);
982 trk->SetTrackletIndex((trackWord >> 0) & 0x3f, 0);
983 trk->SetTrackletIndex((trackWord >> 6) & 0x3f, 1);
984 trk->SetTrackletIndex((trackWord >> 12) & 0x3f, 2);
985 trk->SetTrackletIndex((trackWord >> 18) & 0x3f, 3);
986 trk->SetTrackletIndex((trackWord >> 24) & 0x3f, 4);
987 trk->SetTrackletIndex((trackWord >> 30) & 0x3f, 5);
989 if (trackWord != trk->GetExtendedTrackWord(0)) {
990 AliError(Form("extended track word 0x%016llx does not match the read one 0x%016llx",
991 trk->GetExtendedTrackWord(0), trackWord));
1001 fPayloadCurr += fCurrTrkHeaderSize[stack];
1003 return fCurrTrkHeaderSize[stack];
1006 Int_t AliTRDrawStream::ReadTriggerHeaders()
1008 // read all trigger headers present
1010 AliDebug(1, Form("trigger mask: 0x%03x, fired: 0x%03x\n",
1011 fCurrTriggerEnable, fCurrTriggerFired));
1012 // loop over potential trigger blocks
1013 for (Int_t iTrigger = 0; iTrigger < fgkNtriggers; iTrigger++) {
1014 // check for trigger enable
1015 if (fCurrTriggerEnable & (1 << iTrigger)) {
1016 // check for readout mode and trigger fired
1017 if ((fCurrTrgHeaderReadout == 0) || (fCurrTriggerFired & (1 << iTrigger))) {
1019 AliDebug(1, Form("trigger index word %i: 0x%08x\n", iTrigger, *fPayloadCurr));
1020 fCurrTrgHeaderIndexWord[iTrigger] = *fPayloadCurr;
1021 fCurrTrgHeaderSize[iTrigger] = ((*fPayloadCurr) >> 16) & 0x3ff;
1022 if (iTrigger == 7) {
1023 // timeout trigger, use to extract tracking time
1024 fCurrTrgFlags[fCurrEquipmentId-kDDLOffset] |= (*fPayloadCurr & 0x3ff) << 12;
1029 fPayloadCurr += fCurrTrgHeaderSize[iTrigger];
1037 Int_t AliTRDrawStream::ReadStackHeader(Int_t stack)
1039 // read the stack header
1040 // and store the information in the corresponding variables
1042 fCurrStackIndexWord[stack] = *fPayloadCurr;
1043 fCurrStackHeaderSize[stack] = (((*fPayloadCurr) >> 16) & 0xffff) + 1;
1044 fCurrStackHeaderVersion[stack] = ((*fPayloadCurr) >> 12) & 0xf;
1045 fCurrLinkMask[stack] = (*fPayloadCurr) & 0xfff;
1047 // dumping stack header
1048 AliDebug(1, DumpRaw(Form("stack %i header", stack), fPayloadCurr, fCurrStackHeaderSize[stack]));
1050 if (fPayloadCurr - fPayloadStart >= fPayloadSize - (Int_t) fCurrStackHeaderSize[stack]) {
1051 LinkError(kStackHeaderInvalid, "Stack index header aborted");
1055 switch (fCurrStackHeaderVersion[stack]) {
1057 if (fCurrStackHeaderSize[stack] < 8) {
1058 LinkError(kStackHeaderInvalid, "Stack header smaller than expected!");
1062 fCurrCleanCheckout[stack] = fPayloadCurr[1] & 0x1;
1063 fCurrBoardId[stack] = (fPayloadCurr[1] >> 8) & 0xff;
1064 fCurrHwRevTMU[stack] = (fPayloadCurr[1] >> 16) & 0xffff;
1066 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1068 fCurrLinkMonitorFlags [stack*fgkNlinks + iLayer*2] = fPayloadCurr[iLayer+2] & 0xf;
1069 fCurrLinkDataTypeFlags [stack*fgkNlinks + iLayer*2] = (fPayloadCurr[iLayer+2] >> 4) & 0x3;
1070 fCurrLinkDebugFlags [stack*fgkNlinks + iLayer*2] = (fPayloadCurr[iLayer+2] >> 12) & 0xf;
1072 fCurrLinkMonitorFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 16) & 0xf;
1073 fCurrLinkDataTypeFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 20) & 0x3;
1074 fCurrLinkDebugFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 28) & 0xf;
1079 LinkError(kStackHeaderInvalid, "Invalid Stack Header version %x", fCurrStackHeaderVersion[stack]);
1082 fPayloadCurr += fCurrStackHeaderSize[stack];
1084 return fCurrStackHeaderSize[stack];
1087 Int_t AliTRDrawStream::ReadGTUTrailer()
1089 // read the SM trailer containing CRCs from various stages
1091 UInt_t* trailer = fPayloadStart + fPayloadSize -1;
1093 // look for the trailer index word from the end
1094 for (Int_t iWord = 0; iWord < fPayloadSize; iWord++) {
1095 if ((fPayloadStart[fPayloadSize-1-iWord] & 0xffff) == 0x1f51) {
1096 trailer = fPayloadStart + fPayloadSize - 1 - iWord;
1101 if (((*trailer) & 0xffff) == 0x1f51) {
1102 UInt_t trailerIndexWord = (*trailer);
1103 Int_t trailerSize = (trailerIndexWord >> 16) & 0xffff;
1104 AliDebug(2, DumpRaw("GTU trailer", trailer, trailerSize+1));
1105 // parse the trailer
1108 EquipmentError(kUnknown, "trailer index marker mismatch");
1113 Int_t AliTRDrawStream::ReadLinkData()
1115 // read the data in one link (one HC) until the data endmarker is reached
1116 // returns the number of words read!
1119 UInt_t* startPosLink = fPayloadCurr;
1121 AliDebug(1, DumpRaw(Form("link data from seg %2i slot %i link %2i", fCurrEquipmentId-kDDLOffset, fCurrSlot, fCurrLink),
1122 fPayloadCurr, TMath::Min((Int_t) (fPayloadSize - (fPayloadCurr-fPayloadStart)), 100), 0x00000000));
1125 new ((*fMarkers)[fMarkers->GetEntriesFast()])
1126 AliTRDrawStreamError(-kHCactive, fCurrEquipmentId-kDDLOffset, fCurrStack, fCurrLink);
1128 if (fErrorFlags & kDiscardHC)
1131 //??? add check whether tracklets are enabled
1132 count += ReadTracklets();
1133 if (fErrorFlags & kDiscardHC)
1136 AliDebug(1, DumpRaw("HC header", fPayloadCurr, 4, 0x00000000));
1137 count += ReadHcHeader();
1138 if (fErrorFlags & kDiscardHC)
1141 Int_t det = fCurrSm * 30 + fCurrStack * 6 + fCurrLayer;
1143 if (det > -1 && det < 540) {
1145 // ----- check which kind of data -----
1146 if (fCurrMajor & 0x40) {
1147 if ((fCurrMajor & 0x7) == 0x7) {
1148 AliDebug(1, "This is a config event");
1149 UInt_t *startPos = fPayloadCurr;
1150 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1151 *fPayloadCurr != fgkDataEndmarker)
1153 count += fPayloadCurr - startPos;
1155 // feeding TRAP config
1156 AliTRDtrapConfig *trapcfg = AliTRDtrapConfig::Instance();
1157 trapcfg->ReadPackedConfig(fCurrHC, startPos, fPayloadCurr - startPos);
1160 Int_t tpmode = fCurrMajor & 0x7;
1161 AliDebug(1, Form("Checking testpattern (mode %i) data", tpmode));
1166 // reading real data
1167 if (fDigitsManager) {
1168 if ((fAdcArray = fDigitsManager->GetDigits(det))) {
1169 //fAdcArray->Expand();
1170 if (fAdcArray->GetNtime() != fCurrNtimebins)
1171 fAdcArray->Allocate(16, 144, fCurrNtimebins);
1174 LinkError(kNoDigits);
1177 if (!fDigitsParam) {
1178 fDigitsParam = fDigitsManager->GetDigitsParam();
1181 fDigitsParam->SetPretriggerPhase(det, fCurrPtrgPhase);
1182 fDigitsParam->SetNTimeBins(det, fCurrNtimebins);
1183 fDigitsParam->SetADCbaseline(det, 10);
1186 if (fDigitsManager->UsesDictionaries()) {
1187 fDigitsManager->GetDictionary(det, 0)->Reset();
1188 fDigitsManager->GetDictionary(det, 1)->Reset();
1189 fDigitsManager->GetDictionary(det, 2)->Reset();
1192 if ((fSignalIndex = fDigitsManager->GetIndexes(det))) {
1193 fSignalIndex->SetSM(fCurrSm);
1194 fSignalIndex->SetStack(fCurrStack);
1195 fSignalIndex->SetLayer(fCurrLayer);
1196 fSignalIndex->SetDetNumber(det);
1197 if (!fSignalIndex->IsAllocated())
1198 fSignalIndex->Allocate(16, 144, fCurrNtimebins);
1201 if (fCurrMajor & 0x20) {
1202 AliDebug(1, "This is a zs event");
1203 count += ReadZSData();
1206 AliDebug(1, "This is a nozs event");
1207 count += ReadNonZSData();
1211 // just read until data endmarkers
1212 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1213 *fPayloadCurr != fgkDataEndmarker)
1219 LinkError(kInvalidDetector, "%i", det);
1220 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1221 *fPayloadCurr != fgkDataEndmarker)
1225 if (fCurrSm > -1 && fCurrSm < 18) {
1226 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fBytes += (fPayloadCurr - startPosLink) * sizeof(UInt_t);
1227 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fBytesRead += count * sizeof(UInt_t);
1228 fStats.fStatsSector[fCurrSm].fBytesRead += count * sizeof(UInt_t);
1229 fStats.fBytesRead += count * sizeof(UInt_t);
1235 Int_t AliTRDrawStream::ReadTracklets()
1237 // read the tracklets from one HC
1239 fTrackletArray->Clear();
1241 UInt_t *start = fPayloadCurr;
1242 while (*(fPayloadCurr) != fgkTrackletEndmarker &&
1243 fPayloadCurr - fPayloadStart < fPayloadSize) {
1244 new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(*(fPayloadCurr), fCurrHC);
1249 if (fTrackletArray->GetEntriesFast() > 0) {
1250 AliDebug(1, Form("Found %i tracklets in %i %i %i (ev. %i)", fTrackletArray->GetEntriesFast(),
1251 (fCurrEquipmentId-kDDLOffset), fCurrSlot, fCurrLink, fRawReader->GetEventIndex()));
1252 if (fCurrSm > -1 && fCurrSm < 18) {
1253 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNTracklets += fTrackletArray->GetEntriesFast();
1254 fStats.fStatsSector[fCurrSm].fNTracklets += fTrackletArray->GetEntriesFast();
1257 fTrackletTree->Fill();
1259 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
1260 new ((*fTracklets)[fTracklets->GetEntriesFast()]) AliTRDtrackletWord(*((AliTRDtrackletWord*)(*fTrackletArray)[iTracklet]));
1264 // loop over remaining tracklet endmarkers
1265 while ((*(fPayloadCurr) == fgkTrackletEndmarker &&
1266 fPayloadCurr - fPayloadStart < fPayloadSize))
1269 return fPayloadCurr - start;
1272 Int_t AliTRDrawStream::ReadHcHeader()
1274 // read and parse the HC header of one HC
1275 // and store the information in the corresponding variables
1277 AliDebug(1, Form("HC header: 0x%08x", *fPayloadCurr));
1278 UInt_t *start = fPayloadCurr;
1279 // check not to be at the data endmarker
1280 if (*fPayloadCurr == fgkDataEndmarker)
1283 fCurrSpecial = (*fPayloadCurr >> 31) & 0x1;
1284 fCurrMajor = (*fPayloadCurr >> 24) & 0x7f;
1285 fCurrMinor = (*fPayloadCurr >> 17) & 0x7f;
1286 fCurrAddHcWords = (*fPayloadCurr >> 14) & 0x7;
1287 fCurrSm = (*fPayloadCurr >> 9) & 0x1f;
1288 fCurrLayer = (*fPayloadCurr >> 6) & 0x7;
1289 fCurrStack = (*fPayloadCurr >> 3) & 0x7;
1290 fCurrSide = (*fPayloadCurr >> 2) & 0x1;
1291 fCurrCheck = (*fPayloadCurr) & 0x3;
1293 if ((fCurrSm != (((Int_t) fCurrEquipmentId) - kDDLOffset)) ||
1294 (fCurrStack != fCurrSlot) ||
1295 (fCurrLayer != fCurrLink / 2) ||
1296 (fCurrSide != fCurrLink % 2)) {
1297 LinkError(kHCmismatch,
1298 "HC: %i, %i, %i, %i\n 0x%08x 0x%08x 0x%08x 0x%08x",
1299 fCurrSm, fCurrStack, fCurrLayer, fCurrSide,
1300 fPayloadCurr[0], fPayloadCurr[1], fPayloadCurr[2], fPayloadCurr[3]);
1302 if (fCurrCheck != 0x1) {
1303 LinkError(kHCcheckFailed);
1306 if (fCurrAddHcWords > 0) {
1307 fCurrNtimebins = (fPayloadCurr[1] >> 26) & 0x3f;
1308 fCurrBC = (fPayloadCurr[1] >> 10) & 0xffff;
1309 fCurrPtrgCnt = (fPayloadCurr[1] >> 6) & 0xf;
1310 fCurrPtrgPhase = (fPayloadCurr[1] >> 2) & 0xf;
1313 fPayloadCurr += 1 + fCurrAddHcWords;
1315 return (fPayloadCurr - start);
1318 Int_t AliTRDrawStream::ReadTPData(Int_t mode)
1320 // testing of testpattern 1 to 3 (hardcoded), 0 missing
1321 // evcnt checking missing
1323 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 mcmcount = -1;
1327 Int_t wordcount = 0;
1328 Int_t channelcount = 0;
1330 UInt_t expadcval = 0;
1332 Int_t lastmcmpos = -1;
1333 Int_t lastrobpos = -1;
1335 UInt_t* start = fPayloadCurr;
1337 while (*(fPayloadCurr) != fgkDataEndmarker &&
1338 fPayloadCurr - fPayloadStart < fPayloadSize - 1) {
1340 // ----- Checking MCM Header -----
1341 AliDebug(2, DumpMcmHeader("MCM header: ", *fPayloadCurr));
1342 UInt_t *startPosMCM = fPayloadCurr;
1345 // ----- checking for proper readout order - ROB -----
1346 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1347 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1350 ROBError(kPosUnexp, Form("#%i after #%i in readout order", GetROBReadoutPos(ROB(*fPayloadCurr) / 2), lastrobpos));
1352 fCurrRobPos = ROB(*fPayloadCurr);
1354 // ----- checking for proper readout order - MCM -----
1355 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
1356 lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
1359 MCMError(kPosUnexp, Form("#%i after #%i in readout order", GetMCMReadoutPos(MCM(*fPayloadCurr)), lastmcmpos));
1361 fCurrMcmPos = MCM(*fPayloadCurr);
1366 evcnt = 0x3f & *fPayloadCurr >> 26;
1369 while (channelcount < 21) {
1371 if (cpu != cpufromchannel[channelcount]) {
1372 cpu = cpufromchannel[channelcount];
1373 expadcval = (1 << 9) | (fCurrRobPos << 6) | (fCurrMcmPos << 2) | cpu;
1377 while (count < 10) {
1378 if (*fPayloadCurr == fgkDataEndmarker) {
1379 MCMError(kMissTpData);
1380 return (fPayloadCurr - start);
1383 if (channelcount % 2 == 0)
1390 expword |= expadcval << 2;
1391 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
1392 expword |= expadcval << 12;
1393 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
1394 expword |= expadcval << 22;
1395 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
1397 else if (mode == 2) {
1398 // ----- TP 2 ------
1399 expword = ((0x3f & evcnt) << 26) | ((fCurrSm + 1) << 21) | ((fCurrLayer + 1) << 18) |
1400 ((fCurrStack + 1) << 15) |
1401 (fCurrRobPos << 12) | (fCurrMcmPos << 8) | (cpu << 6) | (wordcount + 1);
1403 else if (mode == 3) {
1405 expword = ((0xfff & evcnt) << 20) | (fCurrSm << 15) | (fCurrLink/2 << 12) | (fCurrStack << 9) |
1406 (fCurrRobPos << 6) | (fCurrMcmPos << 2) | (cpu << 0);
1410 LinkError(kTPmodeInvalid, "Just reading");
1413 diff = *fPayloadCurr ^ expword;
1414 AliDebug(11, Form("Comparing ch %2i, word %2i (cpu %i): 0x%08x <-> 0x%08x",
1415 channelcount, wordcount, cpu, *fPayloadCurr, expword));
1418 MCMError(kTPmismatch,
1419 "Seen 0x%08x, expected 0x%08x, diff: 0x%08x (0x%02x) - word %2i (cpu %i, ch %i)",
1420 *fPayloadCurr, expword, diff, 0xff & (diff | diff >> 8 | diff >> 16 | diff >> 24),
1421 wordcount, cpu, channelcount);;
1429 // continue with next MCM
1431 if (IsDumping() && DumpingMCM(fCurrHC/2, fCurrRobPos, fCurrMcmPos)) {
1432 AliInfo(DumpRaw(Form("Event %i: Det %3i ROB %i MCM %2i", fRawReader->GetEventIndex(), fCurrHC/2, fCurrRobPos, fCurrMcmPos),
1433 startPosMCM, fPayloadCurr - startPosMCM));
1437 return fPayloadCurr - start;
1441 Int_t AliTRDrawStream::ReadZSData()
1443 // read the zs data from one link from the current reading position
1445 UInt_t *start = fPayloadCurr;
1448 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
1449 Int_t channelcount = 0;
1450 Int_t channelcountExp = 0;
1451 Int_t channelcountMax = 0;
1453 Int_t currentTimebin = 0;
1456 Int_t lastmcmpos = -1;
1457 Int_t lastrobpos = -1;
1459 if (fCurrNtimebins != fNtimebins) {
1461 LinkError(kNtimebinsChanged,
1462 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
1463 fNtimebins = fCurrNtimebins;
1466 timebins = fNtimebins;
1468 while (*(fPayloadCurr) != fgkDataEndmarker &&
1469 fPayloadCurr - fPayloadStart < fPayloadSize) {
1471 // ----- Checking MCM Header -----
1472 AliDebug(2, DumpMcmHeader("MCM header: ", *fPayloadCurr));
1473 UInt_t *startPosMCM = fPayloadCurr;
1475 // ----- checking for proper readout order - ROB -----
1476 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1477 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) > lastrobpos)
1479 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1482 ROBError(kPosUnexp, Form("#%i after #%i and #%i in readout order",
1483 GetROBReadoutPos(ROB(*fPayloadCurr) / 2), lastrobpos, GetROBReadoutPos(fCurrRobPos)));
1485 fCurrRobPos = ROB(*fPayloadCurr);
1487 // ----- checking for proper readout order - MCM -----
1488 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) > lastmcmpos) {
1489 lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
1492 MCMError(kPosUnexp, Form("#%i after #%i and #%i in readout order",
1493 GetMCMReadoutPos(MCM(*fPayloadCurr)), lastmcmpos, GetMCMReadoutPos(fCurrMcmPos)));
1495 fCurrMcmPos = MCM(*fPayloadCurr);
1497 if (EvNo(*fPayloadCurr) != evno) {
1499 evno = EvNo(*fPayloadCurr);
1501 MCMError(kPtrgCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
1504 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
1505 Int_t padcoloff = PadColOffset(*fPayloadCurr);
1506 Int_t row = Row(*fPayloadCurr);
1509 // ----- Reading ADC channels -----
1510 AliDebug(2, DumpAdcMask("ADC mask: ", *fPayloadCurr));
1512 // ----- analysing the ADC mask -----
1514 channelcountExp = GetNActiveChannelsFromMask(*fPayloadCurr);
1515 channelcountMax = GetNActiveChannels(*fPayloadCurr);
1516 Int_t channelmask = GetActiveChannels(*fPayloadCurr);
1517 Int_t channelno = -1;
1520 if (channelcountExp != channelcountMax) {
1521 if (channelcountExp > channelcountMax) {
1522 Int_t temp = channelcountExp;
1523 channelcountExp = channelcountMax;
1524 channelcountMax = temp;
1526 while (channelcountExp < channelcountMax && channelcountExp < 21 &&
1527 fPayloadCurr - fPayloadStart < fPayloadSize - 10 * channelcountExp - 1) {
1528 MCMError(kAdcMaskInconsistent,
1529 "Possible MCM-H: 0x%08x, possible ADC-mask: 0x%08x",
1530 *(fPayloadCurr + 10 * channelcountExp),
1531 *(fPayloadCurr + 10 * channelcountExp + 1) );
1532 if (!CouldBeMCMhdr( *(fPayloadCurr + 10 * channelcountExp)) && !CouldBeADCmask( *(fPayloadCurr + 10 * channelcountExp + 1)))
1538 MCMError(kAdcMaskInconsistent,
1539 "Inconsistency in no. of active channels: Counter: %i, Mask: %i, chosen: %i!",
1540 GetNActiveChannels(fPayloadCurr[-1]), GetNActiveChannelsFromMask(fPayloadCurr[-1]), channelcountExp);
1542 AliDebug(2, Form("expecting %i active channels, %i timebins", channelcountExp, fCurrNtimebins));
1544 // ----- reading marked ADC channels -----
1545 while (channelcount < channelcountExp && *(fPayloadCurr) != fgkDataEndmarker) {
1548 while (channelno < 20 && (channelmask & 1 << channelno) == 0)
1551 if (fCurrNtimebins > 30) {
1552 currentTimebin = ((*fPayloadCurr >> 2) & 0x3f);
1553 timebins = ((*fPayloadCurr >> 8) & 0xf) * 3;
1560 Int_t nADCwords = (timebins + 2) / 3;
1561 AliDebug(3, Form("Now reading %i words for channel %2i", nADCwords, channelno));
1562 Int_t adccol = adccoloff - channelno;
1563 Int_t padcol = padcoloff - channelno;
1564 // if (adccol < 3 || adccol > 165)
1565 // AliInfo(Form("writing channel %i of det %3i %i:%2i to adcrow/-col: %i/%i padcol: %i",
1566 // channelno, fCurrHC/2, fCurrRobPos, fCurrMcmPos, row, adccol, padcol));
1568 while ((adcwc < nADCwords) &&
1569 (*(fPayloadCurr) != fgkDataEndmarker) &&
1570 (fPayloadCurr - fPayloadStart < fPayloadSize)) {
1571 int check = 0x3 & *fPayloadCurr;
1572 if (channelno % 2 != 0) { // odd channel
1573 if (check != 0x2 && channelno < 21) {
1574 MCMError(kAdcCheckInvalid,
1575 "%i for %2i. ADC word in odd channel %i",
1576 check, adcwc+1, channelno);
1579 else { // even channel
1580 if (check != 0x3 && channelno < 21) {
1581 MCMError(kAdcCheckInvalid,
1582 "%i for %2i. ADC word in even channel %i",
1583 check, adcwc+1, channelno);
1587 // filling the actual timebin data
1588 int tb2 = 0x3ff & *fPayloadCurr >> 22;
1589 int tb1 = 0x3ff & *fPayloadCurr >> 12;
1590 int tb0 = 0x3ff & *fPayloadCurr >> 2;
1591 if (adcwc != 0 || fCurrNtimebins <= 30)
1592 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
1595 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
1596 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
1602 if (adcwc != nADCwords)
1603 MCMError(kAdcDataAbort);
1606 if (padcol > 0 && padcol < 144) {
1607 fSignalIndex->AddIndexRC(row, padcol);
1613 if (fCurrSm > -1 && fCurrSm < 18) {
1614 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNChannels += channelcount;
1615 fStats.fStatsSector[fCurrSm].fNChannels += channelcount;
1617 if (channelcount != channelcountExp)
1618 MCMError(kAdcChannelsMiss);
1621 if (fCurrSm > -1 && fCurrSm < 18) {
1622 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNMCMs++;
1623 fStats.fStatsSector[fCurrSm].fNMCMs++;
1626 if (IsDumping() && DumpingMCM(fCurrHC/2, fCurrRobPos, fCurrMcmPos)) {
1627 AliInfo(DumpRaw(Form("Event %i: Det %3i ROB %i MCM %2i", fRawReader->GetEventIndex(), fCurrHC/2, fCurrRobPos, fCurrMcmPos),
1628 startPosMCM, fPayloadCurr - startPosMCM));
1631 // continue with next MCM
1634 // check for missing MCMs (if header suppression is inactive)
1635 if (((fCurrMajor & 0x1) == 0) && (mcmcount != mcmcountExp)) {
1636 LinkError(kMissMcmHeaders,
1637 "No. of MCM headers %i not as expected: %i",
1638 mcmcount, mcmcountExp);
1641 return (fPayloadCurr - start);
1644 Int_t AliTRDrawStream::ReadNonZSData()
1646 // read the non-zs data from one link from the current reading position
1648 UInt_t *start = fPayloadCurr;
1651 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
1652 Int_t channelcount = 0;
1653 Int_t channelcountExp = 0;
1655 Int_t currentTimebin = 0;
1658 Int_t lastmcmpos = -1;
1659 Int_t lastrobpos = -1;
1661 if (fCurrNtimebins != fNtimebins) {
1663 LinkError(kNtimebinsChanged,
1664 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
1665 fNtimebins = fCurrNtimebins;
1668 timebins = fNtimebins;
1670 while (*(fPayloadCurr) != fgkDataEndmarker &&
1671 fPayloadCurr - fPayloadStart < fPayloadSize - 2) {
1673 // ----- Checking MCM Header -----
1674 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
1676 // ----- checking for proper readout order - ROB -----
1677 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1678 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1681 ROBError(kPosUnexp, Form("#%i after #%i in readout order", GetROBReadoutPos(ROB(*fPayloadCurr) / 2), lastrobpos));
1683 fCurrRobPos = ROB(*fPayloadCurr);
1685 // ----- checking for proper readout order - MCM -----
1686 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
1687 lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
1690 MCMError(kPosUnexp, Form("#%i after #%i in readout order", GetMCMReadoutPos(MCM(*fPayloadCurr)), lastmcmpos));
1692 fCurrMcmPos = MCM(*fPayloadCurr);
1694 if (EvNo(*fPayloadCurr) != evno) {
1696 evno = EvNo(*fPayloadCurr);
1698 MCMError(kPtrgCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
1703 channelcountExp = 21;
1706 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
1707 Int_t padcoloff = PadColOffset(*fPayloadCurr);
1708 Int_t row = Row(*fPayloadCurr);
1712 // ----- reading marked ADC channels -----
1713 while (channelcount < channelcountExp &&
1714 *(fPayloadCurr) != fgkDataEndmarker) {
1721 Int_t nADCwords = (timebins + 2) / 3;
1722 AliDebug(2, Form("Now looking %i words", nADCwords));
1723 Int_t adccol = adccoloff - channelno;
1724 Int_t padcol = padcoloff - channelno;
1725 while ((adcwc < nADCwords) &&
1726 (*(fPayloadCurr) != fgkDataEndmarker) &&
1727 (fPayloadCurr - fPayloadStart < fPayloadSize)) {
1728 int check = 0x3 & *fPayloadCurr;
1729 if (channelno % 2 != 0) { // odd channel
1730 if (check != 0x2 && channelno < 21) {
1731 MCMError(kAdcCheckInvalid,
1732 "%i for %2i. ADC word in odd channel %i",
1733 check, adcwc+1, channelno);
1736 else { // even channel
1737 if (check != 0x3 && channelno < 21) {
1738 MCMError(kAdcCheckInvalid,
1739 "%i for %2i. ADC word in even channel %i",
1740 check, adcwc+1, channelno);
1744 // filling the actual timebin data
1745 int tb2 = 0x3ff & *fPayloadCurr >> 22;
1746 int tb1 = 0x3ff & *fPayloadCurr >> 12;
1747 int tb0 = 0x3ff & *fPayloadCurr >> 2;
1748 if (adcwc != 0 || fCurrNtimebins <= 30)
1749 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
1752 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
1753 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
1759 if (adcwc != nADCwords)
1760 MCMError(kAdcDataAbort);
1763 if (padcol > 0 && padcol < 144) {
1764 fSignalIndex->AddIndexRC(row, padcol);
1770 if (channelcount != channelcountExp)
1771 MCMError(kAdcChannelsMiss);
1773 // continue with next MCM
1776 // check for missing MCMs (if header suppression is inactive)
1777 if (mcmcount != mcmcountExp) {
1778 LinkError(kMissMcmHeaders,
1779 "%i not as expected: %i", mcmcount, mcmcountExp);
1782 return (fPayloadCurr - start);
1785 Int_t AliTRDrawStream::SeekNextLink()
1787 // proceed in raw data stream till the next link
1789 UInt_t *start = fPayloadCurr;
1791 // read until data endmarkers
1792 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1793 *fPayloadCurr != fgkDataEndmarker)
1796 // read all data endmarkers
1797 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1798 *fPayloadCurr == fgkDataEndmarker)
1801 return (fPayloadCurr - start);
1804 Bool_t AliTRDrawStream::ConnectTracklets(TTree *trklTree)
1806 // connect the tracklet tree used to store the tracklet output
1808 fTrackletTree = trklTree;
1812 if (!fTrackletTree->GetBranch("hc"))
1813 fTrackletTree->Branch("hc", &fCurrHC, "hc/I");
1815 fTrackletTree->SetBranchAddress("hc", &fCurrHC);
1817 if (!fTrackletTree->GetBranch("trkl"))
1818 fTrackletTree->Branch("trkl", &fTrackletArray);
1820 fTrackletTree->SetBranchAddress("trkl", &fTrackletArray);
1826 void AliTRDrawStream::EquipmentError(ErrorCode_t err, const char *const msg, ...)
1828 // register error according to error code on equipment level
1829 // and return the corresponding error message
1831 fLastError.fSector = fCurrEquipmentId - kDDLOffset;
1832 fLastError.fStack = -1;
1833 fLastError.fLink = -1;
1834 fLastError.fRob = -1;
1835 fLastError.fMcm = -1;
1836 fLastError.fError = err;
1837 (this->*fStoreError)();
1840 if (fgErrorDebugLevel[err] > 10)
1841 AliDebug(fgErrorDebugLevel[err],
1842 Form("Event %6i: Eq. %2d - %s : %s",
1843 fRawReader->GetEventIndex(), fCurrEquipmentId, fgkErrorMessages[err],
1844 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1846 AliError(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 fErrorFlags |= fgErrorBehav[err];
1853 void AliTRDrawStream::StackError(ErrorCode_t err, const char *const msg, ...)
1855 // register error according to error code on stack level
1856 // and return the corresponding error message
1858 fLastError.fSector = fCurrEquipmentId - kDDLOffset;
1859 fLastError.fStack = fCurrSlot;
1860 fLastError.fLink = -1;
1861 fLastError.fRob = -1;
1862 fLastError.fMcm = -1;
1863 fLastError.fError = err;
1864 (this->*fStoreError)();
1867 if (fgErrorDebugLevel[err] > 0)
1868 AliDebug(fgErrorDebugLevel[err],
1869 Form("Event %6i: Eq. %2d S %i - %s : %s",
1870 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgkErrorMessages[err],
1871 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1873 AliError(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 fErrorFlags |= fgErrorBehav[err];
1880 void AliTRDrawStream::LinkError(ErrorCode_t err, const char *const msg, ...)
1882 // register error according to error code on link level
1883 // and return the corresponding error message
1885 fLastError.fSector = fCurrEquipmentId - kDDLOffset;
1886 fLastError.fStack = fCurrSlot;
1887 fLastError.fLink = fCurrLink;
1888 fLastError.fRob = -1;
1889 fLastError.fMcm = -1;
1890 fLastError.fError = err;
1891 (this->*fStoreError)();
1894 if (fgErrorDebugLevel[err] > 0)
1895 AliDebug(fgErrorDebugLevel[err],
1896 Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
1897 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgkErrorMessages[err],
1898 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1900 AliError(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 fErrorFlags |= fgErrorBehav[err];
1907 void AliTRDrawStream::ROBError(ErrorCode_t err, const char *const msg, ...)
1909 // register error according to error code on ROB level
1910 // and return the corresponding error message
1912 fLastError.fSector = fCurrEquipmentId - kDDLOffset;
1913 fLastError.fStack = fCurrSlot;
1914 fLastError.fLink = fCurrLink;
1915 fLastError.fRob = fCurrRobPos;
1916 fLastError.fMcm = -1;
1917 fLastError.fError = err;
1918 (this->*fStoreError)();
1921 if (fgErrorDebugLevel[err] > 0)
1922 AliDebug(fgErrorDebugLevel[err],
1923 Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : %s",
1924 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgkErrorMessages[err],
1925 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1927 AliError(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 fErrorFlags |= fgErrorBehav[err];
1934 void AliTRDrawStream::MCMError(ErrorCode_t err, const char *const msg, ...)
1936 // register error according to error code on MCM level
1937 // and return the corresponding error message
1939 fLastError.fSector = fCurrEquipmentId - kDDLOffset;
1940 fLastError.fStack = fCurrSlot;
1941 fLastError.fLink = fCurrLink;
1942 fLastError.fRob = fCurrRobPos;
1943 fLastError.fMcm = fCurrMcmPos;
1944 fLastError.fError = err;
1945 (this->*fStoreError)();
1948 if (fgErrorDebugLevel[err] > 0)
1949 AliDebug(fgErrorDebugLevel[err],
1950 Form("Event %6i: Eq. %2d S %i l %2i ROB %i MCM %2i - %s : %s",
1951 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fCurrMcmPos, fgkErrorMessages[err],
1952 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1954 AliError(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 fErrorFlags |= fgErrorBehav[err];
1960 const char* AliTRDrawStream::GetErrorMessage(ErrorCode_t errCode)
1962 // return the error message for the given error code
1964 if (errCode > 0 && errCode < kLastErrorCode)
1965 return fgkErrorMessages[errCode];
1970 void AliTRDrawStream::AliTRDrawStats::ClearStats()
1972 // clear statistics (includes clearing sector-wise statistics)
1975 for (Int_t iSector = 0; iSector < 18; iSector++) {
1976 fStatsSector[iSector].ClearStats();
1981 void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::ClearStats()
1983 // clear statistics (includes clearing HC-wise statistics)
1991 for (Int_t iHC = 0; iHC < 60; iHC++) {
1992 fStatsHC[iHC].ClearStats();
1996 void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::AliTRDrawStatsHC::ClearStats()
2007 void AliTRDrawStream::SetDumpMCM(Int_t det, Int_t rob, Int_t mcm, Bool_t dump)
2009 // mark MCM for dumping of raw data
2012 fDumpMCM[fNDumpMCMs++] = (det << 7) | (rob << 4) | mcm;
2016 for (iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
2017 if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
2022 for ( ; iMCM < fNDumpMCMs; iMCM++) {
2023 fDumpMCM[iMCM] = fDumpMCM[iMCM+1];
2028 Bool_t AliTRDrawStream::DumpingMCM(Int_t det, Int_t rob, Int_t mcm) const
2030 // check if MCM data should be dumped
2032 for (Int_t iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
2033 if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
2040 TString AliTRDrawStream::DumpRaw(TString title, UInt_t *start, Int_t length, UInt_t endmarker)
2045 for (Int_t pos = 0; pos < length; pos += 4) {
2046 if ((start[pos+0] != endmarker) && pos+0 < length)
2047 if ((start[pos+1] != endmarker && pos+1 < length))
2048 if ((start[pos+2] != endmarker && pos+2 < length))
2049 if ((start[pos+3] != endmarker && pos+3 < length))
2050 title += Form(" 0x%08x 0x%08x 0x%08x 0x%08x\n",
2051 start[pos+0], start[pos+1], start[pos+2], start[pos+3]);
2053 title += Form(" 0x%08x 0x%08x 0x%08x 0x%08x\n",
2054 start[pos+0], start[pos+1], start[pos+2], start[pos+3]);
2058 title += Form(" 0x%08x 0x%08x 0x%08x\n",
2059 start[pos+0], start[pos+1], start[pos+2]);
2063 title += Form(" 0x%08x 0x%08x\n",
2064 start[pos+0], start[pos+1]);
2068 title += Form(" 0x%08x\n",
2076 TString AliTRDrawStream::DumpMcmHeader(TString title, UInt_t word)
2078 title += Form("0x%08x -> ROB: %i, MCM: %2i",
2079 word, ROB(word), MCM(word));
2083 TString AliTRDrawStream::DumpAdcMask(TString title, UInt_t word)
2085 title += Form("0x%08x -> #ch : %2i, 0x%06x (%2i ch)",
2086 word, GetNActiveChannels(word), GetActiveChannels(word), GetNActiveChannelsFromMask(word));
2090 AliTRDrawStream::AliTRDrawStreamError::AliTRDrawStreamError(Int_t error, Int_t sector, Int_t stack, Int_t link, Int_t rob, Int_t mcm) :