]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDrawData.cxx
Add TRD test class for G4 simulations
[u/mrichter/AliRoot.git] / TRD / AliTRDrawData.cxx
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 ClassImp(AliTRDrawData)
52
53 Int_t AliTRDrawData::fgDataSuppressionLevel = 0;
54
55 //_____________________________________________________________________________
56 AliTRDrawData::AliTRDrawData()
57   :TObject()
58   ,fRunLoader(NULL)
59   ,fGeo(NULL)
60   ,fFee(NULL)
61   ,fNumberOfDDLs(0)
62   ,fTrackletTree(NULL)
63   ,fTracklets(new TClonesArray("AliTRDtrackletWord", 1000))
64   ,fTracks(new TClonesArray("AliESDTrdTrack", 500))
65   ,fSMindexPos(0)
66   ,fStackindexPos(0)
67   ,fEventCounter(0)
68   ,fTrgFlags()
69   ,fMcmSim(new AliTRDmcmSim)
70   ,fDigitsParam(NULL)
71 {
72   //
73   // Default constructor
74   //
75
76   fFee = AliTRDfeeParam::Instance();
77   fNumberOfDDLs = AliDAQ::NumberOfDdls("TRD");
78
79 }
80
81 //_____________________________________________________________________________
82 AliTRDrawData::AliTRDrawData(const AliTRDrawData &r)
83   :TObject(r)
84   ,fRunLoader(NULL)
85   ,fGeo(NULL)
86   ,fFee(NULL)
87   ,fNumberOfDDLs(0)
88   ,fTrackletTree(NULL)
89   ,fTracklets(new TClonesArray("AliTRDtrackletWord", 1000))
90   ,fTracks(new TClonesArray("AliESDTrdTrack", 500))
91   ,fSMindexPos(0)
92   ,fStackindexPos(0)
93   ,fEventCounter(0)
94   ,fTrgFlags()
95   ,fMcmSim(new AliTRDmcmSim)
96   ,fDigitsParam(NULL)
97 {
98   //
99   // Copy constructor
100   //
101
102   fFee = AliTRDfeeParam::Instance();
103   fNumberOfDDLs = AliDAQ::NumberOfDdls("TRD");
104
105 }
106
107 //_____________________________________________________________________________
108 AliTRDrawData::~AliTRDrawData()
109 {
110   //
111   // Destructor
112   //
113
114   delete fTracklets;
115   delete fTracks;
116   delete fMcmSim;
117
118 }
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   AliTRDdigitsManager* const digitsManager = new AliTRDdigitsManager();
132
133   if (!digitsManager->ReadDigits(digitsTree)) {
134     delete digitsManager;
135     return kFALSE;
136   }
137
138   if (tracks != NULL) {
139     delete digitsManager;
140     AliError("Tracklet input is not supported yet.");
141     return kFALSE;
142   }
143
144   fGeo = new AliTRDgeometry();
145
146   if (!AliTRDcalibDB::Instance()) {
147     AliError("Could not get calibration object");
148     delete fGeo;
149     delete digitsManager;
150     return kFALSE;
151   }
152
153   Int_t retval = kTRUE;
154   Int_t rv     = fFee->GetRAWversion();
155
156   // Call appropriate Raw Simulator
157   if ( rv > 0 && rv <= 3 ) retval = Digits2Raw(digitsManager); 
158   else {
159     retval = kFALSE;
160     AliWarning(Form("Unsupported raw version (%d).", rv));
161   }
162
163   // Cleanup
164   delete fGeo;
165   delete digitsManager;
166
167   return retval;
168
169 }
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   const Int_t kMaxHcWords = (fGeo->TBmax()/3)
180                           * fGeo->ADCmax()
181                           * fGeo->MCMmax()
182                           * fGeo->ROBmaxC1()/2 
183                           + 100 + 20;
184
185   // Buffer to temporary store half chamber data
186   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   fDigitsParam = digitsManager->GetDigitsParam();
195
196   // sect is same as iDDL, so I use only sect here.
197   for (Int_t sect = 0; sect < fGeo->Nsector(); sect++) { 
198
199     char name[1024];
200     snprintf(name,1024,"TRD_%d.ddl",sect + AliTRDrawStream::kDDLOffset);
201
202     AliFstream* of = new AliFstream(name);
203
204     // Write a dummy data header
205     AliRawDataHeaderSim  header;  // the event header
206     UInt_t hpos = of->Tellp();
207     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     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     iwbuffer[nheader++] = 0x0044b020;   // SM index words | additional SM header:48 = 1 SM header + 47 dummy words(for future use)
222     iwbuffer[nheader++] = 0x10404071;   // SM header
223     for ( Int_t i=0; i<66; i++ ) iwbuffer[nheader++] = 0x00000000;  // dummy words 
224     iwbuffer[nheader++] = 0x10000000;   // end of dummy words
225     
226     for ( Int_t stack= 0; stack < fGeo->Nstack(); stack++) {
227       UInt_t linkMask = 0x0;
228       for( Int_t layer = 0; layer < fGeo->Nlayer(); layer++) {
229         Int_t iDet = fGeo->GetDetector(layer,stack,sect);
230         AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(iDet);
231         if ( fgDataSuppressionLevel==0 || digits->HasData() ) {
232           bStackMask = bStackMask | ( 1 << stack ); // active stack mask for new stack
233           linkMask = linkMask | ( 3 << (2*layer) );    // 3 = 0011
234           bStackHasData = kTRUE;
235           bSMHasData = kTRUE;
236         } // has data
237       } // loop over layer
238       
239       if ( fgDataSuppressionLevel==0 || bStackHasData ){
240         iwbuffer[nheader++] = 0x0007a000 | linkMask;    // stack index word + link masks
241         if (fgDataSuppressionLevel==0) iwbuffer[nheader-1] = 0x0007afff;  // no suppression
242         iwbuffer[nheader++] = 0x04045b01;               // stack header
243         for (Int_t i=0;i<6;i++) iwbuffer[nheader++] = 0x00000000; // 6 dummy words
244         bStackHasData = kFALSE;
245       }
246     } // loop over stack
247     
248     if ( fgDataSuppressionLevel==0 || bSMHasData ){
249       iwbuffer[0] = iwbuffer[0] | bStackMask;  // add stack masks to SM index word
250       if (fgDataSuppressionLevel==0) iwbuffer[0] = 0x0044b03f;    // no suppression : all stacks are active
251       of->WriteBuffer((char *) iwbuffer, nheader*4);
252       AliDebug(11, Form("SM %d index word: %08x", sect, iwbuffer[0]));
253       AliDebug(11, Form("SM %d header: %08x", sect, iwbuffer[1]));
254     }
255     // end of SM & stack header ------------------------------------------------------------------------
256     // -------------------------------------------------------------------------------------------------
257     
258     // Prepare chamber data
259     for( Int_t stack = 0; stack < fGeo->Nstack(); stack++) {
260       for( Int_t layer = 0; layer < fGeo->Nlayer(); layer++) {
261         
262         Int_t iDet = fGeo->GetDetector(layer,stack,sect);
263         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           fEventCounter++;
267           AliDebug(11, Form("New event!! Event counter: %d",fEventCounter));
268         }
269         
270         if ( stack==0 && layer==0 ) newSM = kTRUE;  // new SM flag
271         // Coverity
272         //if ( layer==0 ) newStack = kTRUE;           // new stack flag
273         AliDebug(15, Form("stack : %d, layer : %d, iDec : %d\n",stack,layer,iDet));
274         // Get the digits array
275         AliTRDarrayADC *digits = (AliTRDarrayADC *) digitsManager->GetDigits(iDet);
276         if (fgDataSuppressionLevel==0 || digits->HasData() ) {  // second part is new!! and is for indicating a new event
277           
278           if (digits->HasData()) digits->Expand();
279           
280           Int_t hcwords = 0;
281           
282           // Process A side of the chamber
283           hcwords = ProduceHcData(digits,0,iDet,hcBuffer,kMaxHcWords,newEvent,newSM);
284           if ( newEvent ) newEvent = kFALSE;
285           //AssignLinkMask(hcBuffer, layer);  // active link mask for this layer(2*HC)
286           of->WriteBuffer((char *) hcBuffer, hcwords*4);
287           npayloadbyte += hcwords*4;
288           
289           // Process B side of the chamber
290           hcwords = ProduceHcData(digits,1,iDet,hcBuffer,kMaxHcWords,newEvent,newSM);
291           of->WriteBuffer((char *) hcBuffer, hcwords*4);
292           npayloadbyte += hcwords*4;
293         } // has data
294         
295       } // loop over layer
296     } // loop over stack
297     
298     // Complete header
299     header.fSize = UInt_t(of->Tellp()) - hpos;
300     header.SetAttribute(0);  // Valid data
301     of->Seekp(hpos);         // Rewind to header position
302     of->WriteBuffer((char *) (& header), sizeof(header));
303     delete of;
304
305     delete [] iwbuffer;
306
307   } // loop over sector(SM)
308   
309   delete [] hcBuffer;
310   
311   return kTRUE;
312 }
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     fSMindexPos = nw;       // memorize position of the SM index word for re-allocating stack mask
349     buf[nw++] = 0x0001b020; // SM index words
350     buf[nw++] = 0x10404071; // SM header
351
352     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 }
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     buf[fSMindexPos] = buf[fSMindexPos] | ( 1 << nStack );
368 }
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     if ( ShiftWords(buf, fStackindexPos, 8, nMax)== kFALSE ){
392         AliError("Adding stack header failed.");
393         return 0;
394     }
395
396     buf[fStackindexPos++] = 0x0007a000; // stack index words
397     buf[fStackindexPos++] = 0x04045b01; // stack header
398     for (Int_t i=0;i<6;i++) buf[fStackindexPos++] = 0x00000000; // 6 dummy words 
399     nAddedWords += 8;
400
401     return nAddedWords;
402 }
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     buf[fStackindexPos-8] = buf[fStackindexPos-8] | ( 3 << (2*nLayer) );    // 3 = 0011 
410 }
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     for ( Int_t iw=nMax; iw>nStart-1; iw--){
423         buf[iw+nWords] = buf[iw];
424     }
425     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         Int_t           nw = 0;                       // Number of written    words
437         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         Int_t        layer = fGeo->GetLayer( det );   // Layer
441         Int_t        stack = fGeo->GetStack( det );   // Stack
442         Int_t         sect = fGeo->GetSector( det );  // Sector (=iDDL)
443         const Int_t kCtype = fGeo->GetStack(det) == 2 ? 0 : 1;                       // Chamber type (0:C0, 1:C1)
444
445         Bool_t trackletOn = fFee->GetTracklet();     // tracklet simulation active?
446
447         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         if (trackletOn) {
453           tempBuffer = new UInt_t[maxSize];
454           tempnw     = new Int_t(0);
455           tempof     = new Int_t(0);
456         }
457           
458         WriteIntermediateWords(tempBuffer,*tempnw,*tempof,maxSize,det,side);
459
460         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           for (Int_t iRobRow = 0; iRobRow <= (kCtype + 3)-1; iRobRow++ ) {
464             // ROB number should be increasing
465             Int_t iRob = iRobRow * 2 + side;
466             // MCM on one ROB
467             for (Int_t iMcmRB = 0; iMcmRB < fGeo->MCMmax(); iMcmRB++ ) {
468               Int_t iMcm = 16 - 4*(iMcmRB/4 + 1) + (iMcmRB%4);
469               
470               fMcmSim->Init(det, iRob, iMcm);
471               fMcmSim->SetData(digits);     // no filtering done here (already done in digitizer)
472               if (trackletOn) {
473                 fMcmSim->Tracklet();
474                 Int_t tempNw = fMcmSim->ProduceTrackletStream(&buf[nw], maxSize - nw);
475                 if(  tempNw < 0 ) {
476                   of += tempNw;
477                   nw += maxSize - nw;
478                   AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
479                 } else {
480                   nw += tempNw;
481                 }
482               }
483               fMcmSim->ZSMapping();  // Calculate zero suppression mapping
484               // at the moment it has to be rerun here
485               // Write MCM data to temp. buffer
486               Int_t tempNw = fMcmSim->ProduceRawStream( &tempBuffer[*tempnw], maxSize - *tempnw, fEventCounter );
487               if ( tempNw < 0 ) {
488                 *tempof += tempNw;
489                 *tempnw += maxSize - nw;
490                 AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
491               } else {
492                 *tempnw += tempNw;
493               }
494             }
495           }
496         }
497           
498         // in case of tracklet writing copy temp data to final buffer
499         if (trackletOn) {
500           if (nw + *tempnw < maxSize) {
501             memcpy(&buf[nw], tempBuffer, *tempnw * sizeof(UInt_t));
502             nw += *tempnw;
503           }
504           else {
505             AliError("Buffer overflow detected");
506           }
507
508           delete [] tempBuffer;
509           delete tempof;
510           delete tempnw;
511         }
512     
513         // Write end of raw data marker
514         if (nw+3 < maxSize) {
515           buf[nw++] = fgkEndOfDataMarker;
516           buf[nw++] = fgkEndOfDataMarker;
517           buf[nw++] = fgkEndOfDataMarker;
518           buf[nw++] = fgkEndOfDataMarker;
519         } else {
520           of += 4;
521         }
522         
523         if (of != 0) {
524           AliError("Buffer overflow. Data is truncated. Please increase buffer size and recompile.");
525         }
526
527         return nw;
528 }
529
530 //_____________________________________________________________________________
531 AliTRDdigitsManager *AliTRDrawData::Raw2Digits(AliRawReader *rawReader)
532 {
533   //
534   // Vx of the raw data reading
535   //
536
537   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   AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
547   digitsManager->CreateArrays();
548
549   AliTRDrawStream input(rawReader);
550
551   AliRunLoader *runLoader = AliRunLoader::Instance();
552
553   // ----- preparing tracklet output -----
554   AliDataLoader *trklLoader = runLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets");
555   if (!trklLoader) {
556     trklLoader = new AliDataLoader("TRD.Tracklets.root","tracklets", "tracklets");
557     runLoader->GetLoader("TRDLoader")->AddDataLoader(trklLoader);
558   }
559   if (!trklLoader) {
560     AliError("Could not get the tracklet data loader!");
561     return 0x0;
562   }
563   // trklLoader->SetEvent();
564   trklLoader->Load("update");
565   AliTreeLoader *trklTreeLoader = dynamic_cast<AliTreeLoader*> (trklLoader->GetBaseLoader("tracklets-raw"));
566   if (!trklTreeLoader) {
567     trklTreeLoader = new AliTreeLoader("tracklets-raw", trklLoader);
568     trklLoader->AddBaseLoader(trklTreeLoader);
569   }
570
571   TClonesArray *fTrackletArray = new TClonesArray("AliTRDtrackletWord", 256);
572   Int_t lastHC = -1;
573   if (trklLoader->Tree())
574     trklLoader->MakeTree();
575   TTree *trackletTree = trklTreeLoader->Tree();
576   if (!trackletTree) {
577     trklTreeLoader->MakeTree();
578     trackletTree = trklTreeLoader->Tree();
579   }
580   if (!trackletTree->GetBranch("hc"))
581     trackletTree->Branch("hc", &lastHC, "hc/I");
582   else
583     trackletTree->SetBranchAddress("hc", &lastHC);
584
585   if (!trackletTree->GetBranch("trkl"))
586     trackletTree->Branch("trkl", &fTrackletArray);
587   else
588     trackletTree->SetBranchAddress("trkl", &fTrackletArray);
589
590   // ----- preparing track output -----
591   AliDataLoader *trkLoader = 0x0;
592   trkLoader = runLoader->GetLoader("TRDLoader")->GetDataLoader("gtutracks");
593   if (!trkLoader) {
594     trkLoader = new AliDataLoader("TRD.GtuTracks.root","gtutracks", "gtutracks");
595     runLoader->GetLoader("TRDLoader")->AddDataLoader(trkLoader);
596   }
597   if (!trkLoader) {
598     AliError("Could not get the GTU-track data loader!");
599     return 0x0;
600   }
601   trkLoader->Load("UPDATE");
602
603   TTree *trackTree = trkLoader->Tree();
604   if (!trackTree) {
605     trkLoader->MakeTree();
606     trackTree = trkLoader->Tree();
607   }
608
609   AliTRDtrackGTU trk;
610   if (!trackTree->GetBranch("TRDtrackGTU"))
611     trackTree->Branch("TRDtrackGTU", "AliTRDtrackGTU", &trk, 32000);
612
613   // ----- read the raw data -----
614   // write the tracklets to arrays while reading raw data
615   input.SetTrackletArray(fTracklets);
616   input.SetTrackArray(fTracks);
617
618   // Loop through the digits
619   Int_t det    = 0;
620
621   while (det >= 0)
622     {
623       det = input.NextChamber(digitsManager);
624
625       if (det >= 0)
626         {
627           // get...
628           digits = (AliTRDarrayADC *) digitsManager->GetDigits(det);
629           track0 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,0);
630           track1 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,1);
631           track2 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,2);
632           // and compress
633           if (digits) digits->Compress();  
634           if (track0) track0->Compress();   
635           if (track1) track1->Compress();     
636           if (track2) track2->Compress();
637         }
638     }
639
640   for (Int_t iSector = 0; iSector < fGeo->Nsector(); iSector++) {
641     fTrgFlags[iSector] = input.GetTriggerFlags(iSector);
642   }
643
644   // ----- tracklet output -----
645   fTrackletArray->Clear();
646   Int_t nTracklets = fTracklets ? fTracklets->GetEntries() : 0;
647   AliDebug(1, Form("Writing %i tracklets to loader", nTracklets));
648   for (Int_t iTracklet = 0; iTracklet < nTracklets; ++iTracklet) {
649     AliTRDtrackletWord *trkl = (AliTRDtrackletWord*) fTracklets->At(iTracklet);
650     if (trkl->GetHCId() != lastHC) {
651       if (fTrackletArray->GetEntriesFast() > 0) {
652         trackletTree->Fill();
653         fTrackletArray->Clear();
654       }
655       lastHC = trkl->GetHCId();
656     }
657     new ((*fTrackletArray)[fTrackletArray->GetEntriesFast()]) AliTRDtrackletWord(*trkl);
658   }
659   if (fTrackletArray->GetEntriesFast() > 0) {
660       trackletTree->Fill();
661       fTrackletArray->Clear();
662   }
663   trklTreeLoader->WriteData("OVERWRITE");
664   trklLoader->WriteData("OVERWRITE");
665   trklLoader->UnloadAll();
666   trklLoader->CloseFile();
667
668   // ----- track output -----
669   Int_t nTracks = fTracks ? fTracks->GetEntriesFast() : 0;
670   AliDebug(1, Form("Writing %i tracks to loader", nTracks));
671   for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
672     AliESDTrdTrack *trkEsd = (AliESDTrdTrack*) fTracks->At(iTrack);
673     trk = *trkEsd;
674     trackTree->Fill();
675   }
676   trkLoader->WriteData("OVERWRITE");
677   trkLoader->UnloadAll();
678   trkLoader->CloseFile();
679
680   return digitsManager;
681 }
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     Int_t        layer = fGeo->GetLayer( det );   // Layer
691     Int_t        stack = fGeo->GetStack( det );   // Stack
692     Int_t         sect = fGeo->GetSector( det );  // Sector (=iDDL)
693     Int_t           rv = fFee->GetRAWversion();
694     const Int_t kNTBin = fDigitsParam->GetNTimeBins(det);
695     Bool_t  trackletOn = fFee->GetTracklet();
696     UInt_t           x = 0;
697
698     // Write end of tracklet marker
699     if (nw < maxSize){
700       buf[nw++] = fgkEndOfTrackletMarker;
701       buf[nw++] = fgkEndOfTrackletMarker;     // the number of tracklet end marker should be more than 2
702     }
703     else {
704       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     Int_t zs     = (rv==3) ? 1 : 0;                     // zero suppression
722     Int_t dt     = (trackletOn) ? 0 : 1;        // disable tracklet 
723     
724     majorv = (tp<<6) | (zs<<5) | (dt<<4) | 1;   // major version
725     
726     x = (1<<31) | (majorv<<24) | (minorv<<17) | (add<<14) | (sect<<9) | (layer<<6) | (stack<<3) | (side<<2) | 1;
727     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     x = ((kNTBin)<<26) | (bcCtr<<10) | (ptCtr<<6) | (ptPhase<<2) | 1;
739     if (nw < maxSize) buf[nw++] = x; else of++;
740 }