]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDrawStream.cxx
Removal of old raw reader calls (Raphaelle)
[u/mrichter/AliRoot.git] / TRD / AliTRDrawStream.cxx
CommitLineData
0508ca31 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16////////////////////////////////////////////////////////////////////////////
17// //
18// Decoding data from the TRD raw stream //
19// and translation into ADC values //
20// //
21// Author: J. Klein (jochen.klein@cern.ch) //
22// //
23////////////////////////////////////////////////////////////////////////////
24
92223bf6 25#include <cstdio>
26#include <cstdarg>
27
d60fe037 28#include "TClonesArray.h"
29#include "TTree.h"
30
31#include "AliLog.h"
32#include "AliRawReader.h"
33#include "AliTRDdigitsManager.h"
34#include "AliTRDdigitsParam.h"
35#include "AliTRDtrapConfig.h"
36#include "AliTRDarrayADC.h"
37#include "AliTRDarrayDictionary.h"
38#include "AliTRDSignalIndex.h"
39#include "AliTRDtrackletWord.h"
5fdfc9e4 40#include "AliESDTrdTrack.h"
cc26f39c 41#include "AliTreeLoader.h"
d60fe037 42
43#include "AliTRDrawStream.h"
44
45// temporary
46#include "AliRunLoader.h"
47
48ClassImp(AliTRDrawStream)
49
50// some static information
51const Int_t AliTRDrawStream::fgkMcmOrder[] = {12, 13, 14, 15,
cc26f39c 52 8, 9, 10, 11,
53 4, 5, 6, 7,
54 0, 1, 2, 3};
d60fe037 55const Int_t AliTRDrawStream::fgkRobOrder [] = {0, 1, 2, 3};
56const Int_t AliTRDrawStream::fgkNlinks = 12;
57const Int_t AliTRDrawStream::fgkNstacks = 5;
58const UInt_t AliTRDrawStream::fgkDataEndmarker = 0x00000000;
59const UInt_t AliTRDrawStream::fgkTrackletEndmarker = 0x10001000;
60
cc26f39c 61const char* AliTRDrawStream::fgErrorMessages[] = {
d60fe037 62 "Unknown error",
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",
69 "HC header mismatch",
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",
77 "Missing ADC data",
78 "Missing expected ADC channels",
79 "Missing MCM headers"
80};
81
92223bf6 82const Int_t AliTRDrawStream::fgErrorDebugLevel[] = {
83 0,
84 0,
85 2,
86 1,
87 0,
88 1,
89 1,
90 1,
91 1,
92 2,
93 1,
94 1,
95 1,
96 1,
97 2,
98 1,
99 1,
100 1
101};
102
103AliTRDrawStream::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
122};
123
d60fe037 124AliTRDrawStream::AliTRDrawStream(AliRawReader *rawReader) :
cc26f39c 125 fStats(),
f4b3235e 126 fStoreError(&AliTRDrawStream::StoreErrorTree),
d60fe037 127 fRawReader(rawReader),
128 fDigitsManager(0x0),
129 fDigitsParam(0x0),
130 fErrors(0x0),
131 fLastError(),
92223bf6 132 fErrorFlags(0),
d60fe037 133 fPayloadStart(0x0),
134 fPayloadCurr(0x0),
135 fPayloadSize(0),
136 fNtimebins(-1),
137 fLastEvId(-1),
138 fCurrSlot(-1),
139 fCurrLink(-1),
140 fCurrRobPos(-1),
141 fCurrMcmPos(-1),
142 fCurrEquipmentId(0),
143 fCurrSmuIndexHeaderSize(0),
144 fCurrSmuIndexHeaderVersion(0),
145 fCurrTrackEnable(0),
146 fCurrTrackletEnable(0),
147 fCurrStackMask(0),
148 fCurrStackIndexWord(0x0),
149 fCurrStackHeaderSize(0x0),
150 fCurrStackHeaderVersion(0x0),
151 fCurrLinkMask(0x0),
152 fCurrCleanCheckout(0x0),
153 fCurrBoardId(0x0),
154 fCurrHwRev(0x0),
155 fCurrLinkMonitorFlags(0x0),
156 fCurrLinkDataTypeFlags(0x0),
157 fCurrLinkDebugFlags(0x0),
158 fCurrSpecial(-1),
159 fCurrMajor(-1),
160 fCurrMinor(-1),
161 fCurrAddHcWords(-1),
162 fCurrSm(-1),
163 fCurrStack(-1),
164 fCurrLayer(-1),
165 fCurrSide(-1),
166 fCurrHC(-1),
167 fCurrCheck(-1),
168 fCurrNtimebins(-1),
169 fCurrBC(-1),
170 fCurrPtrgCnt(-1),
171 fCurrPtrgPhase(-1),
cc26f39c 172 fNDumpMCMs(0),
d60fe037 173 fTrackletArray(0x0),
174 fAdcArray(0x0),
175 fSignalIndex(0x0),
5fdfc9e4 176 fTrackletTree(0x0),
177 fTracklets(0x0),
178 fTracks(0x0),
179 fMarkers(0x0)
d60fe037 180{
181 // default constructor
182
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];
5fdfc9e4 193 for (Int_t i = 0; i < 100; i++)
194 fDumpMCM[i] = 0;
d60fe037 195
196 // preparing TClonesArray
197 fTrackletArray = new TClonesArray("AliTRDtrackletWord", 256);
198
199 // setting up the error tree
200 fErrors = new TTree("errorStats", "Error statistics");
201 fErrors->SetDirectory(0x0);
5fdfc9e4 202 fErrors->Branch("error", &fLastError);
d60fe037 203 fErrors->SetCircular(1000);
2f9bdd85 204 for (Int_t i = 0; i < 100; i++) {
205 fErrorBuffer[i] = 0;
206 }
207
d60fe037 208}
209
210AliTRDrawStream::~AliTRDrawStream()
211{
212 // destructor
213
214 delete fErrors;
215
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;
226}
227
228Bool_t AliTRDrawStream::ReadEvent(TTree *trackletTree)
229{
230 // read the current event from the raw reader and fill it to the digits manager
231
232 if (!fRawReader) {
233 AliError("No raw reader available");
234 return kFALSE;
235 }
236
237 // tracklet output
67271412 238 ConnectTracklets(trackletTree);
d60fe037 239
240 // some preparations
241 fDigitsParam = 0x0;
242
243 // loop over all DDLs
244 // data starts with GTU payload, i.e. SMU index word
245 UChar_t *buffer = 0x0;
246
247 while (fRawReader->ReadNextData(buffer)) {
248
249 fCurrEquipmentId = fRawReader->GetEquipmentId();
250 AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
251
252 if (fCurrEquipmentId < 1024 || fCurrEquipmentId > 1041) {
92223bf6 253 EquipmentError(kNonTrdEq, "Skipping");
d60fe037 254 continue;
255 }
256
257 // setting the pointer to data and current reading position
258 fPayloadCurr = fPayloadStart = (UInt_t*) (buffer);
259 fPayloadSize = fRawReader->GetDataSize() / sizeof(UInt_t);
cc26f39c 260 fStats.fStatsSector[fCurrEquipmentId - 1024].fBytes = fRawReader->GetDataSize();
d60fe037 261 AliDebug(2, Form("Read buffer of size: %i", fRawReader->GetDataSize()));
262
263 // read SMU index header
264 if (ReadSmHeader() < 0) {
265 AliError(Form("Reading SMU header failed, skipping this DDL %i", fCurrEquipmentId));
266 continue;
267 }
268
269 // read stack index header
270 for (Int_t iStack = 0; iStack < 5; iStack++) {
67271412 271 if ((fCurrStackMask & (1 << iStack)) != 0)
d60fe037 272 ReadStackIndexHeader(iStack);
273 }
274
275 for (Int_t iStack = 0; iStack < 5; iStack++) {
276 fCurrSlot = iStack;
67271412 277 if ((fCurrStackMask & (1 << fCurrSlot)) == 0)
d60fe037 278 continue;
279
280 AliDebug(2, Form("Stack %i, Link mask: 0x%02x", fCurrSlot, fCurrLinkMask[fCurrSlot]));
281 for (Int_t iLink = 0; iLink < 12; iLink++) {
282 fCurrLink = iLink;
283 fCurrHC = fCurrSm * 60 + fCurrSlot * 12 + iLink;
284 if ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) == 0)
285 continue;
286
92223bf6 287 fErrorFlags = 0;
288 // check for link monitor error flag
d60fe037 289 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0)
290 LinkError(kLinkMonitor);
d60fe037 291
292 // read the data from one HC
293 ReadLinkData();
294
295 // read all data endmarkers
92223bf6 296 SeekNextLink();
d60fe037 297 }
298 }
299 }
300 return kTRUE;
301}
302
303
304Bool_t AliTRDrawStream::NextDDL()
305{
0508ca31 306 // continue reading with the next equipment
307
d60fe037 308 if (!fRawReader)
309 return kFALSE;
310
311 fCurrEquipmentId = 0;
312 fCurrSlot = 0;
313 fCurrLink = 0;
314
315 UChar_t *buffer = 0x0;
316
317 while (fRawReader->ReadNextData(buffer)) {
318
319 fCurrEquipmentId = fRawReader->GetEquipmentId();
320 AliDebug(2, Form("equipment: %i", fCurrEquipmentId));
321
322 if (fCurrEquipmentId < 1024 || fCurrEquipmentId > 1041) {
92223bf6 323 EquipmentError(kNonTrdEq, "Skipping");
d60fe037 324 continue;
325 }
326
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()));
331
332 // read SMU index header
333 if (ReadSmHeader() < 0) {
334 AliError(Form("Reading SMU header failed, skipping this DDL %i", fCurrEquipmentId));
335 continue;
336 }
337
338 // read stack index header
339 for (Int_t iStack = 0; iStack < 5; iStack++) {
67271412 340 if ((fCurrStackMask & (1 << iStack)) != 0) {
d60fe037 341 ReadStackIndexHeader(iStack);
342 }
343 }
344 return kTRUE;
345 }
346
347 return kFALSE;
348}
349
350
351Int_t AliTRDrawStream::NextChamber(AliTRDdigitsManager *digMgr, UInt_t ** /* trackletContainer */, UShort_t ** /* errorContainer */)
352{
0508ca31 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
356
d60fe037 357 fDigitsManager = digMgr;
358 fDigitsParam = 0x0;
359
92223bf6 360 fErrorFlags = 0;
361
d60fe037 362 // tracklet output preparation
67271412 363 TTree *trklTree = 0x0;
364 AliRunLoader *rl = AliRunLoader::Instance();
44ac5cbf 365 AliLoader* trdLoader = rl ? rl->GetLoader("TRDLoader") : NULL;
366 AliDataLoader *trklLoader = trdLoader ? trdLoader->GetDataLoader("tracklets") : NULL;
367 if (trklLoader) {
cc26f39c 368 AliTreeLoader *trklTreeLoader = (AliTreeLoader*) trklLoader->GetBaseLoader("tracklets-raw");
369 if (trklTreeLoader)
370 trklTree = trklTreeLoader->Tree();
371 else
372 trklTree = trklLoader->Tree();
67271412 373 }
374
375 if (fTrackletTree != trklTree)
376 ConnectTracklets(trklTree);
d60fe037 377
378 if (!fRawReader) {
379 AliError("No raw reader available");
380 return -1;
381 }
382
383 if (fCurrSlot < 0 || fCurrSlot >= 5) {
384 if (!NextDDL()) {
385 fCurrSlot = -1;
386 return -1;
387 }
388 }
389
390 AliDebug(2, Form("Stack %i, Link %i, mask: 0x%02x", fCurrSlot, fCurrLink, fCurrLinkMask[fCurrSlot]));
391 fCurrHC = (fCurrEquipmentId - 1024) * 60 + fCurrSlot * 12 + fCurrLink;
392
393 if (fCurrLinkMonitorFlags[fCurrSlot*fgkNlinks + fCurrLink] != 0)
394 LinkError(kLinkMonitor);
395
396 // read the data from one HC
397 ReadLinkData();
398
399 // read all data endmarkers
92223bf6 400 SeekNextLink();
d60fe037 401
402 if (fCurrLink % 2 == 0) {
403 // if we just read the A-side HC then also check the B-side
404 fCurrLink++;
1e7466b5 405 fCurrHC++;
d60fe037 406 if (fCurrLinkMask[fCurrSlot] & (1 << fCurrLink)) {
407 ReadLinkData();
92223bf6 408 SeekNextLink();
d60fe037 409 }
410 }
411
412 //??? to check
413 do {
414 fCurrLink++;
415 if (fCurrLink > 11) {
416 fCurrLink = 0;
417 fCurrSlot++;
418 }
419 } while ((fCurrSlot < 5) &&
67271412 420 (((fCurrStackMask & (1 << fCurrSlot)) == 0) ||
421 ((fCurrLinkMask[fCurrSlot] & (1 << fCurrLink))) == 0));
d60fe037 422
25671e69 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);
427 else
428 return (fCurrSm * 30 + fCurrStack * 6 + fCurrLayer);
d60fe037 429}
430
431
432Int_t AliTRDrawStream::ReadSmHeader()
433{
434 // read the SMU index header at the current reading position
435 // and store the information in the corresponding variables
436
437 if (fPayloadCurr - fPayloadStart >= fPayloadSize - 1) {
92223bf6 438 EquipmentError(kUnknown, "SM Header incomplete");
d60fe037 439 return -1;
440 }
441
442 fCurrSmuIndexHeaderSize = ((*fPayloadCurr) >> 16) & 0xffff;
443 fCurrSmuIndexHeaderVersion = ((*fPayloadCurr) >> 12) & 0xf;
5fdfc9e4 444 // fCurrSmuIndexHeaderTrgAvail = ((*fPayloadCurr) >> 9) & 0x1;
445 // fCurrSmuIndexHeaderEvType = ((*fPayloadCurr) >> 7) & 0x3;
d60fe037 446 fCurrTrackEnable = ((*fPayloadCurr) >> 6) & 0x1;
447 fCurrTrackletEnable = ((*fPayloadCurr) >> 5) & 0x1;
448 fCurrStackMask = ((*fPayloadCurr) ) & 0x1f;
449
450 AliDebug(5, Form("SMU header: size: %i, version: %i, track enable: %i, tracklet enable: %i, stack mask: %2x",
451 fCurrSmuIndexHeaderSize,
452 fCurrSmuIndexHeaderVersion,
453 fCurrTrackEnable,
454 fCurrTrackletEnable,
455 fCurrStackMask));
5fdfc9e4 456
457 // decode GTU track words
8d8be166 458 UInt_t trackWord[2] = { 0, 0 };
5fdfc9e4 459 Int_t stack = 0;
460 Int_t idx = 0;
461 for (UInt_t iWord = 4; iWord < fCurrSmuIndexHeaderSize; iWord++) {
462 if (fPayloadCurr[iWord] == 0x10000000) {
463 stack++;
464 idx = 0;
465 }
466 else {
467 if ((idx == 0) &&
468 ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
469 AliDebug(1,Form("stack %i: fast trigger word: 0x%08x", stack, fPayloadCurr[iWord]));
470 continue;
471 }
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]));
475// if (fTracks)
476// new ((*fTracks)[fTracks->GetEntriesFast()]) AliESDTrdTrack(0, 0, trackWord[0], trackWord[1], fCurrEquipmentId-1024);
477 }
478 else {
479 trackWord[idx&0x1] = fPayloadCurr[iWord];
480 }
481 idx++;
482 }
483 }
484
d60fe037 485 fPayloadCurr += fCurrSmuIndexHeaderSize + 1;
486
487 return fCurrSmuIndexHeaderSize + 1;
488}
489
490Int_t AliTRDrawStream::ReadStackIndexHeader(Int_t stack)
491{
492 // read the stack index header
493 // and store the information in the corresponding variables
494
495 fCurrStackIndexWord[stack] = *fPayloadCurr;
496 fCurrStackHeaderSize[stack] = (((*fPayloadCurr) >> 16) & 0xffff) + 1;
497 fCurrStackHeaderVersion[stack] = ((*fPayloadCurr) >> 12) & 0xf;
498 fCurrLinkMask[stack] = (*fPayloadCurr) & 0xfff;
499
67271412 500 if (fPayloadCurr - fPayloadStart >= fPayloadSize - (Int_t) fCurrStackHeaderSize[stack]) {
92223bf6 501 StackError(kStackHeaderInvalid, "Stack index header aborted");
d60fe037 502 return -1;
503 }
504
505 switch (fCurrStackHeaderVersion[stack]) {
506 case 0xa:
507 if (fCurrStackHeaderSize[stack] < 8) {
92223bf6 508 StackError(kStackHeaderInvalid, "Stack header smaller than expected!");
d60fe037 509 return -1;
510 }
511
512 fCurrCleanCheckout[stack] = fPayloadCurr[1] & 0x1;
513 fCurrBoardId[stack] = (fPayloadCurr[1] >> 8) & 0xff;
514 fCurrHwRev[stack] = (fPayloadCurr[1] >> 16) & 0xffff;
515
516 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
517 // A side
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;
521 // B side
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;
525 }
526 break;
527
528 default:
92223bf6 529 StackError(kStackHeaderInvalid, "Invalid Stack Index Header version %x", fCurrStackHeaderVersion[stack]);
d60fe037 530 }
531
532 fPayloadCurr += fCurrStackHeaderSize[stack];
533
534 return fCurrStackHeaderSize[stack];
535}
536
537Int_t AliTRDrawStream::ReadLinkData()
538{
539 // read the data in one link (one HC) until the data endmarker is reached
cc26f39c 540 // returns the number of words read!
d60fe037 541
542 Int_t count = 0;
cc26f39c 543 UInt_t* startPosLink = fPayloadCurr;
544
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]);
549// }
d60fe037 550
5fdfc9e4 551 if (fMarkers)
552 new ((*fMarkers)[fMarkers->GetEntriesFast()])
553 AliTRDrawStreamError(-kHCactive, fCurrSm, fCurrStack, fCurrLink);
554
92223bf6 555 if (fErrorFlags & kDiscardHC)
556 return count;
557
d60fe037 558 count += ReadTracklets();
92223bf6 559 if (fErrorFlags & kDiscardHC)
560 return count;
d60fe037 561
562 count += ReadHcHeader();
92223bf6 563 if (fErrorFlags & kDiscardHC)
564 return count;
d60fe037 565
566 Int_t det = fCurrSm * 30 + fCurrStack * 6 + fCurrLayer;
567
568 if (det > -1 && det < 540) {
569
570 if ((fAdcArray = fDigitsManager->GetDigits(det))) {
571 //fAdcArray->Expand();
572 if (fAdcArray->GetNtime() != fCurrNtimebins)
573 fAdcArray->Allocate(16, 144, fCurrNtimebins);
574 }
575 else {
92223bf6 576 LinkError(kNoDigits);
d60fe037 577 }
578
579 if (!fDigitsParam) {
580 fDigitsParam = fDigitsManager->GetDigitsParam();
581 }
582 if (fDigitsParam) {
583 fDigitsParam->SetPretriggerPhase(det, fCurrPtrgPhase);
584 fDigitsParam->SetNTimeBins(det, fCurrNtimebins);
585 fDigitsParam->SetADCbaseline(det, 10);
586 }
587
588 if (fDigitsManager->UsesDictionaries()) {
589 fDigitsManager->GetDictionary(det, 0)->Reset();
590 fDigitsManager->GetDictionary(det, 1)->Reset();
591 fDigitsManager->GetDictionary(det, 2)->Reset();
592 }
593
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);
601 }
602
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)
610 fPayloadCurr++;
611 count += fPayloadCurr - startPos;
612
613 // feeding TRAP config
614 AliTRDtrapConfig *trapcfg = AliTRDtrapConfig::Instance();
0508ca31 615 trapcfg->ReadPackedConfig(fCurrHC, startPos, fPayloadCurr - startPos);
d60fe037 616 }
617 else {
618 Int_t tpmode = fCurrMajor & 0x7;
619 AliDebug(1, Form("Checking testpattern (mode %i) data", tpmode));
620 ReadTPData(tpmode);
621 }
622 }
623 else if (fCurrMajor & 0x20) {
624 AliDebug(1, "This is a zs event");
625 count += ReadZSData();
626 }
627 else {
628 AliDebug(1, "This is a nozs event");
629 count += ReadNonZSData();
630 }
631 }
632 else {
92223bf6 633 LinkError(kInvalidDetector, "%i", det);
d60fe037 634 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
635 *fPayloadCurr != fgkDataEndmarker)
636 fPayloadCurr++;
d60fe037 637 }
638
cc26f39c 639 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fBytes += (fPayloadCurr - startPosLink) * sizeof(UInt_t);
640 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fBytesRead += count * sizeof(UInt_t);
641 fStats.fStatsSector[fCurrSm].fBytesRead += count * sizeof(UInt_t);
642 fStats.fBytesRead += count * sizeof(UInt_t);
643
d60fe037 644 return count;
645}
646
647Int_t AliTRDrawStream::ReadTracklets()
648{
649 // read the tracklets from one HC
650
651 fTrackletArray->Clear();
652
653 UInt_t *start = fPayloadCurr;
654 while (*(fPayloadCurr) != fgkTrackletEndmarker &&
655 fPayloadCurr - fPayloadStart < fPayloadSize) {
656
5fdfc9e4 657 new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(*(fPayloadCurr), fCurrHC);
d60fe037 658
659 fPayloadCurr++;
660 }
661
662 if (fTrackletArray->GetEntriesFast() > 0) {
663 AliDebug(1, Form("Found %i tracklets in %i %i %i (ev. %i)", fTrackletArray->GetEntriesFast(),
664 fCurrSm, fCurrSlot, fCurrLink, fRawReader->GetEventIndex()));
5fdfc9e4 665 if (fCurrSm > -1 && fCurrSm < 18) {
666 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNTracklets += fTrackletArray->GetEntriesFast();
667 fStats.fStatsSector[fCurrSm].fNTracklets += fTrackletArray->GetEntriesFast();
668 }
d60fe037 669 if (fTrackletTree)
670 fTrackletTree->Fill();
5fdfc9e4 671 if (fTracklets)
672 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
673 new ((*fTracklets)[fTracklets->GetEntriesFast()]) AliTRDtrackletWord(*((AliTRDtrackletWord*)(*fTrackletArray)[iTracklet]));
674 }
d60fe037 675 }
676
677 // loop over remaining tracklet endmarkers
678 while ((*(fPayloadCurr) == fgkTrackletEndmarker &&
679 fPayloadCurr - fPayloadStart < fPayloadSize))
680 fPayloadCurr++;
681
cc26f39c 682 return fPayloadCurr - start;
d60fe037 683}
684
685Int_t AliTRDrawStream::ReadHcHeader()
686{
687 // read and parse the HC header of one HC
688 // and store the information in the corresponding variables
689
690 UInt_t *start = fPayloadCurr;
691 // check not to be at the data endmarker
692 if (*fPayloadCurr == fgkDataEndmarker)
693 return 0;
694
695 fCurrSpecial = (*fPayloadCurr >> 31) & 0x1;
696 fCurrMajor = (*fPayloadCurr >> 24) & 0x7f;
697 fCurrMinor = (*fPayloadCurr >> 17) & 0x7f;
698 fCurrAddHcWords = (*fPayloadCurr >> 14) & 0x7;
699 fCurrSm = (*fPayloadCurr >> 9) & 0x1f;
700 fCurrLayer = (*fPayloadCurr >> 6) & 0x7;
701 fCurrStack = (*fPayloadCurr >> 3) & 0x7;
702 fCurrSide = (*fPayloadCurr >> 2) & 0x1;
703 fCurrCheck = (*fPayloadCurr) & 0x3;
704
0508ca31 705 if (fCurrSm != (((Int_t) fCurrEquipmentId) - 1024) ||
d60fe037 706 fCurrStack != fCurrSlot ||
707 fCurrLayer != fCurrLink / 2 ||
708 fCurrSide != fCurrLink % 2) {
92223bf6 709 LinkError(kHCmismatch,
710 "HC: %i, %i, %i, %i\n 0x%08x 0x%08x 0x%08x 0x%08x",
711 fCurrSm, fCurrStack, fCurrLayer, fCurrSide,
712 fPayloadCurr[0], fPayloadCurr[1], fPayloadCurr[2], fPayloadCurr[3]);;
d60fe037 713 }
714 if (fCurrCheck != 0x1) {
92223bf6 715 LinkError(kHCcheckFailed);
d60fe037 716 }
717
718 if (fCurrAddHcWords > 0) {
719 fCurrNtimebins = (fPayloadCurr[1] >> 26) & 0x3f;
720 fCurrBC = (fPayloadCurr[1] >> 10) & 0xffff;
721 fCurrPtrgCnt = (fPayloadCurr[1] >> 6) & 0xf;
722 fCurrPtrgPhase = (fPayloadCurr[1] >> 2) & 0xf;
723 }
724
725 fPayloadCurr += 1 + fCurrAddHcWords;
726
5fdfc9e4 727 return (fPayloadCurr - start);
d60fe037 728}
729
730Int_t AliTRDrawStream::ReadTPData(Int_t mode)
731{
732 // testing of testpattern 1 to 3 (hardcoded), 0 missing
733 // evcnt checking missing
734 Int_t cpu = 0;
735 Int_t cpufromchannel[] = {0, 0, 0, 0, 0, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 3, 3, 3, 3, 3, 3};
736 Int_t evcnt = 0;
737 Int_t count = 0;
738 Int_t mcmcount = -1;
739 Int_t wordcount = 0;
740 Int_t channelcount = 0;
741 UInt_t expword = 0;
742 UInt_t expadcval = 0;
743 UInt_t diff = 0;
744 Int_t lastmcmpos = -1;
745 Int_t lastrobpos = -1;
746
747 UInt_t* start = fPayloadCurr;
748
749 while (*(fPayloadCurr) != fgkDataEndmarker &&
750 fPayloadCurr - fPayloadStart < fPayloadSize - 1) {
751
752 // ----- Checking MCM Header -----
753 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
754 mcmcount++;
755
756 // ----- checking for proper readout order - ROB -----
757 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
758 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
759 }
760 else {
92223bf6 761 ROBError(kPosUnexp);
d60fe037 762 }
763 fCurrRobPos = ROB(*fPayloadCurr);
764
765 // ----- checking for proper readout order - MCM -----
766 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
767 lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
768 }
769 else {
92223bf6 770 MCMError(kPosUnexp);
d60fe037 771 }
772 fCurrMcmPos = MCM(*fPayloadCurr);
773
774
775 fPayloadCurr++;
776
777 evcnt = 0x3f & *fPayloadCurr >> 26;
778 cpu = -1;
779 channelcount = 0;
780 while (channelcount < 21) {
781 count = 0;
782 if (cpu != cpufromchannel[channelcount]) {
783 cpu = cpufromchannel[channelcount];
784 expadcval = (1 << 9) | (fCurrRobPos << 6) | (fCurrMcmPos << 2) | cpu;
785 wordcount = 0;
786 }
787
788 while (count < 10) {
789 if (channelcount % 2 == 0)
790 expword = 0x3;
791 else
792 expword = 0x2;
793
794 if (mode == 1) {
795 // ----- TP 1 -----
796 expword |= expadcval << 2;
797 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
798 expword |= expadcval << 12;
799 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
800 expword |= expadcval << 22;
801 expadcval = ( (expadcval << 1) | ( ( (expadcval >> 9) ^ (expadcval >> 6) ) & 1) ) & 0x3FF;
802 }
803 else if (mode == 2) {
804 // ----- TP 2 ------
805 expword = ((0x3f & evcnt) << 26) | ((fCurrSm + 1) << 21) | ((fCurrLayer + 1) << 18) |
806 ((fCurrStack + 1) << 15) |
807 (fCurrRobPos << 12) | (fCurrMcmPos << 8) | (cpu << 6) | (wordcount + 1);
808 }
809 else if (mode == 3) {
810 // ----- TP 3 -----
811 expword = ((0xfff & evcnt) << 20) | (fCurrSm << 15) | (fCurrLink/2 << 12) | (fCurrStack << 9) |
812 (fCurrRobPos << 6) | (fCurrMcmPos << 2) | (cpu << 0);
813 }
814 else {
815 expword = 0;
92223bf6 816 LinkError(kTPmodeInvalid, "Just reading");
d60fe037 817 }
818
819 diff = *fPayloadCurr ^ expword;
820 if (diff != 0) {
92223bf6 821 MCMError(kTPmismatch,
822 "Seen 0x%08x, expected 0x%08x, diff: 0x%08x (0x%02x)",
823 *fPayloadCurr, expword, diff, 0xff & (diff | diff >> 8 | diff >> 16 | diff >> 24));;
d60fe037 824 }
825 fPayloadCurr++;
826 count++;
827 wordcount++;
828 }
829 channelcount++;
830 }
831 // continue with next MCM
832 }
833 return fPayloadCurr - start;
834}
835
836
837Int_t AliTRDrawStream::ReadZSData()
838{
839 // read the zs data from one link from the current reading position
840
841 UInt_t *start = fPayloadCurr;
842
843 Int_t mcmcount = 0;
0508ca31 844 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
d60fe037 845 Int_t channelcount = 0;
0508ca31 846 Int_t channelcountExp = 0;
847 Int_t channelcountMax = 0;
d60fe037 848 Int_t timebins;
849 Int_t currentTimebin = 0;
850 Int_t adcwc = 0;
851 Int_t evno = -1;
852 Int_t lastmcmpos = -1;
853 Int_t lastrobpos = -1;
854
855 if (fCurrNtimebins != fNtimebins) {
856 if (fNtimebins > 0)
92223bf6 857 LinkError(kNtimebinsChanged,
858 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
d60fe037 859 fNtimebins = fCurrNtimebins;
860 }
861
862 timebins = fNtimebins;
863
864 while (*(fPayloadCurr) != fgkDataEndmarker &&
865 fPayloadCurr - fPayloadStart < fPayloadSize) {
866
867 // ----- Checking MCM Header -----
868 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
cc26f39c 869 UInt_t *startPosMCM = fPayloadCurr;
d60fe037 870
871 // ----- checking for proper readout order - ROB -----
872 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
873 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
874 }
875 else {
92223bf6 876 ROBError(kPosUnexp);
d60fe037 877 }
878 fCurrRobPos = ROB(*fPayloadCurr);
879
880 // ----- checking for proper readout order - MCM -----
881 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) > lastmcmpos) {
882 lastmcmpos = GetMCMReadoutPos(lastmcmpos);
883 }
884 else {
92223bf6 885 MCMError(kPosUnexp);
d60fe037 886 }
887 fCurrMcmPos = MCM(*fPayloadCurr);
888
889 if (EvNo(*fPayloadCurr) != evno) {
890 if (evno == -1)
891 evno = EvNo(*fPayloadCurr);
892 else {
92223bf6 893 MCMError(kPtrgCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
d60fe037 894 }
895 }
896 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
897 Int_t padcoloff = PadColOffset(*fPayloadCurr);
898 Int_t row = Row(*fPayloadCurr);
899 fPayloadCurr++;
900
901 // ----- Reading ADC channels -----
902 AliDebug(2, Form("ADC mask: 0x%08x", *fPayloadCurr));
903
904 // ----- analysing the ADC mask -----
905 channelcount = 0;
0508ca31 906 channelcountExp = GetNActiveChannelsFromMask(*fPayloadCurr);
907 channelcountMax = GetNActiveChannels(*fPayloadCurr);
d60fe037 908 Int_t channelmask = GetActiveChannels(*fPayloadCurr);
909 Int_t channelno = -1;
910 fPayloadCurr++;
911
0508ca31 912 if (channelcountExp != channelcountMax) {
913 if (channelcountExp > channelcountMax) {
914 Int_t temp = channelcountExp;
915 channelcountExp = channelcountMax;
916 channelcountMax = temp;
d60fe037 917 }
0508ca31 918 while (channelcountExp < channelcountMax && channelcountExp < 21 &&
919 fPayloadCurr - fPayloadStart < fPayloadSize - 10 * channelcountExp - 1) {
92223bf6 920 MCMError(kAdcMaskInconsistent,
921 "Possible MCM-H: 0x%08x, possible ADC-mask: 0x%08x",
922 *(fPayloadCurr + 10 * channelcountExp),
923 *(fPayloadCurr + 10 * channelcountExp + 1) );
0508ca31 924 if (!CouldBeMCMhdr( *(fPayloadCurr + 10 * channelcountExp)) && !CouldBeADCmask( *(fPayloadCurr + 10 * channelcountExp + 1)))
925 channelcountExp++;
d60fe037 926 else {
927 break;
928 }
929 }
92223bf6 930 MCMError(kAdcMaskInconsistent,
931 "Inconsistency in no. of active channels: Counter: %i, Mask: %i, chosen: %i!",
932 GetNActiveChannels(fPayloadCurr[-1]), GetNActiveChannelsFromMask(fPayloadCurr[-1]), channelcountExp);
d60fe037 933 }
0508ca31 934 AliDebug(2, Form("expecting %i active channels, timebins: %i", channelcountExp, fCurrNtimebins));
d60fe037 935
936 // ----- reading marked ADC channels -----
0508ca31 937 while (channelcount < channelcountExp && *(fPayloadCurr) != fgkDataEndmarker) {
e29e514c 938 if (channelno < 20)
d60fe037 939 channelno++;
e29e514c 940 while (channelno < 20 && (channelmask & 1 << channelno) == 0)
d60fe037 941 channelno++;
942
943 if (fCurrNtimebins > 30) {
944 currentTimebin = ((*fPayloadCurr >> 2) & 0x3f);
945 timebins = ((*fPayloadCurr >> 8) & 0xf) * 3;
946 }
947 else {
948 currentTimebin = 0;
949 }
950
951 adcwc = 0;
952 AliDebug(2, Form("Now looking %i words", timebins / 3));
953 Int_t adccol = adccoloff - channelno;
954 Int_t padcol = padcoloff - channelno;
cc26f39c 955// if (adccol < 3 || adccol > 165)
956// AliInfo(Form("writing channel %i of det %3i %i:%2i to adcrow/-col: %i/%i padcol: %i",
957// channelno, fCurrHC/2, fCurrRobPos, fCurrMcmPos, row, adccol, padcol));
958
d60fe037 959 while (adcwc < timebins / 3 &&
960 *(fPayloadCurr) != fgkDataEndmarker &&
961 fPayloadCurr - fPayloadStart < fPayloadSize) {
962 int check = 0x3 & *fPayloadCurr;
963 if (channelno % 2 != 0) { // odd channel
964 if (check != 0x2 && channelno < 21) {
92223bf6 965 MCMError(kAdcCheckInvalid,
966 "%i for %2i. ADC word in odd channel %i",
967 check, adcwc+1, channelno);
d60fe037 968 }
969 }
970 else { // even channel
971 if (check != 0x3 && channelno < 21) {
92223bf6 972 MCMError(kAdcCheckInvalid,
973 "%i for %2i. ADC word in even channel %i",
974 check, adcwc+1, channelno);
d60fe037 975 }
976 }
977
978 // filling the actual timebin data
979 int tb2 = 0x3ff & *fPayloadCurr >> 22;
980 int tb1 = 0x3ff & *fPayloadCurr >> 12;
981 int tb0 = 0x3ff & *fPayloadCurr >> 2;
982 if (adcwc != 0 || fCurrNtimebins <= 30)
983 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
984 else
985 tb0 = -1;
986 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
987 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
988
989 adcwc++;
990 fPayloadCurr++;
991 }
992
993 if (adcwc != timebins / 3)
92223bf6 994 MCMError(kAdcDataAbort);
d60fe037 995
996 // adding index
997 if (padcol > 0 && padcol < 144) {
998 fSignalIndex->AddIndexRC(row, padcol);
999 }
1000
1001 channelcount++;
1002 }
cc26f39c 1003
1004 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNChannels += channelcount;
1005 fStats.fStatsSector[fCurrSm].fNChannels += channelcount;
0508ca31 1006 if (channelcount != channelcountExp)
92223bf6 1007 MCMError(kAdcChannelsMiss);
d60fe037 1008
1009 mcmcount++;
cc26f39c 1010 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNMCMs++;
1011 fStats.fStatsSector[fCurrSm].fNMCMs++;
1012
1013 if (IsDumping() && DumpingMCM(fCurrHC/2, fCurrRobPos, fCurrMcmPos)) {
92223bf6 1014 DumpRaw(Form("Event %i: Det %3i ROB %i MCM %2i", fRawReader->GetEventIndex(), fCurrHC/2, fCurrRobPos, fCurrMcmPos),
cc26f39c 1015 startPosMCM, fPayloadCurr - startPosMCM);
1016 }
1017
d60fe037 1018 // continue with next MCM
1019 }
1020
1021 // check for missing MCMs (if header suppression is inactive)
67271412 1022 if (((fCurrMajor & 0x1) == 0) && (mcmcount != mcmcountExp)) {
92223bf6 1023 LinkError(kMissMcmHeaders,
1024 "No. of MCM headers %i not as expected: %i",
1025 mcmcount, mcmcountExp);
d60fe037 1026 }
1027
1028 return (fPayloadCurr - start);
1029}
1030
1031Int_t AliTRDrawStream::ReadNonZSData()
1032{
1033 // read the non-zs data from one link from the current reading position
1034
1035 UInt_t *start = fPayloadCurr;
1036
1037 Int_t mcmcount = 0;
0508ca31 1038 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
d60fe037 1039 Int_t channelcount = 0;
0508ca31 1040 Int_t channelcountExp = 0;
d60fe037 1041 Int_t timebins;
1042 Int_t currentTimebin = 0;
1043 Int_t adcwc = 0;
1044 Int_t evno = -1;
1045 Int_t lastmcmpos = -1;
1046 Int_t lastrobpos = -1;
1047
1048 if (fCurrNtimebins != fNtimebins) {
1049 if (fNtimebins > 0)
92223bf6 1050 LinkError(kNtimebinsChanged,
1051 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
d60fe037 1052 fNtimebins = fCurrNtimebins;
1053 }
1054
1055 timebins = fNtimebins;
1056
1057 while (*(fPayloadCurr) != fgkDataEndmarker &&
1058 fPayloadCurr - fPayloadStart < fPayloadSize - 2) {
1059
1060 // ----- Checking MCM Header -----
1061 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
1062
1063 // ----- checking for proper readout order - ROB -----
1064 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1065 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1066 }
1067 else {
92223bf6 1068 ROBError(kPosUnexp);
d60fe037 1069 }
1070 fCurrRobPos = ROB(*fPayloadCurr);
1071
1072 // ----- checking for proper readout order - MCM -----
1073 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
1074 lastmcmpos = GetMCMReadoutPos(*fPayloadCurr);
1075 }
1076 else {
92223bf6 1077 MCMError(kPosUnexp);
d60fe037 1078 }
1079 fCurrMcmPos = MCM(*fPayloadCurr);
1080
1081 if (EvNo(*fPayloadCurr) != evno) {
1082 if (evno == -1)
1083 evno = EvNo(*fPayloadCurr);
1084 else {
92223bf6 1085 MCMError(kPtrgCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
d60fe037 1086 }
1087 }
1088
1089 channelcount = 0;
0508ca31 1090 channelcountExp = 21;
d60fe037 1091 int channelno = -1;
1092
1093 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
1094 Int_t padcoloff = PadColOffset(*fPayloadCurr);
1095 Int_t row = Row(*fPayloadCurr);
1096
1097 fPayloadCurr++;
1098
1099 // ----- reading marked ADC channels -----
0508ca31 1100 while (channelcount < channelcountExp &&
d60fe037 1101 *(fPayloadCurr) != fgkDataEndmarker) {
e29e514c 1102 if (channelno < 20)
d60fe037 1103 channelno++;
1104
1105 currentTimebin = 0;
1106
1107 adcwc = 0;
1108 AliDebug(2, Form("Now looking %i words", timebins / 3));
1109 Int_t adccol = adccoloff - channelno;
1110 Int_t padcol = padcoloff - channelno;
1111 while (adcwc < timebins / 3 &&
1112 *(fPayloadCurr) != fgkDataEndmarker &&
1113 fPayloadCurr - fPayloadStart < fPayloadSize) {
1114 int check = 0x3 & *fPayloadCurr;
1115 if (channelno % 2 != 0) { // odd channel
1116 if (check != 0x2 && channelno < 21) {
92223bf6 1117 MCMError(kAdcCheckInvalid,
1118 "%i for %2i. ADC word in odd channel %i",
1119 check, adcwc+1, channelno);
d60fe037 1120 }
1121 }
1122 else { // even channel
1123 if (check != 0x3 && channelno < 21) {
92223bf6 1124 MCMError(kAdcCheckInvalid,
1125 "%i for %2i. ADC word in even channel %i",
1126 check, adcwc+1, channelno);
d60fe037 1127 }
1128 }
1129
1130 // filling the actual timebin data
1131 int tb2 = 0x3ff & *fPayloadCurr >> 22;
1132 int tb1 = 0x3ff & *fPayloadCurr >> 12;
1133 int tb0 = 0x3ff & *fPayloadCurr >> 2;
1134 if (adcwc != 0 || fCurrNtimebins <= 30)
1135 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb0);
1136 else
1137 tb0 = -1;
1138 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
1139 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
1140
1141 adcwc++;
1142 fPayloadCurr++;
1143 }
1144
1145 if (adcwc != timebins / 3)
92223bf6 1146 MCMError(kAdcDataAbort);
d60fe037 1147
1148 // adding index
1149 if (padcol > 0 && padcol < 144) {
1150 fSignalIndex->AddIndexRC(row, padcol);
1151 }
1152
1153 channelcount++;
1154 }
1155
0508ca31 1156 if (channelcount != channelcountExp)
92223bf6 1157 MCMError(kAdcChannelsMiss);
d60fe037 1158 mcmcount++;
1159 // continue with next MCM
1160 }
1161
1162 // check for missing MCMs (if header suppression is inactive)
0508ca31 1163 if (mcmcount != mcmcountExp) {
92223bf6 1164 LinkError(kMissMcmHeaders,
1165 "%i not as expected: %i", mcmcount, mcmcountExp);
d60fe037 1166 }
1167
1168 return (fPayloadCurr - start);
1169}
1170
92223bf6 1171Int_t AliTRDrawStream::SeekNextLink()
1172{
1173 UInt_t *start = fPayloadCurr;
1174
1175 // read until data endmarkers
1176 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1177 *fPayloadCurr != fgkDataEndmarker)
1178 fPayloadCurr++;
1179
1180 // read all data endmarkers
1181 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1182 *fPayloadCurr == fgkDataEndmarker)
1183 fPayloadCurr++;
1184
1185 return (fPayloadCurr - start);
1186}
1187
67271412 1188Bool_t AliTRDrawStream::ConnectTracklets(TTree *trklTree)
1189{
1190 fTrackletTree = trklTree;
1191 if (!fTrackletTree)
1192 return kTRUE;
1193
cc26f39c 1194 if (!fTrackletTree->GetBranch("hc"))
67271412 1195 fTrackletTree->Branch("hc", &fCurrHC, "hc/I");
cc26f39c 1196 else
67271412 1197 fTrackletTree->SetBranchAddress("hc", &fCurrHC);
cc26f39c 1198
1199 if (!fTrackletTree->GetBranch("trkl"))
1200 fTrackletTree->Branch("trkl", &fTrackletArray);
1201 else
67271412 1202 fTrackletTree->SetBranchAddress("trkl", &fTrackletArray);
cc26f39c 1203
67271412 1204 return kTRUE;
1205}
1206
1207
1d62be37 1208void AliTRDrawStream::EquipmentError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1209{
0508ca31 1210 // register error according to error code on equipment level
1211 // and return the corresponding error message
1212
d60fe037 1213 fLastError.fSector = fCurrEquipmentId - 1024;
1214 fLastError.fStack = -1;
1215 fLastError.fLink = -1;
1216 fLastError.fRob = -1;
1217 fLastError.fMcm = -1;
1218 fLastError.fError = err;
f4b3235e 1219 (this->*fStoreError)();
d60fe037 1220
92223bf6 1221 va_list ap;
1222 if (fgErrorDebugLevel[err] > 10)
1223 AliDebug(fgErrorDebugLevel[err],
1224 Form("Event %6i: Eq. %2d - %s : %s",
1225 fRawReader->GetEventIndex(), fCurrEquipmentId, fgErrorMessages[err],
1d62be37 1226 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1227 else
1228 AliError(Form("Event %6i: Eq. %2d - %s : %s",
1229 fRawReader->GetEventIndex(), fCurrEquipmentId, fgErrorMessages[err],
1d62be37 1230 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1231 fErrorFlags |= fgErrorBehav[err];
d60fe037 1232}
1233
1234
1d62be37 1235void AliTRDrawStream::StackError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1236{
0508ca31 1237 // register error according to error code on stack level
1238 // and return the corresponding error message
1239
d60fe037 1240 fLastError.fSector = fCurrEquipmentId - 1024;
1241 fLastError.fStack = fCurrSlot;
1242 fLastError.fLink = -1;
1243 fLastError.fRob = -1;
1244 fLastError.fMcm = -1;
1245 fLastError.fError = err;
f4b3235e 1246 (this->*fStoreError)();
d60fe037 1247
92223bf6 1248 va_list ap;
e29e514c 1249 if (fgErrorDebugLevel[err] > 0)
92223bf6 1250 AliDebug(fgErrorDebugLevel[err],
1251 Form("Event %6i: Eq. %2d S %i - %s : %s",
1252 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgErrorMessages[err],
1d62be37 1253 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1254 else
1255 AliError(Form("Event %6i: Eq. %2d S %i - %s : %s",
1256 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgErrorMessages[err],
1d62be37 1257 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1258 fErrorFlags |= fgErrorBehav[err];
d60fe037 1259}
1260
1261
1d62be37 1262void AliTRDrawStream::LinkError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1263{
0508ca31 1264 // register error according to error code on link level
1265 // and return the corresponding error message
1266
d60fe037 1267 fLastError.fSector = fCurrEquipmentId - 1024;
1268 fLastError.fStack = fCurrSlot;
1269 fLastError.fLink = fCurrLink;
1270 fLastError.fRob = -1;
1271 fLastError.fMcm = -1;
1272 fLastError.fError = err;
f4b3235e 1273 (this->*fStoreError)();
d60fe037 1274
92223bf6 1275 va_list ap;
e29e514c 1276 if (fgErrorDebugLevel[err] > 0)
92223bf6 1277 AliDebug(fgErrorDebugLevel[err],
1278 Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
1279 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgErrorMessages[err],
1d62be37 1280 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1281 else
1282 AliError(Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
1283 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgErrorMessages[err],
1d62be37 1284 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1285 fErrorFlags |= fgErrorBehav[err];
d60fe037 1286}
1287
1288
1d62be37 1289void AliTRDrawStream::ROBError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1290{
0508ca31 1291 // register error according to error code on ROB level
1292 // and return the corresponding error message
1293
d60fe037 1294 fLastError.fSector = fCurrEquipmentId - 1024;
1295 fLastError.fStack = fCurrSlot;
1296 fLastError.fLink = fCurrLink;
1297 fLastError.fRob = fCurrRobPos;
1298 fLastError.fMcm = -1;
1299 fLastError.fError = err;
f4b3235e 1300 (this->*fStoreError)();
d60fe037 1301
92223bf6 1302 va_list ap;
e29e514c 1303 if (fgErrorDebugLevel[err] > 0)
92223bf6 1304 AliDebug(fgErrorDebugLevel[err],
1305 Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : %s",
1306 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgErrorMessages[err],
1d62be37 1307 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1308 else
1309 AliError(Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : %s",
1310 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgErrorMessages[err],
1d62be37 1311 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1312 fErrorFlags |= fgErrorBehav[err];
d60fe037 1313}
1314
1315
1d62be37 1316void AliTRDrawStream::MCMError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1317{
0508ca31 1318 // register error according to error code on MCM level
1319 // and return the corresponding error message
1320
d60fe037 1321 fLastError.fSector = fCurrEquipmentId - 1024;
1322 fLastError.fStack = fCurrSlot;
1323 fLastError.fLink = fCurrLink;
1324 fLastError.fRob = fCurrRobPos;
1325 fLastError.fMcm = fCurrMcmPos;
1326 fLastError.fError = err;
f4b3235e 1327 (this->*fStoreError)();
d60fe037 1328
92223bf6 1329 va_list ap;
e29e514c 1330 if (fgErrorDebugLevel[err] > 0)
92223bf6 1331 AliDebug(fgErrorDebugLevel[err],
1332 Form("Event %6i: Eq. %2d S %i l %2i ROB %i MCM %2i - %s : %s",
1333 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fCurrMcmPos, fgErrorMessages[err],
1d62be37 1334 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1335 else
1336 AliError(Form("Event %6i: Eq. %2d S %i l %2i ROB %i MCM %2i - %s : %s",
1337 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fCurrMcmPos, fgErrorMessages[err],
1d62be37 1338 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1339 fErrorFlags |= fgErrorBehav[err];
d60fe037 1340}
1341
1342const char* AliTRDrawStream::GetErrorMessage(ErrorCode_t errCode)
1343{
0508ca31 1344 // return the error message for the given error code
1345
d60fe037 1346 if (errCode > 0 && errCode < kLastErrorCode)
1347 return fgErrorMessages[errCode];
1348 else
1349 return "";
1350}
cc26f39c 1351
1352void AliTRDrawStream::AliTRDrawStats::ClearStats()
1353{
1354 // clear statistics (includes clearing sector-wise statistics)
1355
1356 fBytesRead = 0;
1357 for (Int_t iSector = 0; iSector < 18; iSector++) {
1358 fStatsSector[iSector].ClearStats();
1359 }
1360
1361}
1362
1363void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::ClearStats()
1364{
1365 // clear statistics (includes clearing HC-wise statistics)
1366
1367 fBytes = 0;
1368 fBytesRead = 0;
1369 fNTracklets = 0;
1370 fNMCMs = 0;
1371 fNChannels = 0;
1372
1373 for (Int_t iHC = 0; iHC < 60; iHC++) {
1374 fStatsHC[iHC].ClearStats();
1375 }
1376}
1377
1378void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::AliTRDrawStatsHC::ClearStats()
1379{
1380 // clear statistics
1381
1382 fBytes = 0;
1383 fBytesRead = 0;
1384 fNTracklets = 0;
1385 fNMCMs = 0;
1386 fNChannels = 0;
1387}
1388
1389void AliTRDrawStream::SetDumpMCM(Int_t det, Int_t rob, Int_t mcm, Bool_t dump)
1390{
1391 // mark MCM for dumping of raw data
1392
1393 if (dump) {
1394 fDumpMCM[fNDumpMCMs++] = (det << 7) | (rob << 4) | mcm;
1395 }
1396 else {
1397 Int_t iMCM;
1398 for (iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
1399 if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
1400 fNDumpMCMs--;
1401 break;
1402 }
1403 }
1404 for ( ; iMCM < fNDumpMCMs; iMCM++) {
1405 fDumpMCM[iMCM] = fDumpMCM[iMCM+1];
1406 }
1407 }
1408}
1409
1410Bool_t AliTRDrawStream::DumpingMCM(Int_t det, Int_t rob, Int_t mcm)
1411{
1412 // check if MCM data should be dumped
1413
1414 for (Int_t iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
1415 if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
1416 return kTRUE;
1417 }
1418 }
1419 return kFALSE;
1420}
1421
1422void AliTRDrawStream::DumpRaw(TString title, UInt_t *start, Int_t length)
1423{
1424 // dump raw data
1425
1426 title += "\n";
1427 Int_t pos = 0;
1428 for ( ; pos+3 < length; pos += 4) {
1429 title += Form("0x%08x 0x%08x 0x%08x 0x%08x\n",
1430 start[pos+0], start[pos+1], start[pos+2], start[pos+3]);
1431 }
1432 for ( ; pos < length; pos++) {
1433 title += Form("0x%08x ", start[pos]);
1434 }
1435 AliInfo(title);
1436}
f4b3235e 1437
1438AliTRDrawStream::AliTRDrawStreamError::AliTRDrawStreamError(Int_t error, Int_t sector, Int_t stack, Int_t link, Int_t rob, Int_t mcm) :
1439 fError(error),
1440 fSector(sector),
1441 fStack(stack),
1442 fLink(link),
1443 fRob(rob),
1444 fMcm(mcm)
1445{
1446 // ctor
1447
1448}