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 //
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"
43 #include "AliTRDrawStream.h"
46 #include "AliRunLoader.h"
48 ClassImp(AliTRDrawStream)
50 // some static information
51 const Int_t AliTRDrawStream::fgkMcmOrder[] = {12, 13, 14, 15,
55 const Int_t AliTRDrawStream::fgkRobOrder [] = {0, 1, 2, 3};
56 const Int_t AliTRDrawStream::fgkNlinks = 12;
57 const Int_t AliTRDrawStream::fgkNstacks = 5;
58 const UInt_t AliTRDrawStream::fgkDataEndmarker = 0x00000000;
59 const UInt_t AliTRDrawStream::fgkTrackletEndmarker = 0x10001000;
61 const char* AliTRDrawStream::fgErrorMessages[] = {
63 "Link monitor active",
64 "Pretrigger counter mismatch",
65 "not a TRD equipment (1024-1041)",
66 "Invalid Stack header",
67 "Invalid detector number",
68 "No digits could be retrieved from the digitsmanager",
70 "HC check bits wrong",
71 "Unexpected position in readout stream",
72 "Invalid testpattern mode",
73 "Testpattern mismatch",
74 "Number of timebins changed",
75 "ADC mask inconsistent",
76 "ADC check bits invalid",
78 "Missing expected ADC channels",
82 const Int_t AliTRDrawStream::fgErrorDebugLevel[] = {
103 AliTRDrawStream::ErrorBehav_t AliTRDrawStream::fgErrorBehav[] = {
104 AliTRDrawStream::kTolerate,
105 AliTRDrawStream::kDiscardHC,
106 AliTRDrawStream::kTolerate,
107 AliTRDrawStream::kAbort,
108 AliTRDrawStream::kAbort,
109 AliTRDrawStream::kAbort,
110 AliTRDrawStream::kAbort,
111 AliTRDrawStream::kDiscardHC,
112 AliTRDrawStream::kDiscardHC,
113 AliTRDrawStream::kTolerate,
114 AliTRDrawStream::kTolerate,
115 AliTRDrawStream::kTolerate,
116 AliTRDrawStream::kTolerate,
117 AliTRDrawStream::kTolerate,
118 AliTRDrawStream::kTolerate,
119 AliTRDrawStream::kTolerate,
120 AliTRDrawStream::kTolerate,
121 AliTRDrawStream::kTolerate
124 AliTRDrawStream::AliTRDrawStream(AliRawReader *rawReader) :
126 fStoreError(&AliTRDrawStream::StoreErrorTree),
127 fRawReader(rawReader),
143 fCurrSmuIndexHeaderSize(0),
144 fCurrSmuIndexHeaderVersion(0),
146 fCurrTrackletEnable(0),
148 fCurrStackIndexWord(0x0),
149 fCurrStackHeaderSize(0x0),
150 fCurrStackHeaderVersion(0x0),
152 fCurrCleanCheckout(0x0),
155 fCurrLinkMonitorFlags(0x0),
156 fCurrLinkDataTypeFlags(0x0),
157 fCurrLinkDebugFlags(0x0),
181 // default constructor
183 fCurrStackIndexWord = new UInt_t[fgkNstacks];
184 fCurrStackHeaderSize = new UInt_t[fgkNstacks];
185 fCurrStackHeaderVersion = new UInt_t[fgkNstacks];
186 fCurrLinkMask = new UInt_t[fgkNstacks];
187 fCurrCleanCheckout = new UInt_t[fgkNstacks];
188 fCurrBoardId = new UInt_t[fgkNstacks];
189 fCurrHwRev = new UInt_t[fgkNstacks];
190 fCurrLinkMonitorFlags = new UInt_t[fgkNstacks * fgkNlinks];
191 fCurrLinkDataTypeFlags = new UInt_t[fgkNstacks * fgkNlinks];
192 fCurrLinkDebugFlags = new UInt_t[fgkNstacks * fgkNlinks];
193 for (Int_t i = 0; i < 100; i++)
196 // preparing TClonesArray
197 fTrackletArray = new TClonesArray("AliTRDtrackletWord", 256);
199 // setting up the error tree
200 fErrors = new TTree("errorStats", "Error statistics");
201 fErrors->SetDirectory(0x0);
202 fErrors->Branch("error", &fLastError);
203 fErrors->SetCircular(1000);
204 for (Int_t i = 0; i < 100; i++) {
210 AliTRDrawStream::~AliTRDrawStream()
216 delete [] fCurrStackIndexWord;
217 delete [] fCurrStackHeaderSize;
218 delete [] fCurrStackHeaderVersion;
219 delete [] fCurrLinkMask;
220 delete [] fCurrCleanCheckout;
221 delete [] fCurrBoardId;
222 delete [] fCurrHwRev;
223 delete [] fCurrLinkMonitorFlags;
224 delete [] fCurrLinkDataTypeFlags;
225 delete [] fCurrLinkDebugFlags;
228 Bool_t AliTRDrawStream::ReadEvent(TTree *trackletTree)
230 // read the current event from the raw reader and fill it to the digits manager
233 AliError("No raw reader available");
238 ConnectTracklets(trackletTree);
243 // loop over all DDLs
244 // data starts with GTU payload, i.e. SMU index word
245 UChar_t *buffer = 0x0;
247 while (fRawReader->ReadNextData(buffer)) {
249 fCurrEquipmentId = fRawReader->GetEquipmentId();
250 AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
252 if (fCurrEquipmentId < 1024 || fCurrEquipmentId > 1041) {
253 EquipmentError(kNonTrdEq, "Skipping");
257 // setting the pointer to data and current reading position
258 fPayloadCurr = fPayloadStart = (UInt_t*) (buffer);
259 fPayloadSize = fRawReader->GetDataSize() / sizeof(UInt_t);
260 fStats.fStatsSector[fCurrEquipmentId - 1024].fBytes = fRawReader->GetDataSize();
261 AliDebug(2, Form("Read buffer of size: %i", fRawReader->GetDataSize()));
263 // read SMU index header
264 if (ReadSmHeader() < 0) {
265 AliError(Form("Reading SMU header failed, skipping this DDL %i", fCurrEquipmentId));
269 // read stack index header
270 for (Int_t iStack = 0; iStack < 5; iStack++) {
271 if ((fCurrStackMask & (1 << iStack)) != 0)
272 ReadStackIndexHeader(iStack);
275 for (Int_t iStack = 0; iStack < 5; iStack++) {
277 if ((fCurrStackMask & (1 << fCurrSlot)) == 0)
280 AliDebug(2, Form("Stack %i, Link mask: 0x%02x", fCurrSlot, fCurrLinkMask[fCurrSlot]));
281 for (Int_t iLink = 0; iLink < 12; iLink++) {
283 fCurrHC = fCurrSm * 60 + fCurrSlot * 12 + iLink;
284 if ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) == 0)
288 // check for link monitor error flag
289 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0)
290 LinkError(kLinkMonitor);
292 // read the data from one HC
295 // read all data endmarkers
304 Bool_t AliTRDrawStream::NextDDL()
306 // continue reading with the next equipment
311 fCurrEquipmentId = 0;
315 UChar_t *buffer = 0x0;
317 while (fRawReader->ReadNextData(buffer)) {
319 fCurrEquipmentId = fRawReader->GetEquipmentId();
320 AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
322 if (fCurrEquipmentId < 1024 || fCurrEquipmentId > 1041) {
323 EquipmentError(kNonTrdEq, "Skipping");
327 // setting the pointer to data and current reading position
328 fPayloadCurr = fPayloadStart = (UInt_t*) (buffer);
329 fPayloadSize = fRawReader->GetDataSize() / sizeof(UInt_t);
330 AliDebug(2, Form("Read buffer of size: %i", fRawReader->GetDataSize()));
332 // read SMU index header
333 if (ReadSmHeader() < 0) {
334 AliError(Form("Reading SMU header failed, skipping this DDL %i", fCurrEquipmentId));
338 // read stack index header
339 for (Int_t iStack = 0; iStack < 5; iStack++) {
340 if ((fCurrStackMask & (1 << iStack)) != 0) {
341 ReadStackIndexHeader(iStack);
351 Int_t AliTRDrawStream::NextChamber(AliTRDdigitsManager *digMgr, UInt_t ** /* trackletContainer */, UShort_t ** /* errorContainer */)
353 // read the data for the next chamber
354 // in case you only want to read the data of a single chamber
355 // to read all data ReadEvent(...) is recommended
357 fDigitsManager = digMgr;
362 // tracklet output preparation
363 TTree *trklTree = 0x0;
364 AliRunLoader *rl = AliRunLoader::Instance();
365 AliLoader* trdLoader = rl ? rl->GetLoader("TRDLoader") : NULL;
366 AliDataLoader *trklLoader = trdLoader ? trdLoader->GetDataLoader("tracklets") : NULL;
368 AliTreeLoader *trklTreeLoader = (AliTreeLoader*) trklLoader->GetBaseLoader("tracklets-raw");
370 trklTree = trklTreeLoader->Tree();
372 trklTree = trklLoader->Tree();
375 if (fTrackletTree != trklTree)
376 ConnectTracklets(trklTree);
379 AliError("No raw reader available");
383 if (fCurrSlot < 0 || fCurrSlot >= 5) {
390 AliDebug(2, Form("Stack %i, Link %i, mask: 0x%02x", fCurrSlot, fCurrLink, fCurrLinkMask[fCurrSlot]));
391 fCurrHC = (fCurrEquipmentId - 1024) * 60 + fCurrSlot * 12 + fCurrLink;
393 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0)
394 LinkError(kLinkMonitor);
396 // read the data from one HC
399 // read all data endmarkers
402 if (fCurrLink % 2 == 0) {
403 // if we just read the A-side HC then also check the B-side
406 if (fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) {
415 if (fCurrLink > 11) {
419 } while ((fCurrSlot < 5) &&
420 (((fCurrStackMask & (1 << fCurrSlot)) == 0) ||
421 ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink))) == 0));
423 // return chamber information from HC if it is valid
424 // otherwise return information from link position
425 if (fCurrSm < 0 || fCurrSm > 17 || fCurrStack < 0 || fCurrStack > 4 || fCurrLayer < 0 || fCurrLayer > 5)
426 return ((fCurrEquipmentId-1024) + fCurrSlot * 6 + fCurrLink/2);
428 return (fCurrSm * 30 + fCurrStack * 6 + fCurrLayer);
432 Int_t AliTRDrawStream::ReadSmHeader()
434 // read the SMU index header at the current reading position
435 // and store the information in the corresponding variables
437 if (fPayloadCurr - fPayloadStart >= fPayloadSize - 1) {
438 EquipmentError(kUnknown, "SM Header incomplete");
442 fCurrSmuIndexHeaderSize = ((*fPayloadCurr) >> 16) & 0xffff;
443 fCurrSmuIndexHeaderVersion = ((*fPayloadCurr) >> 12) & 0xf;
444 // fCurrSmuIndexHeaderTrgAvail = ((*fPayloadCurr) >> 9) & 0x1;
445 // fCurrSmuIndexHeaderEvType = ((*fPayloadCurr) >> 7) & 0x3;
446 fCurrTrackEnable = ((*fPayloadCurr) >> 6) & 0x1;
447 fCurrTrackletEnable = ((*fPayloadCurr) >> 5) & 0x1;
448 fCurrStackMask = ((*fPayloadCurr) ) & 0x1f;
450 AliDebug(5, Form("SMU header: size: %i, version: %i, track enable: %i, tracklet enable: %i, stack mask: %2x",
451 fCurrSmuIndexHeaderSize,
452 fCurrSmuIndexHeaderVersion,
457 // decode GTU track words
458 UInt_t trackWord[2] = { 0, 0 };
461 for (UInt_t iWord = 4; iWord < fCurrSmuIndexHeaderSize; iWord++) {
462 if (fPayloadCurr[iWord] == 0x10000000) {
468 ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
469 AliDebug(1,Form("stack %i: fast trigger word: 0x%08x", stack, fPayloadCurr[iWord]));
472 else if ((idx & 0x1)==0x1) {
473 trackWord[idx&0x1] = fPayloadCurr[iWord];
474 AliDebug(1,Form("track debug word: 0x%08x%08x", trackWord[1], trackWord[0]));
476 // new ((*fTracks)[fTracks->GetEntriesFast()]) AliESDTrdTrack(0, 0, trackWord[0], trackWord[1], fCurrEquipmentId-1024);
479 trackWord[idx&0x1] = fPayloadCurr[iWord];
485 fPayloadCurr += fCurrSmuIndexHeaderSize + 1;
487 return fCurrSmuIndexHeaderSize + 1;
490 Int_t AliTRDrawStream::ReadStackIndexHeader(Int_t stack)
492 // read the stack index header
493 // and store the information in the corresponding variables
495 fCurrStackIndexWord[stack] = *fPayloadCurr;
496 fCurrStackHeaderSize[stack] = (((*fPayloadCurr) >> 16) & 0xffff) + 1;
497 fCurrStackHeaderVersion[stack] = ((*fPayloadCurr) >> 12) & 0xf;
498 fCurrLinkMask[stack] = (*fPayloadCurr) & 0xfff;
500 if (fPayloadCurr - fPayloadStart >= fPayloadSize - (Int_t) fCurrStackHeaderSize[stack]) {
501 StackError(kStackHeaderInvalid, "Stack index header aborted");
505 switch (fCurrStackHeaderVersion[stack]) {
507 if (fCurrStackHeaderSize[stack] < 8) {
508 StackError(kStackHeaderInvalid, "Stack header smaller than expected!");
512 fCurrCleanCheckout[stack] = fPayloadCurr[1] & 0x1;
513 fCurrBoardId[stack] = (fPayloadCurr[1] >> 8) & 0xff;
514 fCurrHwRev[stack] = (fPayloadCurr[1] >> 16) & 0xffff;
516 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
518 fCurrLinkMonitorFlags [stack*fgkNlinks + iLayer*2] = fPayloadCurr[iLayer+2] & 0xf;
519 fCurrLinkDataTypeFlags [stack*fgkNlinks + iLayer*2] = (fPayloadCurr[iLayer+2] >> 4) & 0x3;
520 fCurrLinkDebugFlags [stack*fgkNlinks + iLayer*2] = (fPayloadCurr[iLayer+2] >> 12) & 0xf;
522 fCurrLinkMonitorFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 16) & 0xf;
523 fCurrLinkDataTypeFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 20) & 0x3;
524 fCurrLinkDebugFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 28) & 0xf;
529 StackError(kStackHeaderInvalid, "Invalid Stack Index Header version %x", fCurrStackHeaderVersion[stack]);
532 fPayloadCurr += fCurrStackHeaderSize[stack];
534 return fCurrStackHeaderSize[stack];
537 Int_t AliTRDrawStream::ReadLinkData()
539 // read the data in one link (one HC) until the data endmarker is reached
540 // returns the number of words read!
543 UInt_t* startPosLink = fPayloadCurr;
545 // printf("----- HC: %i -----\n", fCurrHC);
546 // for (Int_t i = 0; i < 3; i++) {
547 // printf("0x%08x 0x%08x 0x%08x 0x%08x\n",
548 // fPayloadCurr[i*4+0], fPayloadCurr[i*4+1], fPayloadCurr[i*4+2], fPayloadCurr[i*4+3]);
552 new ((*fMarkers)[fMarkers->GetEntriesFast()])
553 AliTRDrawStreamError(-kHCactive, fCurrSm, fCurrStack, fCurrLink);
555 if (fErrorFlags & kDiscardHC)
558 count += ReadTracklets();
559 if (fErrorFlags & kDiscardHC)
562 count += ReadHcHeader();
563 if (fErrorFlags & kDiscardHC)
566 Int_t det = fCurrSm * 30 + fCurrStack * 6 + fCurrLayer;
568 if (det > -1 && det < 540) {
570 if ((fAdcArray = fDigitsManager->GetDigits(det))) {
571 //fAdcArray->Expand();
572 if (fAdcArray->GetNtime() != fCurrNtimebins)
573 fAdcArray->Allocate(16, 144, fCurrNtimebins);
576 LinkError(kNoDigits);
580 fDigitsParam = fDigitsManager->GetDigitsParam();
583 fDigitsParam->SetPretriggerPhase(det, fCurrPtrgPhase);
584 fDigitsParam->SetNTimeBins(det, fCurrNtimebins);
585 fDigitsParam->SetADCbaseline(det, 10);
588 if (fDigitsManager->UsesDictionaries()) {
589 fDigitsManager->GetDictionary(det, 0)->Reset();
590 fDigitsManager->GetDictionary(det, 1)->Reset();
591 fDigitsManager->GetDictionary(det, 2)->Reset();
594 if ((fSignalIndex = fDigitsManager->GetIndexes(det))) {
595 fSignalIndex->SetSM(fCurrSm);
596 fSignalIndex->SetStack(fCurrStack);
597 fSignalIndex->SetLayer(fCurrLayer);
598 fSignalIndex->SetDetNumber(det);
599 if (!fSignalIndex->IsAllocated())
600 fSignalIndex->Allocate(16, 144, fCurrNtimebins);
603 // ----- check which kind of data -----
604 if (fCurrMajor & 0x40) {
605 if ((fCurrMajor & 0x7) == 0x7) {
606 AliDebug(1, "This is a config event");
607 UInt_t *startPos = fPayloadCurr;
608 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
609 *fPayloadCurr != fgkDataEndmarker)
611 count += fPayloadCurr - startPos;
613 // feeding TRAP config
614 AliTRDtrapConfig *trapcfg = AliTRDtrapConfig::Instance();
615 trapcfg->ReadPackedConfig(fCurrHC, startPos, fPayloadCurr - startPos);
618 Int_t tpmode = fCurrMajor & 0x7;
619 AliDebug(1, Form("Checking testpattern (mode %i) data", tpmode));
623 else if (fCurrMajor & 0x20) {
624 AliDebug(1, "This is a zs event");
625 count += ReadZSData();
628 AliDebug(1, "This is a nozs event");
629 count += ReadNonZSData();
633 LinkError(kInvalidDetector, "%i", det);
634 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
635 *fPayloadCurr != fgkDataEndmarker)
639 if (fCurrSm > -1 && fCurrSm < 18) {
640 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fBytes += (fPayloadCurr - startPosLink) * sizeof(UInt_t);
641 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fBytesRead += count * sizeof(UInt_t);
642 fStats.fStatsSector[fCurrSm].fBytesRead += count * sizeof(UInt_t);
643 fStats.fBytesRead += count * sizeof(UInt_t);
649 Int_t AliTRDrawStream::ReadTracklets()
651 // read the tracklets from one HC
653 fTrackletArray->Clear();
655 UInt_t *start = fPayloadCurr;
656 while (*(fPayloadCurr) != fgkTrackletEndmarker &&
657 fPayloadCurr - fPayloadStart < fPayloadSize) {
659 new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(*(fPayloadCurr), fCurrHC);
664 if (fTrackletArray->GetEntriesFast() > 0) {
665 AliDebug(1, Form("Found %i tracklets in %i %i %i (ev. %i)", fTrackletArray->GetEntriesFast(),
666 fCurrSm, fCurrSlot, fCurrLink, fRawReader->GetEventIndex()));
667 if (fCurrSm > -1 && fCurrSm < 18) {
668 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNTracklets += fTrackletArray->GetEntriesFast();
669 fStats.fStatsSector[fCurrSm].fNTracklets += fTrackletArray->GetEntriesFast();
672 fTrackletTree->Fill();
674 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
675 new ((*fTracklets)[fTracklets->GetEntriesFast()]) AliTRDtrackletWord(*((AliTRDtrackletWord*)(*fTrackletArray)[iTracklet]));
679 // loop over remaining tracklet endmarkers
680 while ((*(fPayloadCurr) == fgkTrackletEndmarker &&
681 fPayloadCurr - fPayloadStart < fPayloadSize))
684 return fPayloadCurr - start;
687 Int_t AliTRDrawStream::ReadHcHeader()
689 // read and parse the HC header of one HC
690 // and store the information in the corresponding variables
692 UInt_t *start = fPayloadCurr;
693 // check not to be at the data endmarker
694 if (*fPayloadCurr == fgkDataEndmarker)
697 fCurrSpecial = (*fPayloadCurr >> 31) & 0x1;
698 fCurrMajor = (*fPayloadCurr >> 24) & 0x7f;
699 fCurrMinor = (*fPayloadCurr >> 17) & 0x7f;
700 fCurrAddHcWords = (*fPayloadCurr >> 14) & 0x7;
701 fCurrSm = (*fPayloadCurr >> 9) & 0x1f;
702 fCurrLayer = (*fPayloadCurr >> 6) & 0x7;
703 fCurrStack = (*fPayloadCurr >> 3) & 0x7;
704 fCurrSide = (*fPayloadCurr >> 2) & 0x1;
705 fCurrCheck = (*fPayloadCurr) & 0x3;
707 if (fCurrSm != (((Int_t) fCurrEquipmentId) - 1024) ||
708 fCurrStack != fCurrSlot ||
709 fCurrLayer != fCurrLink / 2 ||
710 fCurrSide != fCurrLink % 2) {
711 LinkError(kHCmismatch,
712 "HC: %i, %i, %i, %i\n 0x%08x 0x%08x 0x%08x 0x%08x",
713 fCurrSm, fCurrStack, fCurrLayer, fCurrSide,
714 fPayloadCurr[0], fPayloadCurr[1], fPayloadCurr[2], fPayloadCurr[3]);;
716 if (fCurrCheck != 0x1) {
717 LinkError(kHCcheckFailed);
720 if (fCurrAddHcWords > 0) {
721 fCurrNtimebins = (fPayloadCurr[1] >> 26) & 0x3f;
722 fCurrBC = (fPayloadCurr[1] >> 10) & 0xffff;
723 fCurrPtrgCnt = (fPayloadCurr[1] >> 6) & 0xf;
724 fCurrPtrgPhase = (fPayloadCurr[1] >> 2) & 0xf;
727 fPayloadCurr += 1 + fCurrAddHcWords;
729 return (fPayloadCurr - start);
732 Int_t AliTRDrawStream::ReadTPData(Int_t mode)
734 // testing of testpattern 1 to 3 (hardcoded), 0 missing
735 // evcnt checking missing
737 Int_t cpufromchannel[] = {0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3};
742 Int_t channelcount = 0;
744 UInt_t expadcval = 0;
746 Int_t lastmcmpos = -1;
747 Int_t lastrobpos = -1;
749 UInt_t* start = fPayloadCurr;
751 while (*(fPayloadCurr) != fgkDataEndmarker &&
752 fPayloadCurr - fPayloadStart < fPayloadSize - 1) {
754 // ----- Checking MCM Header -----
755 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
758 // ----- checking for proper readout order - ROB -----
759 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
760 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
765 fCurrRobPos = ROB(*fPayloadCurr);
767 // ----- checking for proper readout order - MCM -----
768 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
769 lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
774 fCurrMcmPos = MCM(*fPayloadCurr);
779 evcnt = 0x3f & *fPayloadCurr >> 26;
782 while (channelcount < 21) {
784 if (cpu != cpufromchannel[channelcount]) {
785 cpu = cpufromchannel[channelcount];
786 expadcval = (1 << 9) | (fCurrRobPos << 6) | (fCurrMcmPos << 2) | cpu;
791 if (channelcount % 2 == 0)
798 expword |= expadcval << 2;
799 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
800 expword |= expadcval << 12;
801 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
802 expword |= expadcval << 22;
803 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
805 else if (mode == 2) {
807 expword = ((0x3f & evcnt) << 26) | ((fCurrSm + 1) << 21) | ((fCurrLayer + 1) << 18) |
808 ((fCurrStack + 1) << 15) |
809 (fCurrRobPos << 12) | (fCurrMcmPos << 8) | (cpu << 6) | (wordcount + 1);
811 else if (mode == 3) {
813 expword = ((0xfff & evcnt) << 20) | (fCurrSm << 15) | (fCurrLink/2 << 12) | (fCurrStack << 9) |
814 (fCurrRobPos << 6) | (fCurrMcmPos << 2) | (cpu << 0);
818 LinkError(kTPmodeInvalid, "Just reading");
821 diff = *fPayloadCurr ^ expword;
823 MCMError(kTPmismatch,
824 "Seen 0x%08x, expected 0x%08x, diff: 0x%08x (0x%02x)",
825 *fPayloadCurr, expword, diff, 0xff & (diff | diff >> 8 | diff >> 16 | diff >> 24));;
833 // continue with next MCM
835 return fPayloadCurr - start;
839 Int_t AliTRDrawStream::ReadZSData()
841 // read the zs data from one link from the current reading position
843 UInt_t *start = fPayloadCurr;
846 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
847 Int_t channelcount = 0;
848 Int_t channelcountExp = 0;
849 Int_t channelcountMax = 0;
851 Int_t currentTimebin = 0;
854 Int_t lastmcmpos = -1;
855 Int_t lastrobpos = -1;
857 if (fCurrNtimebins != fNtimebins) {
859 LinkError(kNtimebinsChanged,
860 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
861 fNtimebins = fCurrNtimebins;
864 timebins = fNtimebins;
866 while (*(fPayloadCurr) != fgkDataEndmarker &&
867 fPayloadCurr - fPayloadStart < fPayloadSize) {
869 // ----- Checking MCM Header -----
870 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
871 UInt_t *startPosMCM = fPayloadCurr;
873 // ----- checking for proper readout order - ROB -----
874 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
875 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
880 fCurrRobPos = ROB(*fPayloadCurr);
882 // ----- checking for proper readout order - MCM -----
883 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) > lastmcmpos) {
884 lastmcmpos = GetMCMReadoutPos(lastmcmpos);
889 fCurrMcmPos = MCM(*fPayloadCurr);
891 if (EvNo(*fPayloadCurr) != evno) {
893 evno = EvNo(*fPayloadCurr);
895 MCMError(kPtrgCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
898 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
899 Int_t padcoloff = PadColOffset(*fPayloadCurr);
900 Int_t row = Row(*fPayloadCurr);
903 // ----- Reading ADC channels -----
904 AliDebug(2, Form("ADC mask: 0x%08x", *fPayloadCurr));
906 // ----- analysing the ADC mask -----
908 channelcountExp = GetNActiveChannelsFromMask(*fPayloadCurr);
909 channelcountMax = GetNActiveChannels(*fPayloadCurr);
910 Int_t channelmask = GetActiveChannels(*fPayloadCurr);
911 Int_t channelno = -1;
914 if (channelcountExp != channelcountMax) {
915 if (channelcountExp > channelcountMax) {
916 Int_t temp = channelcountExp;
917 channelcountExp = channelcountMax;
918 channelcountMax = temp;
920 while (channelcountExp < channelcountMax && channelcountExp < 21 &&
921 fPayloadCurr - fPayloadStart < fPayloadSize - 10 * channelcountExp - 1) {
922 MCMError(kAdcMaskInconsistent,
923 "Possible MCM-H: 0x%08x, possible ADC-mask: 0x%08x",
924 *(fPayloadCurr + 10 * channelcountExp),
925 *(fPayloadCurr + 10 * channelcountExp + 1) );
926 if (!CouldBeMCMhdr( *(fPayloadCurr + 10 * channelcountExp)) && !CouldBeADCmask( *(fPayloadCurr + 10 * channelcountExp + 1)))
932 MCMError(kAdcMaskInconsistent,
933 "Inconsistency in no. of active channels: Counter: %i, Mask: %i, chosen: %i!",
934 GetNActiveChannels(fPayloadCurr[-1]), GetNActiveChannelsFromMask(fPayloadCurr[-1]), channelcountExp);
936 AliDebug(2, Form("expecting %i active channels, timebins: %i", channelcountExp, fCurrNtimebins));
938 // ----- reading marked ADC channels -----
939 while (channelcount < channelcountExp && *(fPayloadCurr) != fgkDataEndmarker) {
942 while (channelno < 20 && (channelmask & 1 << channelno) == 0)
945 if (fCurrNtimebins > 30) {
946 currentTimebin = ((*fPayloadCurr >> 2) & 0x3f);
947 timebins = ((*fPayloadCurr >> 8) & 0xf) * 3;
954 AliDebug(2, Form("Now looking %i words", timebins / 3));
955 Int_t adccol = adccoloff - channelno;
956 Int_t padcol = padcoloff - channelno;
957 // if (adccol < 3 || adccol > 165)
958 // AliInfo(Form("writing channel %i of det %3i %i:%2i to adcrow/-col: %i/%i padcol: %i",
959 // channelno, fCurrHC/2, fCurrRobPos, fCurrMcmPos, row, adccol, padcol));
961 while (adcwc < timebins / 3 &&
962 *(fPayloadCurr) != fgkDataEndmarker &&
963 fPayloadCurr - fPayloadStart < fPayloadSize) {
964 int check = 0x3 & *fPayloadCurr;
965 if (channelno % 2 != 0) { // odd channel
966 if (check != 0x2 && channelno < 21) {
967 MCMError(kAdcCheckInvalid,
968 "%i for %2i. ADC word in odd channel %i",
969 check, adcwc+1, channelno);
972 else { // even channel
973 if (check != 0x3 && channelno < 21) {
974 MCMError(kAdcCheckInvalid,
975 "%i for %2i. ADC word in even channel %i",
976 check, adcwc+1, channelno);
980 // filling the actual timebin data
981 int tb2 = 0x3ff & *fPayloadCurr >> 22;
982 int tb1 = 0x3ff & *fPayloadCurr >> 12;
983 int tb0 = 0x3ff & *fPayloadCurr >> 2;
984 if (adcwc != 0 || fCurrNtimebins <= 30)
985 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
988 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
989 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
995 if (adcwc != timebins / 3)
996 MCMError(kAdcDataAbort);
999 if (padcol > 0 && padcol < 144) {
1000 fSignalIndex->AddIndexRC(row, padcol);
1006 if (fCurrSm > -1 && fCurrSm < 18) {
1007 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNChannels += channelcount;
1008 fStats.fStatsSector[fCurrSm].fNChannels += channelcount;
1010 if (channelcount != channelcountExp)
1011 MCMError(kAdcChannelsMiss);
1014 if (fCurrSm > -1 && fCurrSm < 18) {
1015 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNMCMs++;
1016 fStats.fStatsSector[fCurrSm].fNMCMs++;
1019 if (IsDumping() && DumpingMCM(fCurrHC/2, fCurrRobPos, fCurrMcmPos)) {
1020 DumpRaw(Form("Event %i: Det %3i ROB %i MCM %2i", fRawReader->GetEventIndex(), fCurrHC/2, fCurrRobPos, fCurrMcmPos),
1021 startPosMCM, fPayloadCurr - startPosMCM);
1024 // continue with next MCM
1027 // check for missing MCMs (if header suppression is inactive)
1028 if (((fCurrMajor & 0x1) == 0) && (mcmcount != mcmcountExp)) {
1029 LinkError(kMissMcmHeaders,
1030 "No. of MCM headers %i not as expected: %i",
1031 mcmcount, mcmcountExp);
1034 return (fPayloadCurr - start);
1037 Int_t AliTRDrawStream::ReadNonZSData()
1039 // read the non-zs data from one link from the current reading position
1041 UInt_t *start = fPayloadCurr;
1044 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
1045 Int_t channelcount = 0;
1046 Int_t channelcountExp = 0;
1048 Int_t currentTimebin = 0;
1051 Int_t lastmcmpos = -1;
1052 Int_t lastrobpos = -1;
1054 if (fCurrNtimebins != fNtimebins) {
1056 LinkError(kNtimebinsChanged,
1057 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
1058 fNtimebins = fCurrNtimebins;
1061 timebins = fNtimebins;
1063 while (*(fPayloadCurr) != fgkDataEndmarker &&
1064 fPayloadCurr - fPayloadStart < fPayloadSize - 2) {
1066 // ----- Checking MCM Header -----
1067 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
1069 // ----- checking for proper readout order - ROB -----
1070 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1071 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1074 ROBError(kPosUnexp);
1076 fCurrRobPos = ROB(*fPayloadCurr);
1078 // ----- checking for proper readout order - MCM -----
1079 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
1080 lastmcmpos = GetMCMReadoutPos(*fPayloadCurr);
1083 MCMError(kPosUnexp);
1085 fCurrMcmPos = MCM(*fPayloadCurr);
1087 if (EvNo(*fPayloadCurr) != evno) {
1089 evno = EvNo(*fPayloadCurr);
1091 MCMError(kPtrgCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
1096 channelcountExp = 21;
1099 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
1100 Int_t padcoloff = PadColOffset(*fPayloadCurr);
1101 Int_t row = Row(*fPayloadCurr);
1105 // ----- reading marked ADC channels -----
1106 while (channelcount < channelcountExp &&
1107 *(fPayloadCurr) != fgkDataEndmarker) {
1114 AliDebug(2, Form("Now looking %i words", timebins / 3));
1115 Int_t adccol = adccoloff - channelno;
1116 Int_t padcol = padcoloff - channelno;
1117 while (adcwc < timebins / 3 &&
1118 *(fPayloadCurr) != fgkDataEndmarker &&
1119 fPayloadCurr - fPayloadStart < fPayloadSize) {
1120 int check = 0x3 & *fPayloadCurr;
1121 if (channelno % 2 != 0) { // odd channel
1122 if (check != 0x2 && channelno < 21) {
1123 MCMError(kAdcCheckInvalid,
1124 "%i for %2i. ADC word in odd channel %i",
1125 check, adcwc+1, channelno);
1128 else { // even channel
1129 if (check != 0x3 && channelno < 21) {
1130 MCMError(kAdcCheckInvalid,
1131 "%i for %2i. ADC word in even channel %i",
1132 check, adcwc+1, channelno);
1136 // filling the actual timebin data
1137 int tb2 = 0x3ff & *fPayloadCurr >> 22;
1138 int tb1 = 0x3ff & *fPayloadCurr >> 12;
1139 int tb0 = 0x3ff & *fPayloadCurr >> 2;
1140 if (adcwc != 0 || fCurrNtimebins <= 30)
1141 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
1144 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
1145 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
1151 if (adcwc != timebins / 3)
1152 MCMError(kAdcDataAbort);
1155 if (padcol > 0 && padcol < 144) {
1156 fSignalIndex->AddIndexRC(row, padcol);
1162 if (channelcount != channelcountExp)
1163 MCMError(kAdcChannelsMiss);
1165 // continue with next MCM
1168 // check for missing MCMs (if header suppression is inactive)
1169 if (mcmcount != mcmcountExp) {
1170 LinkError(kMissMcmHeaders,
1171 "%i not as expected: %i", mcmcount, mcmcountExp);
1174 return (fPayloadCurr - start);
1177 Int_t AliTRDrawStream::SeekNextLink()
1179 UInt_t *start = fPayloadCurr;
1181 // read until data endmarkers
1182 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1183 *fPayloadCurr != fgkDataEndmarker)
1186 // read all data endmarkers
1187 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1188 *fPayloadCurr == fgkDataEndmarker)
1191 return (fPayloadCurr - start);
1194 Bool_t AliTRDrawStream::ConnectTracklets(TTree *trklTree)
1196 fTrackletTree = trklTree;
1200 if (!fTrackletTree->GetBranch("hc"))
1201 fTrackletTree->Branch("hc", &fCurrHC, "hc/I");
1203 fTrackletTree->SetBranchAddress("hc", &fCurrHC);
1205 if (!fTrackletTree->GetBranch("trkl"))
1206 fTrackletTree->Branch("trkl", &fTrackletArray);
1208 fTrackletTree->SetBranchAddress("trkl", &fTrackletArray);
1214 void AliTRDrawStream::EquipmentError(ErrorCode_t err, const char *const msg, ...)
1216 // register error according to error code on equipment level
1217 // and return the corresponding error message
1219 fLastError.fSector = fCurrEquipmentId - 1024;
1220 fLastError.fStack = -1;
1221 fLastError.fLink = -1;
1222 fLastError.fRob = -1;
1223 fLastError.fMcm = -1;
1224 fLastError.fError = err;
1225 (this->*fStoreError)();
1228 if (fgErrorDebugLevel[err] > 10)
1229 AliDebug(fgErrorDebugLevel[err],
1230 Form("Event %6i: Eq. %2d - %s : %s",
1231 fRawReader->GetEventIndex(), fCurrEquipmentId, fgErrorMessages[err],
1232 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1234 AliError(Form("Event %6i: Eq. %2d - %s : %s",
1235 fRawReader->GetEventIndex(), fCurrEquipmentId, fgErrorMessages[err],
1236 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1237 fErrorFlags |= fgErrorBehav[err];
1241 void AliTRDrawStream::StackError(ErrorCode_t err, const char *const msg, ...)
1243 // register error according to error code on stack level
1244 // and return the corresponding error message
1246 fLastError.fSector = fCurrEquipmentId - 1024;
1247 fLastError.fStack = fCurrSlot;
1248 fLastError.fLink = -1;
1249 fLastError.fRob = -1;
1250 fLastError.fMcm = -1;
1251 fLastError.fError = err;
1252 (this->*fStoreError)();
1255 if (fgErrorDebugLevel[err] > 0)
1256 AliDebug(fgErrorDebugLevel[err],
1257 Form("Event %6i: Eq. %2d S %i - %s : %s",
1258 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgErrorMessages[err],
1259 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1261 AliError(Form("Event %6i: Eq. %2d S %i - %s : %s",
1262 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgErrorMessages[err],
1263 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1264 fErrorFlags |= fgErrorBehav[err];
1268 void AliTRDrawStream::LinkError(ErrorCode_t err, const char *const msg, ...)
1270 // register error according to error code on link level
1271 // and return the corresponding error message
1273 fLastError.fSector = fCurrEquipmentId - 1024;
1274 fLastError.fStack = fCurrSlot;
1275 fLastError.fLink = fCurrLink;
1276 fLastError.fRob = -1;
1277 fLastError.fMcm = -1;
1278 fLastError.fError = err;
1279 (this->*fStoreError)();
1282 if (fgErrorDebugLevel[err] > 0)
1283 AliDebug(fgErrorDebugLevel[err],
1284 Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
1285 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgErrorMessages[err],
1286 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1288 AliError(Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
1289 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgErrorMessages[err],
1290 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1291 fErrorFlags |= fgErrorBehav[err];
1295 void AliTRDrawStream::ROBError(ErrorCode_t err, const char *const msg, ...)
1297 // register error according to error code on ROB level
1298 // and return the corresponding error message
1300 fLastError.fSector = fCurrEquipmentId - 1024;
1301 fLastError.fStack = fCurrSlot;
1302 fLastError.fLink = fCurrLink;
1303 fLastError.fRob = fCurrRobPos;
1304 fLastError.fMcm = -1;
1305 fLastError.fError = err;
1306 (this->*fStoreError)();
1309 if (fgErrorDebugLevel[err] > 0)
1310 AliDebug(fgErrorDebugLevel[err],
1311 Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : %s",
1312 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgErrorMessages[err],
1313 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1315 AliError(Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : %s",
1316 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgErrorMessages[err],
1317 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1318 fErrorFlags |= fgErrorBehav[err];
1322 void AliTRDrawStream::MCMError(ErrorCode_t err, const char *const msg, ...)
1324 // register error according to error code on MCM level
1325 // and return the corresponding error message
1327 fLastError.fSector = fCurrEquipmentId - 1024;
1328 fLastError.fStack = fCurrSlot;
1329 fLastError.fLink = fCurrLink;
1330 fLastError.fRob = fCurrRobPos;
1331 fLastError.fMcm = fCurrMcmPos;
1332 fLastError.fError = err;
1333 (this->*fStoreError)();
1336 if (fgErrorDebugLevel[err] > 0)
1337 AliDebug(fgErrorDebugLevel[err],
1338 Form("Event %6i: Eq. %2d S %i l %2i ROB %i MCM %2i - %s : %s",
1339 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fCurrMcmPos, fgErrorMessages[err],
1340 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1342 AliError(Form("Event %6i: Eq. %2d S %i l %2i ROB %i MCM %2i - %s : %s",
1343 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fCurrMcmPos, fgErrorMessages[err],
1344 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
1345 fErrorFlags |= fgErrorBehav[err];
1348 const char* AliTRDrawStream::GetErrorMessage(ErrorCode_t errCode)
1350 // return the error message for the given error code
1352 if (errCode > 0 && errCode < kLastErrorCode)
1353 return fgErrorMessages[errCode];
1358 void AliTRDrawStream::AliTRDrawStats::ClearStats()
1360 // clear statistics (includes clearing sector-wise statistics)
1363 for (Int_t iSector = 0; iSector < 18; iSector++) {
1364 fStatsSector[iSector].ClearStats();
1369 void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::ClearStats()
1371 // clear statistics (includes clearing HC-wise statistics)
1379 for (Int_t iHC = 0; iHC < 60; iHC++) {
1380 fStatsHC[iHC].ClearStats();
1384 void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::AliTRDrawStatsHC::ClearStats()
1395 void AliTRDrawStream::SetDumpMCM(Int_t det, Int_t rob, Int_t mcm, Bool_t dump)
1397 // mark MCM for dumping of raw data
1400 fDumpMCM[fNDumpMCMs++] = (det << 7) | (rob << 4) | mcm;
1404 for (iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
1405 if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
1410 for ( ; iMCM < fNDumpMCMs; iMCM++) {
1411 fDumpMCM[iMCM] = fDumpMCM[iMCM+1];
1416 Bool_t AliTRDrawStream::DumpingMCM(Int_t det, Int_t rob, Int_t mcm)
1418 // check if MCM data should be dumped
1420 for (Int_t iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
1421 if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
1428 void AliTRDrawStream::DumpRaw(TString title, UInt_t *start, Int_t length)
1434 for ( ; pos+3 < length; pos += 4) {
1435 title += Form("0x%08x 0x%08x 0x%08x 0x%08x\n",
1436 start[pos+0], start[pos+1], start[pos+2], start[pos+3]);
1438 for ( ; pos < length; pos++) {
1439 title += Form("0x%08x ", start[pos]);
1444 AliTRDrawStream::AliTRDrawStreamError::AliTRDrawStreamError(Int_t error, Int_t sector, Int_t stack, Int_t link, Int_t rob, Int_t mcm) :