Line data Source code
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 : /* $Id$ */
17 :
18 : ///////////////////////////////////////////////////////////////////////////////
19 : // //
20 : // TRD raw data conversion class //
21 : // //
22 : ///////////////////////////////////////////////////////////////////////////////
23 :
24 :
25 : #include <TMath.h>
26 : #include "TClass.h"
27 :
28 : #include "AliDAQ.h"
29 : #include "AliRawDataHeaderSim.h"
30 : #include "AliRawReader.h"
31 : #include "AliLog.h"
32 : #include "AliFstream.h"
33 : #include "AliLoader.h"
34 : #include "AliTreeLoader.h"
35 :
36 : #include "AliTRDrawData.h"
37 : #include "AliTRDdigitsManager.h"
38 : #include "AliTRDgeometry.h"
39 : #include "AliTRDarrayDictionary.h"
40 : #include "AliTRDarrayADC.h"
41 : #include "AliTRDrawStream.h"
42 : #include "AliTRDcalibDB.h"
43 : #include "AliTRDSignalIndex.h"
44 : #include "AliTRDfeeParam.h"
45 : #include "AliTRDmcmSim.h"
46 : #include "AliTRDtrackletWord.h"
47 : #include "AliTRDtrackGTU.h"
48 : #include "AliESDTrdTrack.h"
49 : #include "AliTRDdigitsParam.h"
50 :
51 48 : ClassImp(AliTRDrawData)
52 :
53 : Int_t AliTRDrawData::fgDataSuppressionLevel = 0;
54 :
55 : //_____________________________________________________________________________
56 : AliTRDrawData::AliTRDrawData()
57 4 : :TObject()
58 4 : ,fRunLoader(NULL)
59 4 : ,fGeo(NULL)
60 4 : ,fFee(NULL)
61 4 : ,fNumberOfDDLs(0)
62 4 : ,fTrackletTree(NULL)
63 12 : ,fTracklets(new TClonesArray("AliTRDtrackletWord", 1000))
64 12 : ,fTracks(new TClonesArray("AliESDTrdTrack", 500))
65 4 : ,fSMindexPos(0)
66 4 : ,fStackindexPos(0)
67 4 : ,fEventCounter(0)
68 4 : ,fTrgFlags()
69 12 : ,fMcmSim(new AliTRDmcmSim)
70 4 : ,fDigitsParam(NULL)
71 20 : {
72 : //
73 : // Default constructor
74 : //
75 :
76 8 : fFee = AliTRDfeeParam::Instance();
77 8 : fNumberOfDDLs = AliDAQ::NumberOfDdls("TRD");
78 :
79 8 : }
80 :
81 : //_____________________________________________________________________________
82 : AliTRDrawData::AliTRDrawData(const AliTRDrawData &r)
83 0 : :TObject(r)
84 0 : ,fRunLoader(NULL)
85 0 : ,fGeo(NULL)
86 0 : ,fFee(NULL)
87 0 : ,fNumberOfDDLs(0)
88 0 : ,fTrackletTree(NULL)
89 0 : ,fTracklets(new TClonesArray("AliTRDtrackletWord", 1000))
90 0 : ,fTracks(new TClonesArray("AliESDTrdTrack", 500))
91 0 : ,fSMindexPos(0)
92 0 : ,fStackindexPos(0)
93 0 : ,fEventCounter(0)
94 0 : ,fTrgFlags()
95 0 : ,fMcmSim(new AliTRDmcmSim)
96 0 : ,fDigitsParam(NULL)
97 0 : {
98 : //
99 : // Copy constructor
100 : //
101 :
102 0 : fFee = AliTRDfeeParam::Instance();
103 0 : fNumberOfDDLs = AliDAQ::NumberOfDdls("TRD");
104 :
105 0 : }
106 :
107 : //_____________________________________________________________________________
108 : AliTRDrawData::~AliTRDrawData()
109 16 : {
110 : //
111 : // Destructor
112 : //
113 :
114 8 : delete fTracklets;
115 8 : delete fTracks;
116 8 : delete fMcmSim;
117 :
118 8 : }
119 :
120 : //_____________________________________________________________________________
121 : Bool_t AliTRDrawData::Digits2Raw(TTree *digitsTree, const TTree *tracks )
122 : {
123 : //
124 : // Initialize necessary parameters and call one
125 : // of the raw data simulator selected by SetRawVersion.
126 : //
127 : // Currently tracklet output is not spported yet and it
128 : // will be supported in higher version simulator.
129 : //
130 :
131 8 : AliTRDdigitsManager* const digitsManager = new AliTRDdigitsManager();
132 :
133 4 : if (!digitsManager->ReadDigits(digitsTree)) {
134 0 : delete digitsManager;
135 0 : return kFALSE;
136 : }
137 :
138 4 : if (tracks != NULL) {
139 0 : delete digitsManager;
140 0 : AliError("Tracklet input is not supported yet.");
141 0 : return kFALSE;
142 : }
143 :
144 8 : fGeo = new AliTRDgeometry();
145 :
146 4 : if (!AliTRDcalibDB::Instance()) {
147 0 : AliError("Could not get calibration object");
148 0 : delete fGeo;
149 0 : delete digitsManager;
150 0 : return kFALSE;
151 : }
152 :
153 : Int_t retval = kTRUE;
154 4 : Int_t rv = fFee->GetRAWversion();
155 :
156 : // Call appropriate Raw Simulator
157 8 : if ( rv > 0 && rv <= 3 ) retval = Digits2Raw(digitsManager);
158 : else {
159 : retval = kFALSE;
160 0 : AliWarning(Form("Unsupported raw version (%d).", rv));
161 : }
162 :
163 : // Cleanup
164 8 : delete fGeo;
165 8 : delete digitsManager;
166 :
167 4 : return retval;
168 :
169 4 : }
170 :
171 : //_____________________________________________________________________________
172 : Bool_t AliTRDrawData::Digits2Raw(AliTRDdigitsManager *digitsManager)
173 : {
174 : //
175 : // Raw data simulator for all versions > 0. This is prepared for real data.
176 : // This version simulate only raw data with ADC data and not with tracklet.
177 : //
178 :
179 8 : const Int_t kMaxHcWords = (fGeo->TBmax()/3)
180 4 : * fGeo->ADCmax()
181 4 : * fGeo->MCMmax()
182 4 : * fGeo->ROBmaxC1()/2
183 4 : + 100 + 20;
184 :
185 : // Buffer to temporary store half chamber data
186 4 : UInt_t *hcBuffer = new UInt_t[kMaxHcWords];
187 :
188 : Bool_t newEvent = kFALSE; // only for correct readout tree
189 : Bool_t newSM = kFALSE; // new SM flag, for writing SM index words
190 : // Coverity
191 : //Bool_t newStack = kFALSE; // new stack flag, for writing stack index words
192 :
193 : // Get digits parameter
194 4 : fDigitsParam = digitsManager->GetDigitsParam();
195 :
196 : // sect is same as iDDL, so I use only sect here.
197 152 : for (Int_t sect = 0; sect < fGeo->Nsector(); sect++) {
198 :
199 72 : char name[1024];
200 72 : snprintf(name,1024,"TRD_%d.ddl",sect + AliTRDrawStream::kDDLOffset);
201 :
202 72 : AliFstream* of = new AliFstream(name);
203 :
204 : // Write a dummy data header
205 72 : AliRawDataHeaderSim header; // the event header
206 72 : UInt_t hpos = of->Tellp();
207 72 : of->WriteBuffer((char *) (& header), sizeof(header));
208 :
209 : // Reset payload byte size (payload does not include header).
210 : Int_t npayloadbyte = 0;
211 :
212 : // check the existance of the data
213 : // SM index word and Stack index word
214 72 : UInt_t *iwbuffer = new UInt_t[109]; // index word buffer; max 109 = 2 SM headers + 67 dummy headers + 5*8 stack headers
215 : Int_t nheader = 0;
216 : UInt_t bStackMask = 0x0;
217 : Bool_t bStackHasData = kFALSE;
218 : Bool_t bSMHasData = kFALSE;
219 :
220 : //iwbuffer[nheader++] = 0x0001a020; // SM index words
221 72 : iwbuffer[nheader++] = 0x0044b020; // SM index words | additional SM header:48 = 1 SM header + 47 dummy words(for future use)
222 72 : iwbuffer[nheader++] = 0x10404071; // SM header
223 9648 : for ( Int_t i=0; i<66; i++ ) iwbuffer[nheader++] = 0x00000000; // dummy words
224 72 : iwbuffer[nheader++] = 0x10000000; // end of dummy words
225 :
226 864 : for ( Int_t stack= 0; stack < fGeo->Nstack(); stack++) {
227 : UInt_t linkMask = 0x0;
228 5040 : for( Int_t layer = 0; layer < fGeo->Nlayer(); layer++) {
229 2160 : Int_t iDet = fGeo->GetDetector(layer,stack,sect);
230 2160 : AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(iDet);
231 2160 : if ( fgDataSuppressionLevel==0 || digits->HasData() ) {
232 2160 : bStackMask = bStackMask | ( 1 << stack ); // active stack mask for new stack
233 2160 : linkMask = linkMask | ( 3 << (2*layer) ); // 3 = 0011
234 : bStackHasData = kTRUE;
235 : bSMHasData = kTRUE;
236 2160 : } // has data
237 : } // loop over layer
238 :
239 360 : if ( fgDataSuppressionLevel==0 || bStackHasData ){
240 360 : iwbuffer[nheader++] = 0x0007a000 | linkMask; // stack index word + link masks
241 720 : if (fgDataSuppressionLevel==0) iwbuffer[nheader-1] = 0x0007afff; // no suppression
242 360 : iwbuffer[nheader++] = 0x04045b01; // stack header
243 5040 : for (Int_t i=0;i<6;i++) iwbuffer[nheader++] = 0x00000000; // 6 dummy words
244 : bStackHasData = kFALSE;
245 360 : }
246 : } // loop over stack
247 :
248 72 : if ( fgDataSuppressionLevel==0 || bSMHasData ){
249 72 : iwbuffer[0] = iwbuffer[0] | bStackMask; // add stack masks to SM index word
250 144 : if (fgDataSuppressionLevel==0) iwbuffer[0] = 0x0044b03f; // no suppression : all stacks are active
251 72 : of->WriteBuffer((char *) iwbuffer, nheader*4);
252 216 : AliDebug(11, Form("SM %d index word: %08x", sect, iwbuffer[0]));
253 216 : AliDebug(11, Form("SM %d header: %08x", sect, iwbuffer[1]));
254 : }
255 : // end of SM & stack header ------------------------------------------------------------------------
256 : // -------------------------------------------------------------------------------------------------
257 :
258 : // Prepare chamber data
259 864 : for( Int_t stack = 0; stack < fGeo->Nstack(); stack++) {
260 5040 : for( Int_t layer = 0; layer < fGeo->Nlayer(); layer++) {
261 :
262 2160 : Int_t iDet = fGeo->GetDetector(layer,stack,sect);
263 2160 : if (iDet == 0){
264 : newEvent = kTRUE; // it is expected that each event has at least one tracklet;
265 : // this is only needed for correct readout tree
266 4 : fEventCounter++;
267 12 : AliDebug(11, Form("New event!! Event counter: %d",fEventCounter));
268 : }
269 :
270 2232 : if ( stack==0 && layer==0 ) newSM = kTRUE; // new SM flag
271 : // Coverity
272 : //if ( layer==0 ) newStack = kTRUE; // new stack flag
273 6480 : AliDebug(15, Form("stack : %d, layer : %d, iDec : %d\n",stack,layer,iDet));
274 : // Get the digits array
275 2160 : AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(iDet);
276 2160 : if (fgDataSuppressionLevel==0 || digits->HasData() ) { // second part is new!! and is for indicating a new event
277 :
278 2348 : if (digits->HasData()) digits->Expand();
279 :
280 : Int_t hcwords = 0;
281 :
282 : // Process A side of the chamber
283 2160 : hcwords = ProduceHcData(digits,0,iDet,hcBuffer,kMaxHcWords,newEvent,newSM);
284 2164 : if ( newEvent ) newEvent = kFALSE;
285 : //AssignLinkMask(hcBuffer, layer); // active link mask for this layer(2*HC)
286 2160 : of->WriteBuffer((char *) hcBuffer, hcwords*4);
287 2160 : npayloadbyte += hcwords*4;
288 :
289 : // Process B side of the chamber
290 2160 : hcwords = ProduceHcData(digits,1,iDet,hcBuffer,kMaxHcWords,newEvent,newSM);
291 2160 : of->WriteBuffer((char *) hcBuffer, hcwords*4);
292 2160 : npayloadbyte += hcwords*4;
293 2160 : } // has data
294 :
295 : } // loop over layer
296 : } // loop over stack
297 :
298 : // Complete header
299 72 : header.fSize = UInt_t(of->Tellp()) - hpos;
300 72 : header.SetAttribute(0); // Valid data
301 72 : of->Seekp(hpos); // Rewind to header position
302 72 : of->WriteBuffer((char *) (& header), sizeof(header));
303 144 : delete of;
304 :
305 144 : delete [] iwbuffer;
306 :
307 72 : } // loop over sector(SM)
308 :
309 8 : delete [] hcBuffer;
310 :
311 4 : return kTRUE;
312 0 : }
313 :
314 : //_____________________________________________________________________________
315 : void AliTRDrawData::ProduceSMIndexData(UInt_t *buf, Int_t& nw){
316 : //
317 : // This function generates
318 : // 1) SM index words : ssssssss ssssssss vvvv rrrr r d t mmmmm
319 : // - s : size of SM header (number of header, default = 0x0001)
320 : // - v : SM header version (default = 0xa)
321 : // - r : reserved for future use (default = 00000)
322 : // - d : track data enabled bit (default = 0)
323 : // - t : tracklet data enabled bit (default = 1)
324 : // - m : stack mask (each bit corresponds a stack, default = 11111)
325 : //
326 : // 2) SM header : rrr c vvvv vvvvvvvv vvvv rrrr bbbbbbbb
327 : // - r : reserved for future use (default = 000)
328 : // - c : clean check out flag (default = 1)
329 : // - v : hardware design revision (default = 0x0404)
330 : // - r : reserved for future use (default = 0x0)
331 : // - b : physical board ID (default = 0x71)
332 : //
333 : // 3) stack index words : ssssssss ssssssss vvvv mmmm mmmmmmmm
334 : // - s : size of stack header (number of header, (default = 0x0007)
335 : // - v : header version (default = 0xa)
336 : // - m : link mask (default = 0xfff)
337 : //
338 : // 4) stack header : vvvvvvvv vvvvvvvv bbbbbbbb rrrr rrr c
339 : // - v : hardware design revision (default = 0x0404)
340 : // - b : physical board ID (default = 0x5b)
341 : // - r : reserved for future use (default = 0000 000)
342 : // - c : clean checkout flag (default = 1)
343 : //
344 : // and 6 dummy words(0x00000000)
345 : //
346 :
347 : //buf[nw++] = 0x0001a03f; // SM index words
348 0 : fSMindexPos = nw; // memorize position of the SM index word for re-allocating stack mask
349 0 : buf[nw++] = 0x0001b020; // SM index words
350 0 : buf[nw++] = 0x10404071; // SM header
351 :
352 0 : fStackindexPos = nw; // memorize position of the stack index word for future adding
353 : /*
354 : for (Int_t istack=0; istack<5; istack++){
355 : buf[nw++] = 0x0007afff; // stack index words
356 : buf[nw++] = 0x04045b01; // stack header
357 : for (Int_t i=0;i<6;i++) buf[nw++] = 0x00000000; // 6 dummy words
358 : } // loop over 5 stacks
359 : */
360 0 : }
361 :
362 : //_____________________________________________________________________________
363 : void AliTRDrawData::AssignStackMask(UInt_t *buf, Int_t nStack){
364 : //
365 : // This function re-assign stack mask active(from 0 to 1) in the SM index word
366 : //
367 0 : buf[fSMindexPos] = buf[fSMindexPos] | ( 1 << nStack );
368 0 : }
369 :
370 : //_____________________________________________________________________________
371 : Int_t AliTRDrawData::AddStackIndexWords(UInt_t *buf, Int_t /*nStack*/, Int_t nMax){
372 : //
373 : // This function add stack index words and stack header when there is data for the stack
374 : //
375 : // 1) stack index words : ssssssss ssssssss vvvv mmmm mmmmmmmm
376 : // - s : size of stack header (number of header, (default = 0x0007)
377 : // - v : header version (default = 0xa)
378 : // - m : link mask (default = 0xfff)
379 : // - m : link mask (starting value = 0x000)
380 : //
381 : // 2) stack header : vvvvvvvv vvvvvvvv bbbbbbbb rrrr rrr c
382 : // - v : hardware design revision (default = 0x0404)
383 : // - b : physical board ID (default = 0x5b)
384 : // - r : reserved for future use (default = 0000 000)
385 : // - c : clean checkout flag (default = 1)
386 : //
387 : // and 6 dummy words(0x00000000)
388 : //
389 :
390 : Int_t nAddedWords = 0; // Number of added words
391 0 : if ( ShiftWords(buf, fStackindexPos, 8, nMax)== kFALSE ){
392 0 : AliError("Adding stack header failed.");
393 0 : return 0;
394 : }
395 :
396 0 : buf[fStackindexPos++] = 0x0007a000; // stack index words
397 0 : buf[fStackindexPos++] = 0x04045b01; // stack header
398 0 : for (Int_t i=0;i<6;i++) buf[fStackindexPos++] = 0x00000000; // 6 dummy words
399 : nAddedWords += 8;
400 :
401 0 : return nAddedWords;
402 0 : }
403 :
404 : //_____________________________________________________________________________
405 : void AliTRDrawData::AssignLinkMask(UInt_t *buf, Int_t nLayer){
406 : //
407 : // This function re-assign link mask active(from 0 to 1) in the stack index word
408 : //
409 0 : buf[fStackindexPos-8] = buf[fStackindexPos-8] | ( 3 << (2*nLayer) ); // 3 = 0011
410 0 : }
411 :
412 : //_____________________________________________________________________________
413 : Bool_t AliTRDrawData::ShiftWords(UInt_t *buf, Int_t nStart, Int_t nWords, Int_t nMax){
414 : //
415 : // This function shifts n words
416 : //
417 : //if ( nStart+nWords > sizeof(buf)/sizeof(UInt_t) ){
418 : // AliError("Words shift failed. No more buffer space.");
419 : // return kFALSE;
420 : //}
421 :
422 0 : for ( Int_t iw=nMax; iw>nStart-1; iw--){
423 0 : buf[iw+nWords] = buf[iw];
424 : }
425 0 : return kTRUE;
426 : }
427 :
428 : //_____________________________________________________________________________
429 : Int_t AliTRDrawData::ProduceHcData(AliTRDarrayADC *digits, Int_t side, Int_t det, UInt_t *buf, Int_t maxSize, Bool_t /*newEvent = kFALSE*/, Bool_t /*newSM = kFALSE*/){
430 : //
431 : // This function produces the raw data for one HC, i.e. tracklets, tracklet endmarkers,
432 : // raw data, raw data endmarkers.
433 : // This function can be used for both ZS and NZS data
434 : //
435 :
436 4320 : Int_t nw = 0; // Number of written words
437 4320 : Int_t of = 0; // Number of overflowed words
438 : Int_t *tempnw = &nw; // Number of written words for temp. buffer
439 : Int_t *tempof = &of; // Number of overflowed words for temp. buffer
440 4320 : Int_t layer = fGeo->GetLayer( det ); // Layer
441 4320 : Int_t stack = fGeo->GetStack( det ); // Stack
442 4320 : Int_t sect = fGeo->GetSector( det ); // Sector (=iDDL)
443 4320 : const Int_t kCtype = fGeo->GetStack(det) == 2 ? 0 : 1; // Chamber type (0:C0, 1:C1)
444 :
445 4320 : Bool_t trackletOn = fFee->GetTracklet(); // tracklet simulation active?
446 :
447 12960 : AliDebug(1,Form("Producing raw data for sect=%d layer=%d stack=%d side=%d",sect,layer,stack,side));
448 :
449 : UInt_t *tempBuffer = buf; // tempBuffer used to write ADC data
450 : // different in case of tracklet writing
451 :
452 4320 : if (trackletOn) {
453 4320 : tempBuffer = new UInt_t[maxSize];
454 4320 : tempnw = new Int_t(0);
455 4320 : tempof = new Int_t(0);
456 4320 : }
457 :
458 4320 : WriteIntermediateWords(tempBuffer,*tempnw,*tempof,maxSize,det,side);
459 :
460 4320 : if (digits->HasData()) {
461 : // scanning direction such, that tracklet-words are sorted in ascending z and then in ascending y order
462 : // ROB numbering on chamber and MCM numbering on ROB increase with decreasing z and increasing y
463 3352 : for (Int_t iRobRow = 0; iRobRow <= (kCtype + 3)-1; iRobRow++ ) {
464 : // ROB number should be increasing
465 1300 : Int_t iRob = iRobRow * 2 + side;
466 : // MCM on one ROB
467 44200 : for (Int_t iMcmRB = 0; iMcmRB < fGeo->MCMmax(); iMcmRB++ ) {
468 20800 : Int_t iMcm = 16 - 4*(iMcmRB/4 + 1) + (iMcmRB%4);
469 :
470 20800 : fMcmSim->Init(det, iRob, iMcm);
471 20800 : fMcmSim->SetData(digits); // no filtering done here (already done in digitizer)
472 20800 : if (trackletOn) {
473 20800 : fMcmSim->Tracklet();
474 20800 : Int_t tempNw = fMcmSim->ProduceTrackletStream(&buf[nw], maxSize - nw);
475 20800 : if( tempNw < 0 ) {
476 0 : of += tempNw;
477 0 : nw += maxSize - nw;
478 0 : AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
479 0 : } else {
480 20800 : nw += tempNw;
481 : }
482 20800 : }
483 20800 : fMcmSim->ZSMapping(); // Calculate zero suppression mapping
484 : // at the moment it has to be rerun here
485 : // Write MCM data to temp. buffer
486 20800 : Int_t tempNw = fMcmSim->ProduceRawStream( &tempBuffer[*tempnw], maxSize - *tempnw, fEventCounter );
487 20800 : if ( tempNw < 0 ) {
488 0 : *tempof += tempNw;
489 0 : *tempnw += maxSize - nw;
490 0 : AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
491 0 : } else {
492 20800 : *tempnw += tempNw;
493 : }
494 : }
495 : }
496 376 : }
497 :
498 : // in case of tracklet writing copy temp data to final buffer
499 4320 : if (trackletOn) {
500 4320 : if (nw + *tempnw < maxSize) {
501 4320 : memcpy(&buf[nw], tempBuffer, *tempnw * sizeof(UInt_t));
502 4320 : nw += *tempnw;
503 4320 : }
504 : else {
505 0 : AliError("Buffer overflow detected");
506 : }
507 :
508 8640 : delete [] tempBuffer;
509 8640 : delete tempof;
510 8640 : delete tempnw;
511 : }
512 :
513 : // Write end of raw data marker
514 4320 : if (nw+3 < maxSize) {
515 4320 : buf[nw++] = fgkEndOfDataMarker;
516 4320 : buf[nw++] = fgkEndOfDataMarker;
517 4320 : buf[nw++] = fgkEndOfDataMarker;
518 4320 : buf[nw++] = fgkEndOfDataMarker;
519 4320 : } else {
520 0 : of += 4;
521 : }
522 :
523 4320 : if (of != 0) {
524 0 : AliError("Buffer overflow. Data is truncated. Please increase buffer size and recompile.");
525 0 : }
526 :
527 8640 : return nw;
528 4320 : }
529 :
530 : //_____________________________________________________________________________
531 : AliTRDdigitsManager *AliTRDrawData::Raw2Digits(AliRawReader *rawReader)
532 : {
533 : //
534 : // Vx of the raw data reading
535 : //
536 :
537 0 : rawReader->Select("TRD"); //[mj]
538 :
539 : AliTRDarrayADC *digits = 0;
540 : AliTRDarrayDictionary *track0 = 0;
541 : AliTRDarrayDictionary *track1 = 0;
542 : AliTRDarrayDictionary *track2 = 0;
543 :
544 : //AliTRDSignalIndex *indexes = 0;
545 : // Create the digits manager
546 0 : AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
547 0 : digitsManager->CreateArrays();
548 :
549 0 : AliTRDrawStream input(rawReader);
550 :
551 0 : AliRunLoader *runLoader = AliRunLoader::Instance();
552 :
553 : // ----- preparing tracklet output -----
554 0 : AliDataLoader *trklLoader = runLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets");
555 0 : if (!trklLoader) {
556 0 : trklLoader = new AliDataLoader("TRD.Tracklets.root","tracklets", "tracklets");
557 0 : runLoader->GetLoader("TRDLoader")->AddDataLoader(trklLoader);
558 : }
559 0 : if (!trklLoader) {
560 0 : AliError("Could not get the tracklet data loader!");
561 0 : return 0x0;
562 : }
563 : // trklLoader->SetEvent();
564 0 : trklLoader->Load("update");
565 0 : AliTreeLoader *trklTreeLoader = dynamic_cast<AliTreeLoader*> (trklLoader->GetBaseLoader("tracklets-raw"));
566 0 : if (!trklTreeLoader) {
567 0 : trklTreeLoader = new AliTreeLoader("tracklets-raw", trklLoader);
568 0 : trklLoader->AddBaseLoader(trklTreeLoader);
569 : }
570 :
571 0 : TClonesArray *fTrackletArray = new TClonesArray("AliTRDtrackletWord", 256);
572 0 : Int_t lastHC = -1;
573 0 : if (trklLoader->Tree())
574 0 : trklLoader->MakeTree();
575 0 : TTree *trackletTree = trklTreeLoader->Tree();
576 0 : if (!trackletTree) {
577 0 : trklTreeLoader->MakeTree();
578 0 : trackletTree = trklTreeLoader->Tree();
579 0 : }
580 0 : if (!trackletTree->GetBranch("hc"))
581 0 : trackletTree->Branch("hc", &lastHC, "hc/I");
582 : else
583 0 : trackletTree->SetBranchAddress("hc", &lastHC);
584 :
585 0 : if (!trackletTree->GetBranch("trkl"))
586 0 : trackletTree->Branch("trkl", &fTrackletArray);
587 : else
588 0 : trackletTree->SetBranchAddress("trkl", &fTrackletArray);
589 :
590 : // ----- preparing track output -----
591 : AliDataLoader *trkLoader = 0x0;
592 0 : trkLoader = runLoader->GetLoader("TRDLoader")->GetDataLoader("gtutracks");
593 0 : if (!trkLoader) {
594 0 : trkLoader = new AliDataLoader("TRD.GtuTracks.root","gtutracks", "gtutracks");
595 0 : runLoader->GetLoader("TRDLoader")->AddDataLoader(trkLoader);
596 : }
597 0 : if (!trkLoader) {
598 0 : AliError("Could not get the GTU-track data loader!");
599 0 : return 0x0;
600 : }
601 0 : trkLoader->Load("UPDATE");
602 :
603 0 : TTree *trackTree = trkLoader->Tree();
604 0 : if (!trackTree) {
605 0 : trkLoader->MakeTree();
606 0 : trackTree = trkLoader->Tree();
607 0 : }
608 :
609 0 : AliTRDtrackGTU trk;
610 0 : if (!trackTree->GetBranch("TRDtrackGTU"))
611 0 : trackTree->Branch("TRDtrackGTU", "AliTRDtrackGTU", &trk, 32000);
612 :
613 : // ----- read the raw data -----
614 : // write the tracklets to arrays while reading raw data
615 0 : input.SetTrackletArray(fTracklets);
616 0 : input.SetTrackArray(fTracks);
617 :
618 : // Loop through the digits
619 : Int_t det = 0;
620 :
621 0 : while (det >= 0)
622 : {
623 0 : det = input.NextChamber(digitsManager);
624 :
625 0 : if (det >= 0)
626 : {
627 : // get...
628 0 : digits = (AliTRDarrayADC *) digitsManager->GetDigits(det);
629 0 : track0 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,0);
630 0 : track1 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,1);
631 0 : track2 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,2);
632 : // and compress
633 0 : if (digits) digits->Compress();
634 0 : if (track0) track0->Compress();
635 0 : if (track1) track1->Compress();
636 0 : if (track2) track2->Compress();
637 : }
638 : }
639 :
640 0 : for (Int_t iSector = 0; iSector < fGeo->Nsector(); iSector++) {
641 0 : fTrgFlags[iSector] = input.GetTriggerFlags(iSector);
642 : }
643 :
644 : // ----- tracklet output -----
645 0 : fTrackletArray->Clear();
646 0 : Int_t nTracklets = fTracklets ? fTracklets->GetEntries() : 0;
647 0 : AliDebug(1, Form("Writing %i tracklets to loader", nTracklets));
648 0 : for (Int_t iTracklet = 0; iTracklet < nTracklets; ++iTracklet) {
649 0 : AliTRDtrackletWord *trkl = (AliTRDtrackletWord*) fTracklets->At(iTracklet);
650 0 : if (trkl->GetHCId() != lastHC) {
651 0 : if (fTrackletArray->GetEntriesFast() > 0) {
652 0 : trackletTree->Fill();
653 0 : fTrackletArray->Clear();
654 : }
655 0 : lastHC = trkl->GetHCId();
656 0 : }
657 0 : new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(*trkl);
658 : }
659 0 : if (fTrackletArray->GetEntriesFast() > 0) {
660 0 : trackletTree->Fill();
661 0 : fTrackletArray->Clear();
662 : }
663 0 : trklTreeLoader->WriteData("OVERWRITE");
664 0 : trklLoader->WriteData("OVERWRITE");
665 0 : trklLoader->UnloadAll();
666 0 : trklLoader->CloseFile();
667 :
668 : // ----- track output -----
669 0 : Int_t nTracks = fTracks ? fTracks->GetEntriesFast() : 0;
670 0 : AliDebug(1, Form("Writing %i tracks to loader", nTracks));
671 0 : for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
672 0 : AliESDTrdTrack *trkEsd = (AliESDTrdTrack*) fTracks->At(iTrack);
673 0 : trk = *trkEsd;
674 0 : trackTree->Fill();
675 : }
676 0 : trkLoader->WriteData("OVERWRITE");
677 0 : trkLoader->UnloadAll();
678 0 : trkLoader->CloseFile();
679 :
680 : return digitsManager;
681 0 : }
682 :
683 : //_____________________________________________________________________________
684 : void AliTRDrawData::WriteIntermediateWords(UInt_t* buf, Int_t& nw, Int_t& of, const Int_t& maxSize, const Int_t& det, const Int_t& side) {
685 : //
686 : // write tracklet end marker(0x10001000)
687 : // and half chamber headers(H[0] and H[1])
688 : //
689 :
690 8640 : Int_t layer = fGeo->GetLayer( det ); // Layer
691 4320 : Int_t stack = fGeo->GetStack( det ); // Stack
692 4320 : Int_t sect = fGeo->GetSector( det ); // Sector (=iDDL)
693 4320 : Int_t rv = fFee->GetRAWversion();
694 4320 : const Int_t kNTBin = fDigitsParam->GetNTimeBins(det);
695 4320 : Bool_t trackletOn = fFee->GetTracklet();
696 : UInt_t x = 0;
697 :
698 : // Write end of tracklet marker
699 4320 : if (nw < maxSize){
700 4320 : buf[nw++] = fgkEndOfTrackletMarker;
701 4320 : buf[nw++] = fgkEndOfTrackletMarker; // the number of tracklet end marker should be more than 2
702 4320 : }
703 : else {
704 0 : of++;
705 : }
706 :
707 : // Half Chamber header
708 : // h[0] (there are 2 HC headers) xmmm mmmm nnnn nnnq qqss sssp ppcc ci01
709 : // , where x : Raw version speacial number (=1)
710 : // m : Raw version major number (test pattern, ZS, disable tracklet, 0, options)
711 : // n : Raw version minor number
712 : // q : number of addtional header words (default = 1)
713 : // s : SM sector number (ALICE numbering)
714 : // p : plane(layer) number
715 : // c : chamber(stack) number
716 : // i : side number (0:A, 1:B)
717 : Int_t majorv = 0; // The major version number
718 : Int_t minorv = 0; // The minor version number
719 : Int_t add = 1; // The number of additional header words to follow : now 1, previous 2
720 : Int_t tp = 0; // test pattern (default=0)
721 4320 : Int_t zs = (rv==3) ? 1 : 0; // zero suppression
722 4320 : Int_t dt = (trackletOn) ? 0 : 1; // disable tracklet
723 :
724 4320 : majorv = (tp<<6) | (zs<<5) | (dt<<4) | 1; // major version
725 :
726 4320 : x = (1<<31) | (majorv<<24) | (minorv<<17) | (add<<14) | (sect<<9) | (layer<<6) | (stack<<3) | (side<<2) | 1;
727 8640 : if (nw < maxSize) buf[nw++] = x; else of++;
728 :
729 : // h[1] tttt ttbb bbbb bbbb bbbb bbpp pphh hh01
730 : // , where t : number of time bins
731 : // b : bunch crossing number
732 : // p : pretrigger counter
733 : // h : pretrigger phase
734 : Int_t bcCtr = 99; // bunch crossing counter. Here it is set to 99 always for no reason
735 : Int_t ptCtr = 15; // pretrigger counter. Here it is set to 15 always for no reason
736 : Int_t ptPhase = 11; // pretrigger phase. Here it is set to 11 always for no reason
737 : //x = (bcCtr<<16) | (ptCtr<<12) | (ptPhase<<8) | ((kNTBin-1)<<2) | 1; // old format
738 4320 : x = ((kNTBin)<<26) | (bcCtr<<10) | (ptCtr<<6) | (ptPhase<<2) | 1;
739 8640 : if (nw < maxSize) buf[nw++] = x; else of++;
740 4320 : }
|