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 ////////////////////////////////////////////////////////////////////////////
25 #include "TClonesArray.h"
29 #include "AliRawReader.h"
30 #include "AliTRDdigitsManager.h"
31 #include "AliTRDdigitsParam.h"
32 #include "AliTRDtrapConfig.h"
33 #include "AliTRDarrayADC.h"
34 #include "AliTRDarrayDictionary.h"
35 #include "AliTRDSignalIndex.h"
36 #include "AliTRDtrackletWord.h"
38 #include "AliTRDrawStream.h"
41 #include "AliRunLoader.h"
43 ClassImp(AliTRDrawStream)
45 // some static information
46 const Int_t AliTRDrawStream::fgkMcmOrder[] = {12, 13, 14, 15,
50 const Int_t AliTRDrawStream::fgkRobOrder [] = {0, 1, 2, 3};
51 const Int_t AliTRDrawStream::fgkNlinks = 12;
52 const Int_t AliTRDrawStream::fgkNstacks = 5;
53 const UInt_t AliTRDrawStream::fgkDataEndmarker = 0x00000000;
54 const UInt_t AliTRDrawStream::fgkTrackletEndmarker = 0x10001000;
56 char* AliTRDrawStream::fgErrorMessages[] = {
58 "Link monitor active",
59 "Pretrigger counter mismatch",
60 "not a TRD equipment (1024-1041)",
61 "Invalid Stack header",
62 "Invalid detector number",
63 "No digits could be retrieved from the digitsmanager",
65 "HC check bits wrong",
66 "Unexpected position in readout stream",
67 "Invalid testpattern mode",
68 "Testpattern mismatch",
69 "Number of timebins changed",
70 "ADC mask inconsistent",
71 "ADC check bits invalid",
73 "Missing expected ADC channels",
77 AliTRDrawStream::AliTRDrawStream(AliRawReader *rawReader) :
78 fRawReader(rawReader),
93 fCurrSmuIndexHeaderSize(0),
94 fCurrSmuIndexHeaderVersion(0),
96 fCurrTrackletEnable(0),
98 fCurrStackIndexWord(0x0),
99 fCurrStackHeaderSize(0x0),
100 fCurrStackHeaderVersion(0x0),
102 fCurrCleanCheckout(0x0),
105 fCurrLinkMonitorFlags(0x0),
106 fCurrLinkDataTypeFlags(0x0),
107 fCurrLinkDebugFlags(0x0),
127 // default constructor
129 fCurrStackIndexWord = new UInt_t[fgkNstacks];
130 fCurrStackHeaderSize = new UInt_t[fgkNstacks];
131 fCurrStackHeaderVersion = new UInt_t[fgkNstacks];
132 fCurrLinkMask = new UInt_t[fgkNstacks];
133 fCurrCleanCheckout = new UInt_t[fgkNstacks];
134 fCurrBoardId = new UInt_t[fgkNstacks];
135 fCurrHwRev = new UInt_t[fgkNstacks];
136 fCurrLinkMonitorFlags = new UInt_t[fgkNstacks * fgkNlinks];
137 fCurrLinkDataTypeFlags = new UInt_t[fgkNstacks * fgkNlinks];
138 fCurrLinkDebugFlags = new UInt_t[fgkNstacks * fgkNlinks];
140 // preparing TClonesArray
141 fTrackletArray = new TClonesArray("AliTRDtrackletWord", 256);
143 // setting up the error tree
144 fErrors = new TTree("errorStats", "Error statistics");
145 fErrors->SetDirectory(0x0);
146 fErrors->Branch("error", &fLastError, "sector/I:stack:link:error:rob:mcm");
147 fErrors->SetCircular(1000);
150 AliTRDrawStream::~AliTRDrawStream()
156 delete [] fCurrStackIndexWord;
157 delete [] fCurrStackHeaderSize;
158 delete [] fCurrStackHeaderVersion;
159 delete [] fCurrLinkMask;
160 delete [] fCurrCleanCheckout;
161 delete [] fCurrBoardId;
162 delete [] fCurrHwRev;
163 delete [] fCurrLinkMonitorFlags;
164 delete [] fCurrLinkDataTypeFlags;
165 delete [] fCurrLinkDebugFlags;
168 Bool_t AliTRDrawStream::ReadEvent(TTree *trackletTree)
170 // read the current event from the raw reader and fill it to the digits manager
173 AliError("No raw reader available");
178 fTrackletTree = trackletTree;
180 if (!fTrackletTree->GetBranch("trkl")) {
181 fTrackletTree->Branch("hc", &fCurrHC, "hc/I");
182 fTrackletTree->Branch("trkl", &fTrackletArray);
185 fTrackletTree->SetBranchAddress("hc", &fCurrHC);
186 fTrackletTree->SetBranchAddress("trkl", &fTrackletArray);
193 // loop over all DDLs
194 // data starts with GTU payload, i.e. SMU index word
195 UChar_t *buffer = 0x0;
197 while (fRawReader->ReadNextData(buffer)) {
199 fCurrEquipmentId = fRawReader->GetEquipmentId();
200 AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
202 if (fCurrEquipmentId < 1024 || fCurrEquipmentId > 1041) {
203 AliError(EquipmentError(kNonTrdEq, "Skipping"));
207 // setting the pointer to data and current reading position
208 fPayloadCurr = fPayloadStart = (UInt_t*) (buffer);
209 fPayloadSize = fRawReader->GetDataSize() / sizeof(UInt_t);
210 AliDebug(2, Form("Read buffer of size: %i", fRawReader->GetDataSize()));
212 // read SMU index header
213 if (ReadSmHeader() < 0) {
214 AliError(Form("Reading SMU header failed, skipping this DDL %i", fCurrEquipmentId));
218 // read stack index header
219 for (Int_t iStack = 0; iStack < 5; iStack++) {
220 if (fCurrStackMask & (1 << iStack) != 0)
221 ReadStackIndexHeader(iStack);
224 for (Int_t iStack = 0; iStack < 5; iStack++) {
226 if (fCurrStackMask & (1 << fCurrSlot) == 0)
229 AliDebug(2, Form("Stack %i, Link mask: 0x%02x", fCurrSlot, fCurrLinkMask[fCurrSlot]));
230 for (Int_t iLink = 0; iLink < 12; iLink++) {
232 fCurrHC = fCurrSm * 60 + fCurrSlot * 12 + iLink;
233 if ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) == 0)
236 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0)
237 LinkError(kLinkMonitor);
238 AliDebug(2, Form("Payload for S%il%i", fCurrSlot, fCurrLink));
239 for (Int_t i = 0; i < 10; i++) //???
240 AliDebug(5, Form("%3i: 0x%08x 0x%08x 0x%08x 0x%08x", i*4,
241 fPayloadCurr[4*i], fPayloadCurr[4*i+1],
242 fPayloadCurr[4*i+2], fPayloadCurr[4*i+3]));
244 // read the data from one HC
247 // read all data endmarkers
248 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
249 *fPayloadCurr == fgkDataEndmarker)
259 Bool_t AliTRDrawStream::NextDDL()
261 // continue reading with the next equipment
266 fCurrEquipmentId = 0;
270 UChar_t *buffer = 0x0;
272 while (fRawReader->ReadNextData(buffer)) {
274 fCurrEquipmentId = fRawReader->GetEquipmentId();
275 AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
277 if (fCurrEquipmentId < 1024 || fCurrEquipmentId > 1041) {
278 AliError(EquipmentError(kNonTrdEq, "Skipping"));
282 // setting the pointer to data and current reading position
283 fPayloadCurr = fPayloadStart = (UInt_t*) (buffer);
284 fPayloadSize = fRawReader->GetDataSize() / sizeof(UInt_t);
285 AliDebug(2, Form("Read buffer of size: %i", fRawReader->GetDataSize()));
287 // read SMU index header
288 if (ReadSmHeader() < 0) {
289 AliError(Form("Reading SMU header failed, skipping this DDL %i", fCurrEquipmentId));
293 // read stack index header
294 for (Int_t iStack = 0; iStack < 5; iStack++) {
295 if (fCurrStackMask & (1 << iStack) != 0) {
296 ReadStackIndexHeader(iStack);
306 Int_t AliTRDrawStream::NextChamber(AliTRDdigitsManager *digMgr, UInt_t ** /* trackletContainer */, UShort_t ** /* errorContainer */)
308 // read the data for the next chamber
309 // in case you only want to read the data of a single chamber
310 // to read all data ReadEvent(...) is recommended
312 fDigitsManager = digMgr;
315 // tracklet output preparation
317 //AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets");
319 // fTrackletTree = trklLoader->Tree();
322 AliError("No raw reader available");
326 if (fCurrSlot < 0 || fCurrSlot >= 5) {
333 AliDebug(2, Form("Stack %i, Link %i, mask: 0x%02x", fCurrSlot, fCurrLink, fCurrLinkMask[fCurrSlot]));
334 fCurrHC = (fCurrEquipmentId - 1024) * 60 + fCurrSlot * 12 + fCurrLink;
336 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0)
337 LinkError(kLinkMonitor);
339 // read the data from one HC
342 // read all data endmarkers
343 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
344 *fPayloadCurr == fgkDataEndmarker)
347 if (fCurrLink % 2 == 0) {
348 // if we just read the A-side HC then also check the B-side
350 if (fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) {
352 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
353 *fPayloadCurr == fgkDataEndmarker)
361 if (fCurrLink > 11) {
365 } while ((fCurrSlot < 5) &&
366 ((fCurrStackMask & (1 << fCurrSlot) == 0) ||
367 (fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) == 0));
369 return (fCurrSm * 30 + fCurrStack * 6 + fCurrLayer);
373 Int_t AliTRDrawStream::ReadSmHeader()
375 // read the SMU index header at the current reading position
376 // and store the information in the corresponding variables
378 if (fPayloadCurr - fPayloadStart >= fPayloadSize - 1) {
379 AliError(EquipmentError(kUnknown, "SM Header incomplete"));
383 fCurrSmuIndexHeaderSize = ((*fPayloadCurr) >> 16) & 0xffff;
384 fCurrSmuIndexHeaderVersion = ((*fPayloadCurr) >> 12) & 0xf;
385 fCurrTrackEnable = ((*fPayloadCurr) >> 6) & 0x1;
386 fCurrTrackletEnable = ((*fPayloadCurr) >> 5) & 0x1;
387 fCurrStackMask = ((*fPayloadCurr) ) & 0x1f;
389 AliDebug(5, Form("SMU header: size: %i, version: %i, track enable: %i, tracklet enable: %i, stack mask: %2x",
390 fCurrSmuIndexHeaderSize,
391 fCurrSmuIndexHeaderVersion,
396 fPayloadCurr += fCurrSmuIndexHeaderSize + 1;
398 return fCurrSmuIndexHeaderSize + 1;
401 Int_t AliTRDrawStream::ReadStackIndexHeader(Int_t stack)
403 // read the stack index header
404 // and store the information in the corresponding variables
406 fCurrStackIndexWord[stack] = *fPayloadCurr;
407 fCurrStackHeaderSize[stack] = (((*fPayloadCurr) >> 16) & 0xffff) + 1;
408 fCurrStackHeaderVersion[stack] = ((*fPayloadCurr) >> 12) & 0xf;
409 fCurrLinkMask[stack] = (*fPayloadCurr) & 0xfff;
411 if (fPayloadCurr - fPayloadStart >= fPayloadSize - fCurrStackHeaderSize[stack]) {
412 AliError(StackError(kStackHeaderInvalid, "Stack index header aborted"));
416 switch (fCurrStackHeaderVersion[stack]) {
418 if (fCurrStackHeaderSize[stack] < 8) {
419 AliError(StackError(kStackHeaderInvalid, "Stack header smaller than expected!"));
423 fCurrCleanCheckout[stack] = fPayloadCurr[1] & 0x1;
424 fCurrBoardId[stack] = (fPayloadCurr[1] >> 8) & 0xff;
425 fCurrHwRev[stack] = (fPayloadCurr[1] >> 16) & 0xffff;
427 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
429 fCurrLinkMonitorFlags [stack*fgkNlinks + iLayer*2] = fPayloadCurr[iLayer+2] & 0xf;
430 fCurrLinkDataTypeFlags [stack*fgkNlinks + iLayer*2] = (fPayloadCurr[iLayer+2] >> 4) & 0x3;
431 fCurrLinkDebugFlags [stack*fgkNlinks + iLayer*2] = (fPayloadCurr[iLayer+2] >> 12) & 0xf;
433 fCurrLinkMonitorFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 16) & 0xf;
434 fCurrLinkDataTypeFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 20) & 0x3;
435 fCurrLinkDebugFlags [stack*fgkNlinks + iLayer*2 + 1] = (fPayloadCurr[iLayer+2] >> 28) & 0xf;
440 AliError(StackError(kStackHeaderInvalid, Form("Invalid Stack Index Header version %x",
441 fCurrStackHeaderVersion[stack])));
444 fPayloadCurr += fCurrStackHeaderSize[stack];
446 return fCurrStackHeaderSize[stack];
449 Int_t AliTRDrawStream::ReadLinkData()
451 // read the data in one link (one HC) until the data endmarker is reached
455 count += ReadTracklets();
457 count += ReadHcHeader();
459 Int_t det = fCurrSm * 30 + fCurrStack * 6 + fCurrLayer;
461 if (det > -1 && det < 540) {
463 if ((fAdcArray = fDigitsManager->GetDigits(det))) {
464 //fAdcArray->Expand();
465 if (fAdcArray->GetNtime() != fCurrNtimebins)
466 fAdcArray->Allocate(16, 144, fCurrNtimebins);
469 AliError(LinkError(kNoDigits));
473 fDigitsParam = fDigitsManager->GetDigitsParam();
476 fDigitsParam->SetPretriggerPhase(det, fCurrPtrgPhase);
477 fDigitsParam->SetNTimeBins(det, fCurrNtimebins);
478 fDigitsParam->SetADCbaseline(det, 10);
481 if (fDigitsManager->UsesDictionaries()) {
482 fDigitsManager->GetDictionary(det, 0)->Reset();
483 fDigitsManager->GetDictionary(det, 1)->Reset();
484 fDigitsManager->GetDictionary(det, 2)->Reset();
487 if ((fSignalIndex = fDigitsManager->GetIndexes(det))) {
488 fSignalIndex->SetSM(fCurrSm);
489 fSignalIndex->SetStack(fCurrStack);
490 fSignalIndex->SetLayer(fCurrLayer);
491 fSignalIndex->SetDetNumber(det);
492 if (!fSignalIndex->IsAllocated())
493 fSignalIndex->Allocate(16, 144, fCurrNtimebins);
496 // ----- check which kind of data -----
497 if (fCurrMajor & 0x40) {
498 if ((fCurrMajor & 0x7) == 0x7) {
499 AliDebug(1, "This is a config event");
500 UInt_t *startPos = fPayloadCurr;
501 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
502 *fPayloadCurr != fgkDataEndmarker)
504 count += fPayloadCurr - startPos;
506 // feeding TRAP config
507 AliTRDtrapConfig *trapcfg = AliTRDtrapConfig::Instance();
508 trapcfg->ReadPackedConfig(fCurrHC, startPos, fPayloadCurr - startPos);
511 Int_t tpmode = fCurrMajor & 0x7;
512 AliDebug(1, Form("Checking testpattern (mode %i) data", tpmode));
516 else if (fCurrMajor & 0x20) {
517 AliDebug(1, "This is a zs event");
518 count += ReadZSData();
521 AliDebug(1, "This is a nozs event");
522 count += ReadNonZSData();
526 AliError(LinkError(kInvalidDetector, Form("%i", det)));
527 UInt_t *startPos = fPayloadCurr;
528 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
529 *fPayloadCurr != fgkDataEndmarker)
531 count += fPayloadCurr - startPos;
537 Int_t AliTRDrawStream::ReadTracklets()
539 // read the tracklets from one HC
541 fTrackletArray->Clear();
543 UInt_t *start = fPayloadCurr;
544 while (*(fPayloadCurr) != fgkTrackletEndmarker &&
545 fPayloadCurr - fPayloadStart < fPayloadSize) {
547 new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(*(fPayloadCurr));
552 if (fTrackletArray->GetEntriesFast() > 0) {
553 AliDebug(1, Form("Found %i tracklets in %i %i %i (ev. %i)", fTrackletArray->GetEntriesFast(),
554 fCurrSm, fCurrSlot, fCurrLink, fRawReader->GetEventIndex()));
556 fTrackletTree->Fill();
559 // loop over remaining tracklet endmarkers
560 while ((*(fPayloadCurr) == fgkTrackletEndmarker &&
561 fPayloadCurr - fPayloadStart < fPayloadSize))
564 return (fPayloadCurr - start) / sizeof(UInt_t);
567 Int_t AliTRDrawStream::ReadHcHeader()
569 // read and parse the HC header of one HC
570 // and store the information in the corresponding variables
572 UInt_t *start = fPayloadCurr;
573 // check not to be at the data endmarker
574 if (*fPayloadCurr == fgkDataEndmarker)
577 fCurrSpecial = (*fPayloadCurr >> 31) & 0x1;
578 fCurrMajor = (*fPayloadCurr >> 24) & 0x7f;
579 fCurrMinor = (*fPayloadCurr >> 17) & 0x7f;
580 fCurrAddHcWords = (*fPayloadCurr >> 14) & 0x7;
581 fCurrSm = (*fPayloadCurr >> 9) & 0x1f;
582 fCurrLayer = (*fPayloadCurr >> 6) & 0x7;
583 fCurrStack = (*fPayloadCurr >> 3) & 0x7;
584 fCurrSide = (*fPayloadCurr >> 2) & 0x1;
585 fCurrCheck = (*fPayloadCurr) & 0x3;
587 if (fCurrSm != (((Int_t) fCurrEquipmentId) - 1024) ||
588 fCurrStack != fCurrSlot ||
589 fCurrLayer != fCurrLink / 2 ||
590 fCurrSide != fCurrLink % 2) {
591 AliError(LinkError(kHCmismatch,
592 Form("HC: %i, %i, %i, %i\n 0x%08x 0x%08x 0x%08x 0x%08x",
593 fCurrSm, fCurrStack, fCurrLayer, fCurrSide,
594 fPayloadCurr[0], fPayloadCurr[1], fPayloadCurr[2], fPayloadCurr[3]) ));
596 if (fCurrCheck != 0x1) {
597 AliError(LinkError(kHCcheckFailed));
600 if (fCurrAddHcWords > 0) {
601 fCurrNtimebins = (fPayloadCurr[1] >> 26) & 0x3f;
602 fCurrBC = (fPayloadCurr[1] >> 10) & 0xffff;
603 fCurrPtrgCnt = (fPayloadCurr[1] >> 6) & 0xf;
604 fCurrPtrgPhase = (fPayloadCurr[1] >> 2) & 0xf;
607 fPayloadCurr += 1 + fCurrAddHcWords;
609 return (fPayloadCurr - start) / sizeof(UInt_t);
612 Int_t AliTRDrawStream::ReadTPData(Int_t mode)
614 // testing of testpattern 1 to 3 (hardcoded), 0 missing
615 // evcnt checking missing
617 Int_t cpufromchannel[] = {0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3};
622 Int_t channelcount = 0;
624 UInt_t expadcval = 0;
626 Int_t lastmcmpos = -1;
627 Int_t lastrobpos = -1;
629 UInt_t* start = fPayloadCurr;
631 while (*(fPayloadCurr) != fgkDataEndmarker &&
632 fPayloadCurr - fPayloadStart < fPayloadSize - 1) {
634 // ----- Checking MCM Header -----
635 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
638 // ----- checking for proper readout order - ROB -----
639 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
640 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
643 AliError(ROBError(kPosUnexp));
645 fCurrRobPos = ROB(*fPayloadCurr);
647 // ----- checking for proper readout order - MCM -----
648 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
649 lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
652 AliError(MCMError(kPosUnexp));
654 fCurrMcmPos = MCM(*fPayloadCurr);
659 evcnt = 0x3f & *fPayloadCurr >> 26;
662 while (channelcount < 21) {
664 if (cpu != cpufromchannel[channelcount]) {
665 cpu = cpufromchannel[channelcount];
666 expadcval = (1 << 9) | (fCurrRobPos << 6) | (fCurrMcmPos << 2) | cpu;
671 if (channelcount % 2 == 0)
678 expword |= expadcval << 2;
679 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
680 expword |= expadcval << 12;
681 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
682 expword |= expadcval << 22;
683 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
685 else if (mode == 2) {
687 expword = ((0x3f & evcnt) << 26) | ((fCurrSm + 1) << 21) | ((fCurrLayer + 1) << 18) |
688 ((fCurrStack + 1) << 15) |
689 (fCurrRobPos << 12) | (fCurrMcmPos << 8) | (cpu << 6) | (wordcount + 1);
691 else if (mode == 3) {
693 expword = ((0xfff & evcnt) << 20) | (fCurrSm << 15) | (fCurrLink/2 << 12) | (fCurrStack << 9) |
694 (fCurrRobPos << 6) | (fCurrMcmPos << 2) | (cpu << 0);
698 AliError(LinkError(kTPmodeInvalid, "Just reading"));
701 diff = *fPayloadCurr ^ expword;
703 AliError(MCMError(kTPmismatch,
704 Form("Seen 0x%08x, expected 0x%08x, diff: 0x%08x (0x%02x)",
705 *fPayloadCurr, expword, diff, 0xff & (diff | diff >> 8 | diff >> 16 | diff >> 24)) ));
713 // continue with next MCM
715 return fPayloadCurr - start;
719 Int_t AliTRDrawStream::ReadZSData()
721 // read the zs data from one link from the current reading position
723 UInt_t *start = fPayloadCurr;
726 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
727 Int_t channelcount = 0;
728 Int_t channelcountExp = 0;
729 Int_t channelcountMax = 0;
731 Int_t currentTimebin = 0;
734 Int_t lastmcmpos = -1;
735 Int_t lastrobpos = -1;
737 if (fCurrNtimebins != fNtimebins) {
739 AliError(LinkError(kNtimebinsChanged,
740 Form("No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins) ));
741 fNtimebins = fCurrNtimebins;
744 timebins = fNtimebins;
746 while (*(fPayloadCurr) != fgkDataEndmarker &&
747 fPayloadCurr - fPayloadStart < fPayloadSize) {
749 // ----- Checking MCM Header -----
750 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
752 // ----- checking for proper readout order - ROB -----
753 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
754 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
757 AliError(ROBError(kPosUnexp));
759 fCurrRobPos = ROB(*fPayloadCurr);
761 // ----- checking for proper readout order - MCM -----
762 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) > lastmcmpos) {
763 lastmcmpos = GetMCMReadoutPos(lastmcmpos);
766 AliError(MCMError(kPosUnexp));
768 fCurrMcmPos = MCM(*fPayloadCurr);
770 if (EvNo(*fPayloadCurr) != evno) {
772 evno = EvNo(*fPayloadCurr);
774 AliError(MCMError(kPtrgCntMismatch,
775 Form("%i <-> %i", evno, EvNo(*fPayloadCurr)) ));
778 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
779 Int_t padcoloff = PadColOffset(*fPayloadCurr);
780 Int_t row = Row(*fPayloadCurr);
783 // ----- Reading ADC channels -----
784 AliDebug(2, Form("ADC mask: 0x%08x", *fPayloadCurr));
786 // ----- analysing the ADC mask -----
788 channelcountExp = GetNActiveChannelsFromMask(*fPayloadCurr);
789 channelcountMax = GetNActiveChannels(*fPayloadCurr);
790 Int_t channelmask = GetActiveChannels(*fPayloadCurr);
791 Int_t channelno = -1;
794 if (channelcountExp != channelcountMax) {
795 if (channelcountExp > channelcountMax) {
796 Int_t temp = channelcountExp;
797 channelcountExp = channelcountMax;
798 channelcountMax = temp;
800 while (channelcountExp < channelcountMax && channelcountExp < 21 &&
801 fPayloadCurr - fPayloadStart < fPayloadSize - 10 * channelcountExp - 1) {
802 AliError(MCMError(kAdcMaskInconsistent,
803 Form("Possible MCM-H: 0x%08x, possible ADC-mask: 0x%08x",
804 *(fPayloadCurr + 10 * channelcountExp),
805 *(fPayloadCurr + 10 * channelcountExp + 1) ) ));
806 if (!CouldBeMCMhdr( *(fPayloadCurr + 10 * channelcountExp)) && !CouldBeADCmask( *(fPayloadCurr + 10 * channelcountExp + 1)))
812 AliError(MCMError(kAdcMaskInconsistent,
813 Form("Inconsistency in no. of active channels: Counter: %i, Mask: %i, chosen: %i!",
814 GetNActiveChannels(fPayloadCurr[-1]), GetNActiveChannelsFromMask(fPayloadCurr[-1]), channelcountExp) ));
816 AliDebug(2, Form("expecting %i active channels, timebins: %i", channelcountExp, fCurrNtimebins));
818 // ----- reading marked ADC channels -----
819 while (channelcount < channelcountExp && *(fPayloadCurr) != fgkDataEndmarker) {
822 while (channelno < 21 && (channelmask & 1 << channelno) == 0)
825 if (fCurrNtimebins > 30) {
826 currentTimebin = ((*fPayloadCurr >> 2) & 0x3f);
827 timebins = ((*fPayloadCurr >> 8) & 0xf) * 3;
834 AliDebug(2, Form("Now looking %i words", timebins / 3));
835 Int_t adccol = adccoloff - channelno;
836 Int_t padcol = padcoloff - channelno;
837 while (adcwc < timebins / 3 &&
838 *(fPayloadCurr) != fgkDataEndmarker &&
839 fPayloadCurr - fPayloadStart < fPayloadSize) {
840 int check = 0x3 & *fPayloadCurr;
841 if (channelno % 2 != 0) { // odd channel
842 if (check != 0x2 && channelno < 21) {
843 AliError(MCMError(kAdcCheckInvalid,
844 Form("Invalid check bits: %i for %2i. ADC word in odd channel: %i",
845 check, adcwc+1, channelno) ));
848 else { // even channel
849 if (check != 0x3 && channelno < 21) {
850 AliError(MCMError(kAdcCheckInvalid,
851 Form("Invalid check bits: %i for %2i. ADC word in even channel: %i",
852 check, adcwc+1, channelno) ));
856 // filling the actual timebin data
857 int tb2 = 0x3ff & *fPayloadCurr >> 22;
858 int tb1 = 0x3ff & *fPayloadCurr >> 12;
859 int tb0 = 0x3ff & *fPayloadCurr >> 2;
860 if (adcwc != 0 || fCurrNtimebins <= 30)
861 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
864 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
865 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
871 if (adcwc != timebins / 3)
872 AliError(MCMError(kAdcDataAbort));
875 if (padcol > 0 && padcol < 144) {
876 fSignalIndex->AddIndexRC(row, padcol);
882 if (channelcount != channelcountExp)
883 AliError(MCMError(kAdcChannelsMiss));
886 // continue with next MCM
889 // check for missing MCMs (if header suppression is inactive)
890 if (fCurrMajor & 0x1 == 0 && mcmcount != mcmcountExp) {
891 AliError(LinkError(kMissMcmHeaders,
892 Form("No. of MCM headers %i not as expected: %i",
893 mcmcount, mcmcountExp) ));
896 return (fPayloadCurr - start);
899 Int_t AliTRDrawStream::ReadNonZSData()
901 // read the non-zs data from one link from the current reading position
903 UInt_t *start = fPayloadCurr;
906 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
907 Int_t channelcount = 0;
908 Int_t channelcountExp = 0;
910 Int_t currentTimebin = 0;
913 Int_t lastmcmpos = -1;
914 Int_t lastrobpos = -1;
916 if (fCurrNtimebins != fNtimebins) {
918 AliError(LinkError(kNtimebinsChanged,
919 Form("No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins) ));
920 fNtimebins = fCurrNtimebins;
923 timebins = fNtimebins;
925 while (*(fPayloadCurr) != fgkDataEndmarker &&
926 fPayloadCurr - fPayloadStart < fPayloadSize - 2) {
928 // ----- Checking MCM Header -----
929 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
931 // ----- checking for proper readout order - ROB -----
932 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
933 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
936 AliError(ROBError(kPosUnexp));
938 fCurrRobPos = ROB(*fPayloadCurr);
940 // ----- checking for proper readout order - MCM -----
941 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
942 lastmcmpos = GetMCMReadoutPos(*fPayloadCurr);
945 AliError(MCMError(kPosUnexp));
947 fCurrMcmPos = MCM(*fPayloadCurr);
949 if (EvNo(*fPayloadCurr) != evno) {
951 evno = EvNo(*fPayloadCurr);
953 AliError(MCMError(kPtrgCntMismatch,
954 Form("%i <-> %i", evno, EvNo(*fPayloadCurr)) ));
959 channelcountExp = 21;
962 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
963 Int_t padcoloff = PadColOffset(*fPayloadCurr);
964 Int_t row = Row(*fPayloadCurr);
968 // ----- reading marked ADC channels -----
969 while (channelcount < channelcountExp &&
970 *(fPayloadCurr) != fgkDataEndmarker) {
977 AliDebug(2, Form("Now looking %i words", timebins / 3));
978 Int_t adccol = adccoloff - channelno;
979 Int_t padcol = padcoloff - channelno;
980 while (adcwc < timebins / 3 &&
981 *(fPayloadCurr) != fgkDataEndmarker &&
982 fPayloadCurr - fPayloadStart < fPayloadSize) {
983 int check = 0x3 & *fPayloadCurr;
984 if (channelno % 2 != 0) { // odd channel
985 if (check != 0x2 && channelno < 21) {
986 AliError(MCMError(kAdcCheckInvalid,
987 Form("Invalid check bits: %i for %2i. ADC word in odd channel: %i",
988 check, adcwc+1, channelno) ));
991 else { // even channel
992 if (check != 0x3 && channelno < 21) {
993 AliError(MCMError(kAdcCheckInvalid,
994 Form("Invalid check bits: %i for %2i. ADC word in even channel: %i",
995 check, adcwc+1, channelno) ));
999 // filling the actual timebin data
1000 int tb2 = 0x3ff & *fPayloadCurr >> 22;
1001 int tb1 = 0x3ff & *fPayloadCurr >> 12;
1002 int tb0 = 0x3ff & *fPayloadCurr >> 2;
1003 if (adcwc != 0 || fCurrNtimebins <= 30)
1004 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
1007 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
1008 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
1014 if (adcwc != timebins / 3)
1015 AliError(MCMError(kAdcDataAbort));
1018 if (padcol > 0 && padcol < 144) {
1019 fSignalIndex->AddIndexRC(row, padcol);
1025 if (channelcount != channelcountExp)
1026 AliError(MCMError(kAdcChannelsMiss));
1028 // continue with next MCM
1031 // check for missing MCMs (if header suppression is inactive)
1032 if (mcmcount != mcmcountExp) {
1033 AliError(LinkError(kMissMcmHeaders,
1034 Form("%i not as expected: %i", mcmcount, mcmcountExp) ));
1037 return (fPayloadCurr - start);
1041 Int_t AliTRDrawStream::GetNActiveChannelsFromMask(UInt_t adcmask) const
1043 // return number of active bits in the ADC mask
1045 adcmask = GetActiveChannels(adcmask);
1046 adcmask = adcmask - ((adcmask >> 1) & 0x55555555);
1047 adcmask = (adcmask & 0x33333333) + ((adcmask >> 2) & 0x33333333);
1048 return (((adcmask + (adcmask >> 4)) & 0xF0F0F0F) * 0x1010101) >> 24;
1052 TString AliTRDrawStream::EquipmentError(ErrorCode_t err, TString msg)
1054 // register error according to error code on equipment level
1055 // and return the corresponding error message
1057 fLastError.fSector = fCurrEquipmentId - 1024;
1058 fLastError.fStack = -1;
1059 fLastError.fLink = -1;
1060 fLastError.fRob = -1;
1061 fLastError.fMcm = -1;
1062 fLastError.fError = err;
1065 return (Form("Event %6i: Eq. %2d - %s : ",
1066 fRawReader->GetEventIndex(), fCurrEquipmentId, fgErrorMessages[err]) + msg);
1070 TString AliTRDrawStream::StackError(ErrorCode_t err, TString msg)
1072 // register error according to error code on stack level
1073 // and return the corresponding error message
1075 fLastError.fSector = fCurrEquipmentId - 1024;
1076 fLastError.fStack = fCurrSlot;
1077 fLastError.fLink = -1;
1078 fLastError.fRob = -1;
1079 fLastError.fMcm = -1;
1080 fLastError.fError = err;
1083 return (Form("Event %6i: Eq. %2d S %i - %s : ",
1084 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgErrorMessages[err]) + msg);
1088 TString AliTRDrawStream::LinkError(ErrorCode_t err, TString msg)
1090 // register error according to error code on link level
1091 // and return the corresponding error message
1093 fLastError.fSector = fCurrEquipmentId - 1024;
1094 fLastError.fStack = fCurrSlot;
1095 fLastError.fLink = fCurrLink;
1096 fLastError.fRob = -1;
1097 fLastError.fMcm = -1;
1098 fLastError.fError = err;
1101 return (Form("Event %6i: Eq. %2d S %i l %2i - %s : ",
1102 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgErrorMessages[err]) + msg);
1106 TString AliTRDrawStream::ROBError(ErrorCode_t err, TString msg)
1108 // register error according to error code on ROB level
1109 // and return the corresponding error message
1111 fLastError.fSector = fCurrEquipmentId - 1024;
1112 fLastError.fStack = fCurrSlot;
1113 fLastError.fLink = fCurrLink;
1114 fLastError.fRob = fCurrRobPos;
1115 fLastError.fMcm = -1;
1116 fLastError.fError = err;
1119 return (Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : ",
1120 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgErrorMessages[err]) + msg);
1124 TString AliTRDrawStream::MCMError(ErrorCode_t err, TString msg)
1126 // register error according to error code on MCM level
1127 // and return the corresponding error message
1129 fLastError.fSector = fCurrEquipmentId - 1024;
1130 fLastError.fStack = fCurrSlot;
1131 fLastError.fLink = fCurrLink;
1132 fLastError.fRob = fCurrRobPos;
1133 fLastError.fMcm = fCurrMcmPos;
1134 fLastError.fError = err;
1137 return (Form("Event %6i: Eq. %2d S %i l %2i ROB %i MCM %2i - %s : ",
1138 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fCurrMcmPos, fgErrorMessages[err]) + msg);
1141 const char* AliTRDrawStream::GetErrorMessage(ErrorCode_t errCode)
1143 // return the error message for the given error code
1145 if (errCode > 0 && errCode < kLastErrorCode)
1146 return fgErrorMessages[errCode];