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