]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDrawStream.cxx
leak fix
[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 }
637666cd 638
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);
644 }
cc26f39c 645
d60fe037 646 return count;
647}
648
649Int_t AliTRDrawStream::ReadTracklets()
650{
651 // read the tracklets from one HC
652
653 fTrackletArray->Clear();
654
655 UInt_t *start = fPayloadCurr;
656 while (*(fPayloadCurr) != fgkTrackletEndmarker &&
657 fPayloadCurr - fPayloadStart < fPayloadSize) {
658
5fdfc9e4 659 new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(*(fPayloadCurr), fCurrHC);
d60fe037 660
661 fPayloadCurr++;
662 }
663
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()));
5fdfc9e4 667 if (fCurrSm > -1 && fCurrSm < 18) {
668 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNTracklets += fTrackletArray->GetEntriesFast();
669 fStats.fStatsSector[fCurrSm].fNTracklets += fTrackletArray->GetEntriesFast();
670 }
d60fe037 671 if (fTrackletTree)
672 fTrackletTree->Fill();
5fdfc9e4 673 if (fTracklets)
674 for (Int_t iTracklet = 0; iTracklet < fTrackletArray->GetEntriesFast(); iTracklet++) {
675 new ((*fTracklets)[fTracklets->GetEntriesFast()]) AliTRDtrackletWord(*((AliTRDtrackletWord*)(*fTrackletArray)[iTracklet]));
676 }
d60fe037 677 }
678
679 // loop over remaining tracklet endmarkers
680 while ((*(fPayloadCurr) == fgkTrackletEndmarker &&
681 fPayloadCurr - fPayloadStart < fPayloadSize))
682 fPayloadCurr++;
683
cc26f39c 684 return fPayloadCurr - start;
d60fe037 685}
686
687Int_t AliTRDrawStream::ReadHcHeader()
688{
689 // read and parse the HC header of one HC
690 // and store the information in the corresponding variables
691
692 UInt_t *start = fPayloadCurr;
693 // check not to be at the data endmarker
694 if (*fPayloadCurr == fgkDataEndmarker)
695 return 0;
696
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;
706
0508ca31 707 if (fCurrSm != (((Int_t) fCurrEquipmentId) - 1024) ||
d60fe037 708 fCurrStack != fCurrSlot ||
709 fCurrLayer != fCurrLink / 2 ||
710 fCurrSide != fCurrLink % 2) {
92223bf6 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]);;
d60fe037 715 }
716 if (fCurrCheck != 0x1) {
92223bf6 717 LinkError(kHCcheckFailed);
d60fe037 718 }
719
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;
725 }
726
727 fPayloadCurr += 1 + fCurrAddHcWords;
728
5fdfc9e4 729 return (fPayloadCurr - start);
d60fe037 730}
731
732Int_t AliTRDrawStream::ReadTPData(Int_t mode)
733{
734 // testing of testpattern 1 to 3 (hardcoded), 0 missing
735 // evcnt checking missing
736 Int_t cpu = 0;
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};
738 Int_t evcnt = 0;
739 Int_t count = 0;
740 Int_t mcmcount = -1;
741 Int_t wordcount = 0;
742 Int_t channelcount = 0;
743 UInt_t expword = 0;
744 UInt_t expadcval = 0;
745 UInt_t diff = 0;
746 Int_t lastmcmpos = -1;
747 Int_t lastrobpos = -1;
748
749 UInt_t* start = fPayloadCurr;
750
751 while (*(fPayloadCurr) != fgkDataEndmarker &&
752 fPayloadCurr - fPayloadStart < fPayloadSize - 1) {
753
754 // ----- Checking MCM Header -----
755 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
756 mcmcount++;
757
758 // ----- checking for proper readout order - ROB -----
759 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
760 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
761 }
762 else {
92223bf6 763 ROBError(kPosUnexp);
d60fe037 764 }
765 fCurrRobPos = ROB(*fPayloadCurr);
766
767 // ----- checking for proper readout order - MCM -----
768 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
769 lastmcmpos = GetMCMReadoutPos(MCM(*fPayloadCurr));
770 }
771 else {
92223bf6 772 MCMError(kPosUnexp);
d60fe037 773 }
774 fCurrMcmPos = MCM(*fPayloadCurr);
775
776
777 fPayloadCurr++;
778
779 evcnt = 0x3f & *fPayloadCurr >> 26;
780 cpu = -1;
781 channelcount = 0;
782 while (channelcount < 21) {
783 count = 0;
784 if (cpu != cpufromchannel[channelcount]) {
785 cpu = cpufromchannel[channelcount];
786 expadcval = (1 << 9) | (fCurrRobPos << 6) | (fCurrMcmPos << 2) | cpu;
787 wordcount = 0;
788 }
789
790 while (count < 10) {
791 if (channelcount % 2 == 0)
792 expword = 0x3;
793 else
794 expword = 0x2;
795
796 if (mode == 1) {
797 // ----- TP 1 -----
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;
804 }
805 else if (mode == 2) {
806 // ----- TP 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);
810 }
811 else if (mode == 3) {
812 // ----- TP 3 -----
813 expword = ((0xfff & evcnt) << 20) | (fCurrSm << 15) | (fCurrLink/2 << 12) | (fCurrStack << 9) |
814 (fCurrRobPos << 6) | (fCurrMcmPos << 2) | (cpu << 0);
815 }
816 else {
817 expword = 0;
92223bf6 818 LinkError(kTPmodeInvalid, "Just reading");
d60fe037 819 }
820
821 diff = *fPayloadCurr ^ expword;
822 if (diff != 0) {
92223bf6 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));;
d60fe037 826 }
827 fPayloadCurr++;
828 count++;
829 wordcount++;
830 }
831 channelcount++;
832 }
833 // continue with next MCM
834 }
835 return fPayloadCurr - start;
836}
837
838
839Int_t AliTRDrawStream::ReadZSData()
840{
841 // read the zs data from one link from the current reading position
842
843 UInt_t *start = fPayloadCurr;
844
845 Int_t mcmcount = 0;
0508ca31 846 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
d60fe037 847 Int_t channelcount = 0;
0508ca31 848 Int_t channelcountExp = 0;
849 Int_t channelcountMax = 0;
d60fe037 850 Int_t timebins;
851 Int_t currentTimebin = 0;
852 Int_t adcwc = 0;
853 Int_t evno = -1;
854 Int_t lastmcmpos = -1;
855 Int_t lastrobpos = -1;
856
857 if (fCurrNtimebins != fNtimebins) {
858 if (fNtimebins > 0)
92223bf6 859 LinkError(kNtimebinsChanged,
860 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
d60fe037 861 fNtimebins = fCurrNtimebins;
862 }
863
864 timebins = fNtimebins;
865
866 while (*(fPayloadCurr) != fgkDataEndmarker &&
867 fPayloadCurr - fPayloadStart < fPayloadSize) {
868
869 // ----- Checking MCM Header -----
870 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
cc26f39c 871 UInt_t *startPosMCM = fPayloadCurr;
d60fe037 872
873 // ----- checking for proper readout order - ROB -----
874 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
875 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
876 }
877 else {
92223bf6 878 ROBError(kPosUnexp);
d60fe037 879 }
880 fCurrRobPos = ROB(*fPayloadCurr);
881
882 // ----- checking for proper readout order - MCM -----
883 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) > lastmcmpos) {
884 lastmcmpos = GetMCMReadoutPos(lastmcmpos);
885 }
886 else {
92223bf6 887 MCMError(kPosUnexp);
d60fe037 888 }
889 fCurrMcmPos = MCM(*fPayloadCurr);
890
891 if (EvNo(*fPayloadCurr) != evno) {
892 if (evno == -1)
893 evno = EvNo(*fPayloadCurr);
894 else {
92223bf6 895 MCMError(kPtrgCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
d60fe037 896 }
897 }
898 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
899 Int_t padcoloff = PadColOffset(*fPayloadCurr);
900 Int_t row = Row(*fPayloadCurr);
901 fPayloadCurr++;
902
903 // ----- Reading ADC channels -----
904 AliDebug(2, Form("ADC mask: 0x%08x", *fPayloadCurr));
905
906 // ----- analysing the ADC mask -----
907 channelcount = 0;
0508ca31 908 channelcountExp = GetNActiveChannelsFromMask(*fPayloadCurr);
909 channelcountMax = GetNActiveChannels(*fPayloadCurr);
d60fe037 910 Int_t channelmask = GetActiveChannels(*fPayloadCurr);
911 Int_t channelno = -1;
912 fPayloadCurr++;
913
0508ca31 914 if (channelcountExp != channelcountMax) {
915 if (channelcountExp > channelcountMax) {
916 Int_t temp = channelcountExp;
917 channelcountExp = channelcountMax;
918 channelcountMax = temp;
d60fe037 919 }
0508ca31 920 while (channelcountExp < channelcountMax && channelcountExp < 21 &&
921 fPayloadCurr - fPayloadStart < fPayloadSize - 10 * channelcountExp - 1) {
92223bf6 922 MCMError(kAdcMaskInconsistent,
923 "Possible MCM-H: 0x%08x, possible ADC-mask: 0x%08x",
924 *(fPayloadCurr + 10 * channelcountExp),
925 *(fPayloadCurr + 10 * channelcountExp + 1) );
0508ca31 926 if (!CouldBeMCMhdr( *(fPayloadCurr + 10 * channelcountExp)) && !CouldBeADCmask( *(fPayloadCurr + 10 * channelcountExp + 1)))
927 channelcountExp++;
d60fe037 928 else {
929 break;
930 }
931 }
92223bf6 932 MCMError(kAdcMaskInconsistent,
933 "Inconsistency in no. of active channels: Counter: %i, Mask: %i, chosen: %i!",
934 GetNActiveChannels(fPayloadCurr[-1]), GetNActiveChannelsFromMask(fPayloadCurr[-1]), channelcountExp);
d60fe037 935 }
0508ca31 936 AliDebug(2, Form("expecting %i active channels, timebins: %i", channelcountExp, fCurrNtimebins));
d60fe037 937
938 // ----- reading marked ADC channels -----
0508ca31 939 while (channelcount < channelcountExp && *(fPayloadCurr) != fgkDataEndmarker) {
e29e514c 940 if (channelno < 20)
d60fe037 941 channelno++;
e29e514c 942 while (channelno < 20 && (channelmask & 1 << channelno) == 0)
d60fe037 943 channelno++;
944
945 if (fCurrNtimebins > 30) {
946 currentTimebin = ((*fPayloadCurr >> 2) & 0x3f);
947 timebins = ((*fPayloadCurr >> 8) & 0xf) * 3;
948 }
949 else {
950 currentTimebin = 0;
951 }
952
953 adcwc = 0;
954 AliDebug(2, Form("Now looking %i words", timebins / 3));
955 Int_t adccol = adccoloff - channelno;
956 Int_t padcol = padcoloff - channelno;
cc26f39c 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));
960
d60fe037 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) {
92223bf6 967 MCMError(kAdcCheckInvalid,
968 "%i for %2i. ADC word in odd channel %i",
969 check, adcwc+1, channelno);
d60fe037 970 }
971 }
972 else { // even channel
973 if (check != 0x3 && channelno < 21) {
92223bf6 974 MCMError(kAdcCheckInvalid,
975 "%i for %2i. ADC word in even channel %i",
976 check, adcwc+1, channelno);
d60fe037 977 }
978 }
979
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);
986 else
987 tb0 = -1;
988 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
989 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
990
991 adcwc++;
992 fPayloadCurr++;
993 }
994
995 if (adcwc != timebins / 3)
92223bf6 996 MCMError(kAdcDataAbort);
d60fe037 997
998 // adding index
999 if (padcol > 0 && padcol < 144) {
1000 fSignalIndex->AddIndexRC(row, padcol);
1001 }
1002
1003 channelcount++;
1004 }
cc26f39c 1005
637666cd 1006 if (fCurrSm > -1 && fCurrSm < 18) {
1007 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNChannels += channelcount;
1008 fStats.fStatsSector[fCurrSm].fNChannels += channelcount;
1009 }
0508ca31 1010 if (channelcount != channelcountExp)
92223bf6 1011 MCMError(kAdcChannelsMiss);
d60fe037 1012
1013 mcmcount++;
637666cd 1014 if (fCurrSm > -1 && fCurrSm < 18) {
1015 fStats.fStatsSector[fCurrSm].fStatsHC[fCurrHC%60].fNMCMs++;
1016 fStats.fStatsSector[fCurrSm].fNMCMs++;
1017 }
cc26f39c 1018
1019 if (IsDumping() && DumpingMCM(fCurrHC/2, fCurrRobPos, fCurrMcmPos)) {
92223bf6 1020 DumpRaw(Form("Event %i: Det %3i ROB %i MCM %2i", fRawReader->GetEventIndex(), fCurrHC/2, fCurrRobPos, fCurrMcmPos),
cc26f39c 1021 startPosMCM, fPayloadCurr - startPosMCM);
1022 }
1023
d60fe037 1024 // continue with next MCM
1025 }
1026
1027 // check for missing MCMs (if header suppression is inactive)
67271412 1028 if (((fCurrMajor & 0x1) == 0) && (mcmcount != mcmcountExp)) {
92223bf6 1029 LinkError(kMissMcmHeaders,
1030 "No. of MCM headers %i not as expected: %i",
1031 mcmcount, mcmcountExp);
d60fe037 1032 }
1033
1034 return (fPayloadCurr - start);
1035}
1036
1037Int_t AliTRDrawStream::ReadNonZSData()
1038{
1039 // read the non-zs data from one link from the current reading position
1040
1041 UInt_t *start = fPayloadCurr;
1042
1043 Int_t mcmcount = 0;
0508ca31 1044 Int_t mcmcountExp = fCurrStack == 2 ? 48 : 64;
d60fe037 1045 Int_t channelcount = 0;
0508ca31 1046 Int_t channelcountExp = 0;
d60fe037 1047 Int_t timebins;
1048 Int_t currentTimebin = 0;
1049 Int_t adcwc = 0;
1050 Int_t evno = -1;
1051 Int_t lastmcmpos = -1;
1052 Int_t lastrobpos = -1;
1053
1054 if (fCurrNtimebins != fNtimebins) {
1055 if (fNtimebins > 0)
92223bf6 1056 LinkError(kNtimebinsChanged,
1057 "No. of timebins changed from %i to %i", fNtimebins, fCurrNtimebins);
d60fe037 1058 fNtimebins = fCurrNtimebins;
1059 }
1060
1061 timebins = fNtimebins;
1062
1063 while (*(fPayloadCurr) != fgkDataEndmarker &&
1064 fPayloadCurr - fPayloadStart < fPayloadSize - 2) {
1065
1066 // ----- Checking MCM Header -----
1067 AliDebug(2, Form("MCM header: 0x%08x", *fPayloadCurr));
1068
1069 // ----- checking for proper readout order - ROB -----
1070 if (GetROBReadoutPos(ROB(*fPayloadCurr) / 2) >= lastrobpos) {
1071 lastrobpos = GetROBReadoutPos(ROB(*fPayloadCurr) / 2);
1072 }
1073 else {
92223bf6 1074 ROBError(kPosUnexp);
d60fe037 1075 }
1076 fCurrRobPos = ROB(*fPayloadCurr);
1077
1078 // ----- checking for proper readout order - MCM -----
1079 if (GetMCMReadoutPos(MCM(*fPayloadCurr)) >= (lastmcmpos + 1) % 16) {
1080 lastmcmpos = GetMCMReadoutPos(*fPayloadCurr);
1081 }
1082 else {
92223bf6 1083 MCMError(kPosUnexp);
d60fe037 1084 }
1085 fCurrMcmPos = MCM(*fPayloadCurr);
1086
1087 if (EvNo(*fPayloadCurr) != evno) {
1088 if (evno == -1)
1089 evno = EvNo(*fPayloadCurr);
1090 else {
92223bf6 1091 MCMError(kPtrgCntMismatch, "%i <-> %i", evno, EvNo(*fPayloadCurr));
d60fe037 1092 }
1093 }
1094
1095 channelcount = 0;
0508ca31 1096 channelcountExp = 21;
d60fe037 1097 int channelno = -1;
1098
1099 Int_t adccoloff = AdcColOffset(*fPayloadCurr);
1100 Int_t padcoloff = PadColOffset(*fPayloadCurr);
1101 Int_t row = Row(*fPayloadCurr);
1102
1103 fPayloadCurr++;
1104
1105 // ----- reading marked ADC channels -----
0508ca31 1106 while (channelcount < channelcountExp &&
d60fe037 1107 *(fPayloadCurr) != fgkDataEndmarker) {
e29e514c 1108 if (channelno < 20)
d60fe037 1109 channelno++;
1110
1111 currentTimebin = 0;
1112
1113 adcwc = 0;
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) {
92223bf6 1123 MCMError(kAdcCheckInvalid,
1124 "%i for %2i. ADC word in odd channel %i",
1125 check, adcwc+1, channelno);
d60fe037 1126 }
1127 }
1128 else { // even channel
1129 if (check != 0x3 && channelno < 21) {
92223bf6 1130 MCMError(kAdcCheckInvalid,
1131 "%i for %2i. ADC word in even channel %i",
1132 check, adcwc+1, channelno);
d60fe037 1133 }
1134 }
1135
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);
1142 else
1143 tb0 = -1;
1144 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb1);
1145 fAdcArray->SetDataByAdcCol(row, adccol, currentTimebin++, tb2);
1146
1147 adcwc++;
1148 fPayloadCurr++;
1149 }
1150
1151 if (adcwc != timebins / 3)
92223bf6 1152 MCMError(kAdcDataAbort);
d60fe037 1153
1154 // adding index
1155 if (padcol > 0 && padcol < 144) {
1156 fSignalIndex->AddIndexRC(row, padcol);
1157 }
1158
1159 channelcount++;
1160 }
1161
0508ca31 1162 if (channelcount != channelcountExp)
92223bf6 1163 MCMError(kAdcChannelsMiss);
d60fe037 1164 mcmcount++;
1165 // continue with next MCM
1166 }
1167
1168 // check for missing MCMs (if header suppression is inactive)
0508ca31 1169 if (mcmcount != mcmcountExp) {
92223bf6 1170 LinkError(kMissMcmHeaders,
1171 "%i not as expected: %i", mcmcount, mcmcountExp);
d60fe037 1172 }
1173
1174 return (fPayloadCurr - start);
1175}
1176
92223bf6 1177Int_t AliTRDrawStream::SeekNextLink()
1178{
1179 UInt_t *start = fPayloadCurr;
1180
1181 // read until data endmarkers
1182 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1183 *fPayloadCurr != fgkDataEndmarker)
1184 fPayloadCurr++;
1185
1186 // read all data endmarkers
1187 while (fPayloadCurr - fPayloadStart < fPayloadSize &&
1188 *fPayloadCurr == fgkDataEndmarker)
1189 fPayloadCurr++;
1190
1191 return (fPayloadCurr - start);
1192}
1193
67271412 1194Bool_t AliTRDrawStream::ConnectTracklets(TTree *trklTree)
1195{
1196 fTrackletTree = trklTree;
1197 if (!fTrackletTree)
1198 return kTRUE;
1199
cc26f39c 1200 if (!fTrackletTree->GetBranch("hc"))
67271412 1201 fTrackletTree->Branch("hc", &fCurrHC, "hc/I");
cc26f39c 1202 else
67271412 1203 fTrackletTree->SetBranchAddress("hc", &fCurrHC);
cc26f39c 1204
1205 if (!fTrackletTree->GetBranch("trkl"))
1206 fTrackletTree->Branch("trkl", &fTrackletArray);
1207 else
67271412 1208 fTrackletTree->SetBranchAddress("trkl", &fTrackletArray);
cc26f39c 1209
67271412 1210 return kTRUE;
1211}
1212
1213
1d62be37 1214void AliTRDrawStream::EquipmentError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1215{
0508ca31 1216 // register error according to error code on equipment level
1217 // and return the corresponding error message
1218
d60fe037 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;
f4b3235e 1225 (this->*fStoreError)();
d60fe037 1226
92223bf6 1227 va_list ap;
1228 if (fgErrorDebugLevel[err] > 10)
1229 AliDebug(fgErrorDebugLevel[err],
1230 Form("Event %6i: Eq. %2d - %s : %s",
1231 fRawReader->GetEventIndex(), fCurrEquipmentId, fgErrorMessages[err],
1d62be37 1232 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1233 else
1234 AliError(Form("Event %6i: Eq. %2d - %s : %s",
1235 fRawReader->GetEventIndex(), fCurrEquipmentId, fgErrorMessages[err],
1d62be37 1236 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1237 fErrorFlags |= fgErrorBehav[err];
d60fe037 1238}
1239
1240
1d62be37 1241void AliTRDrawStream::StackError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1242{
0508ca31 1243 // register error according to error code on stack level
1244 // and return the corresponding error message
1245
d60fe037 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;
f4b3235e 1252 (this->*fStoreError)();
d60fe037 1253
92223bf6 1254 va_list ap;
e29e514c 1255 if (fgErrorDebugLevel[err] > 0)
92223bf6 1256 AliDebug(fgErrorDebugLevel[err],
1257 Form("Event %6i: Eq. %2d S %i - %s : %s",
1258 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgErrorMessages[err],
1d62be37 1259 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1260 else
1261 AliError(Form("Event %6i: Eq. %2d S %i - %s : %s",
1262 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fgErrorMessages[err],
1d62be37 1263 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1264 fErrorFlags |= fgErrorBehav[err];
d60fe037 1265}
1266
1267
1d62be37 1268void AliTRDrawStream::LinkError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1269{
0508ca31 1270 // register error according to error code on link level
1271 // and return the corresponding error message
1272
d60fe037 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;
f4b3235e 1279 (this->*fStoreError)();
d60fe037 1280
92223bf6 1281 va_list ap;
e29e514c 1282 if (fgErrorDebugLevel[err] > 0)
92223bf6 1283 AliDebug(fgErrorDebugLevel[err],
1284 Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
1285 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgErrorMessages[err],
1d62be37 1286 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1287 else
1288 AliError(Form("Event %6i: Eq. %2d S %i l %2i - %s : %s",
1289 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fgErrorMessages[err],
1d62be37 1290 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1291 fErrorFlags |= fgErrorBehav[err];
d60fe037 1292}
1293
1294
1d62be37 1295void AliTRDrawStream::ROBError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1296{
0508ca31 1297 // register error according to error code on ROB level
1298 // and return the corresponding error message
1299
d60fe037 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;
f4b3235e 1306 (this->*fStoreError)();
d60fe037 1307
92223bf6 1308 va_list ap;
e29e514c 1309 if (fgErrorDebugLevel[err] > 0)
92223bf6 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],
1d62be37 1313 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1314 else
1315 AliError(Form("Event %6i: Eq. %2d S %i l %2i ROB %i - %s : %s",
1316 fRawReader->GetEventIndex(), fCurrEquipmentId, fCurrSlot, fCurrLink, fCurrRobPos, fgErrorMessages[err],
1d62be37 1317 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1318 fErrorFlags |= fgErrorBehav[err];
d60fe037 1319}
1320
1321
1d62be37 1322void AliTRDrawStream::MCMError(ErrorCode_t err, const char *const msg, ...)
d60fe037 1323{
0508ca31 1324 // register error according to error code on MCM level
1325 // and return the corresponding error message
1326
d60fe037 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;
f4b3235e 1333 (this->*fStoreError)();
d60fe037 1334
92223bf6 1335 va_list ap;
e29e514c 1336 if (fgErrorDebugLevel[err] > 0)
92223bf6 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],
1d62be37 1340 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1341 else
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],
1d62be37 1344 (va_start(ap, msg), vsprintf(fErrorBuffer, msg, ap), va_end(ap), fErrorBuffer) ));
92223bf6 1345 fErrorFlags |= fgErrorBehav[err];
d60fe037 1346}
1347
1348const char* AliTRDrawStream::GetErrorMessage(ErrorCode_t errCode)
1349{
0508ca31 1350 // return the error message for the given error code
1351
d60fe037 1352 if (errCode > 0 && errCode < kLastErrorCode)
1353 return fgErrorMessages[errCode];
1354 else
1355 return "";
1356}
cc26f39c 1357
1358void AliTRDrawStream::AliTRDrawStats::ClearStats()
1359{
1360 // clear statistics (includes clearing sector-wise statistics)
1361
1362 fBytesRead = 0;
1363 for (Int_t iSector = 0; iSector < 18; iSector++) {
1364 fStatsSector[iSector].ClearStats();
1365 }
1366
1367}
1368
1369void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::ClearStats()
1370{
1371 // clear statistics (includes clearing HC-wise statistics)
1372
1373 fBytes = 0;
1374 fBytesRead = 0;
1375 fNTracklets = 0;
1376 fNMCMs = 0;
1377 fNChannels = 0;
1378
1379 for (Int_t iHC = 0; iHC < 60; iHC++) {
1380 fStatsHC[iHC].ClearStats();
1381 }
1382}
1383
1384void AliTRDrawStream::AliTRDrawStats::AliTRDrawStatsSector::AliTRDrawStatsHC::ClearStats()
1385{
1386 // clear statistics
1387
1388 fBytes = 0;
1389 fBytesRead = 0;
1390 fNTracklets = 0;
1391 fNMCMs = 0;
1392 fNChannels = 0;
1393}
1394
1395void AliTRDrawStream::SetDumpMCM(Int_t det, Int_t rob, Int_t mcm, Bool_t dump)
1396{
1397 // mark MCM for dumping of raw data
1398
1399 if (dump) {
1400 fDumpMCM[fNDumpMCMs++] = (det << 7) | (rob << 4) | mcm;
1401 }
1402 else {
1403 Int_t iMCM;
1404 for (iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
1405 if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
1406 fNDumpMCMs--;
1407 break;
1408 }
1409 }
1410 for ( ; iMCM < fNDumpMCMs; iMCM++) {
1411 fDumpMCM[iMCM] = fDumpMCM[iMCM+1];
1412 }
1413 }
1414}
1415
1416Bool_t AliTRDrawStream::DumpingMCM(Int_t det, Int_t rob, Int_t mcm)
1417{
1418 // check if MCM data should be dumped
1419
1420 for (Int_t iMCM = 0; iMCM < fNDumpMCMs; iMCM++) {
1421 if (fDumpMCM[iMCM] == ((det << 7) | (rob << 4) | mcm)) {
1422 return kTRUE;
1423 }
1424 }
1425 return kFALSE;
1426}
1427
1428void AliTRDrawStream::DumpRaw(TString title, UInt_t *start, Int_t length)
1429{
1430 // dump raw data
1431
1432 title += "\n";
1433 Int_t pos = 0;
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]);
1437 }
1438 for ( ; pos < length; pos++) {
1439 title += Form("0x%08x ", start[pos]);
1440 }
1441 AliInfo(title);
1442}
f4b3235e 1443
1444AliTRDrawStream::AliTRDrawStreamError::AliTRDrawStreamError(Int_t error, Int_t sector, Int_t stack, Int_t link, Int_t rob, Int_t mcm) :
1445 fError(error),
1446 fSector(sector),
1447 fStack(stack),
1448 fLink(link),
1449 fRob(rob),
1450 fMcm(mcm)
1451{
1452 // ctor
1453
1454}