]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDrawData.cxx
Introduce smaller AliTRDCalDCSFEEv2 object (Frederick)
[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     snprintf(name,1024,"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
299     delete [] iwbuffer;
300
301   } // loop over sector(SM)
302   
303   delete [] hcBuffer;
304   
305   return kTRUE;
306 }
307
308 //_____________________________________________________________________________
309 void AliTRDrawData::ProduceSMIndexData(UInt_t *buf, Int_t& nw){
310         // 
311         // This function generates 
312         //       1) SM index words : ssssssss ssssssss vvvv rrrr r d t mmmmm
313         //          - s : size of SM header (number of header, default = 0x0001)
314         //          - v : SM header version (default = 0xa)
315         //          - r : reserved for future use (default = 00000)
316         //          - d : track data enabled bit (default = 0)
317         //          - t : tracklet data enabled bit (default = 1)
318         //          - m : stack mask (each bit corresponds a stack, default = 11111)
319         //
320         //       2) SM header : rrr c vvvv vvvvvvvv vvvv rrrr bbbbbbbb
321         //          - r : reserved for future use (default = 000)
322         //          - c : clean check out flag (default = 1)
323         //          - v : hardware design revision (default = 0x0404)
324         //          - r : reserved for future use (default = 0x0)
325         //          - b : physical board ID (default = 0x71)
326         //
327         //       3) stack index words : ssssssss ssssssss vvvv mmmm mmmmmmmm
328         //          - s : size of stack header (number of header, (default = 0x0007)
329         //          - v : header version (default = 0xa)
330         //          - m : link mask (default = 0xfff)
331         //
332         //       4) stack header : vvvvvvvv vvvvvvvv bbbbbbbb rrrr rrr c
333         //          - v : hardware design revision (default = 0x0404)
334         //          - b : physical board ID (default = 0x5b)
335         //          - r : reserved for future use (default = 0000 000)
336         //          - c : clean checkout flag (default = 1)
337         //      
338         //       and 6 dummy words(0x00000000)
339         //
340         
341     //buf[nw++] = 0x0001a03f;   // SM index words
342     fSMindexPos = nw;       // memorize position of the SM index word for re-allocating stack mask
343     buf[nw++] = 0x0001a020; // SM index words
344     buf[nw++] = 0x10404071; // SM header
345
346     fStackindexPos = nw;    // memorize position of the stack index word for future adding
347         /*  
348     for (Int_t istack=0; istack<5; istack++){
349         buf[nw++] = 0x0007afff; // stack index words
350         buf[nw++] = 0x04045b01; // stack header
351         for (Int_t i=0;i<6;i++) buf[nw++] = 0x00000000; // 6 dummy words
352     } // loop over 5 stacks
353         */
354 }
355
356 //_____________________________________________________________________________
357 void AliTRDrawData::AssignStackMask(UInt_t *buf, Int_t nStack){
358     //
359     // This function re-assign stack mask active(from 0 to 1) in the SM index word
360     //   
361     buf[fSMindexPos] = buf[fSMindexPos] | ( 1 << nStack );
362 }
363
364 //_____________________________________________________________________________  
365 Int_t AliTRDrawData::AddStackIndexWords(UInt_t *buf, Int_t /*nStack*/, Int_t nMax){
366     // 
367     // This function add stack index words and stack header when there is data for the stack
368     //
369     //   1) stack index words : ssssssss ssssssss vvvv mmmm mmmmmmmm 
370     //      - s : size of stack header (number of header, (default = 0x0007)       
371     //      - v : header version (default = 0xa)
372     //      - m : link mask (default = 0xfff)
373     //      - m : link mask (starting value = 0x000)
374     //
375     //   2) stack header : vvvvvvvv vvvvvvvv bbbbbbbb rrrr rrr c
376     //      - v : hardware design revision (default = 0x0404)
377     //      - b : physical board ID (default = 0x5b)
378     //      - r : reserved for future use (default = 0000 000)
379     //      - c : clean checkout flag (default = 1)
380     //  
381     //   and 6 dummy words(0x00000000)
382     //
383
384     Int_t nAddedWords = 0;  // Number of added words
385     if ( ShiftWords(buf, fStackindexPos, 8, nMax)== kFALSE ){
386         AliError("Adding stack header failed.");
387         return 0;
388     }
389
390     buf[fStackindexPos++] = 0x0007a000; // stack index words
391     buf[fStackindexPos++] = 0x04045b01; // stack header
392     for (Int_t i=0;i<6;i++) buf[fStackindexPos++] = 0x00000000; // 6 dummy words 
393     nAddedWords += 8;
394
395     return nAddedWords;
396 }
397
398 //_____________________________________________________________________________
399 void AliTRDrawData::AssignLinkMask(UInt_t *buf, Int_t nLayer){
400     //
401     // This function re-assign link mask active(from 0 to 1) in the stack index word
402     //   
403     buf[fStackindexPos-8] = buf[fStackindexPos-8] | ( 3 << (2*nLayer) );    // 3 = 0011 
404 }
405
406 //_____________________________________________________________________________ 
407 Bool_t AliTRDrawData::ShiftWords(UInt_t *buf, Int_t nStart, Int_t nWords, Int_t nMax){
408     //  
409     // This function shifts n words
410     //
411     //if ( nStart+nWords > sizeof(buf)/sizeof(UInt_t) ){
412     //  AliError("Words shift failed. No more buffer space.");
413     //  return kFALSE;
414     //}
415
416     for ( Int_t iw=nMax; iw>nStart-1; iw--){
417         buf[iw+nWords] = buf[iw];
418     }
419     return kTRUE;
420 }
421
422 //_____________________________________________________________________________
423 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*/){
424         //
425         // This function produces the raw data for one HC, i.e. tracklets, tracklet endmarkers, 
426         // raw data, raw data endmarkers. 
427         // This function can be used for both ZS and NZS data
428         //
429
430         Int_t           nw = 0;                       // Number of written    words
431         Int_t           of = 0;                       // Number of overflowed words
432         Int_t      *tempnw = &nw;                     // Number of written    words for temp. buffer
433         Int_t      *tempof = &of;                     // Number of overflowed words for temp. buffer
434         Int_t        layer = fGeo->GetLayer( det );   // Layer
435         Int_t        stack = fGeo->GetStack( det );   // Stack
436         Int_t         sect = fGeo->GetSector( det );  // Sector (=iDDL)
437         const Int_t kCtype = fGeo->GetStack(det) == 2 ? 0 : 1;                       // Chamber type (0:C0, 1:C1)
438
439         Bool_t trackletOn = fFee->GetTracklet();     // tracklet simulation active?
440
441         AliDebug(1,Form("Producing raw data for sect=%d layer=%d stack=%d side=%d",sect,layer,stack,side));
442         
443         UInt_t *tempBuffer = buf; // tempBuffer used to write ADC data
444                                   // different in case of tracklet writing
445         
446         if (trackletOn) {
447           tempBuffer = new UInt_t[maxSize];
448           tempnw     = new Int_t(0);
449           tempof     = new Int_t(0);
450         }
451           
452         WriteIntermediateWords(tempBuffer,*tempnw,*tempof,maxSize,det,side);
453
454         if (digits->HasData()) {
455           // scanning direction such, that tracklet-words are sorted in ascending z and then in ascending y order
456           // ROB numbering on chamber and MCM numbering on ROB increase with decreasing z and increasing y
457           for (Int_t iRobRow = 0; iRobRow <= (kCtype + 3)-1; iRobRow++ ) {
458             // ROB number should be increasing
459             Int_t iRob = iRobRow * 2 + side;
460             // MCM on one ROB
461             for (Int_t iMcmRB = 0; iMcmRB < fGeo->MCMmax(); iMcmRB++ ) {
462               Int_t iMcm = 16 - 4*(iMcmRB/4 + 1) + (iMcmRB%4);
463               
464               fMcmSim->Init(det, iRob, iMcm);
465               fMcmSim->SetData(digits);     // no filtering done here (already done in digitizer)
466               if (trackletOn) {
467                 fMcmSim->Tracklet();
468                 Int_t tempNw = fMcmSim->ProduceTrackletStream(&buf[nw], maxSize - nw);
469                 if(  tempNw < 0 ) {
470                   of += tempNw;
471                   nw += maxSize - nw;
472                   AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
473                 } else {
474                   nw += tempNw;
475                 }
476               }
477               fMcmSim->ZSMapping();  // Calculate zero suppression mapping
478               // at the moment it has to be rerun here
479               // Write MCM data to temp. buffer
480               Int_t tempNw = fMcmSim->ProduceRawStream( &tempBuffer[*tempnw], maxSize - *tempnw, fEventCounter );
481               if ( tempNw < 0 ) {
482                 *tempof += tempNw;
483                 *tempnw += maxSize - nw;
484                 AliError(Form("Buffer overflow detected. Please increase the buffer size and recompile."));
485               } else {
486                 *tempnw += tempNw;
487               }
488             }
489           }
490           
491           // in case of tracklet writing copy temp data to final buffer
492           if (trackletOn) {
493             if (nw + *tempnw < maxSize) {
494               memcpy(&buf[nw], tempBuffer, *tempnw * sizeof(UInt_t));
495               nw += *tempnw;
496             }
497             else {
498               AliError("Buffer overflow detected");
499             }
500           }
501         }
502
503         if (trackletOn) {
504           delete [] tempBuffer;
505           delete tempof;
506           delete tempnw;
507         }
508     
509         // Write end of raw data marker
510         if (nw+3 < maxSize) {
511           buf[nw++] = fgkEndOfDataMarker;
512           buf[nw++] = fgkEndOfDataMarker;
513           buf[nw++] = fgkEndOfDataMarker;
514           buf[nw++] = fgkEndOfDataMarker;
515         } else {
516           of += 4;
517         }
518         
519         if (of != 0) {
520           AliError("Buffer overflow. Data is truncated. Please increase buffer size and recompile.");
521         }
522
523         return nw;
524 }
525
526 //_____________________________________________________________________________
527 AliTRDdigitsManager *AliTRDrawData::Raw2Digits(AliRawReader *rawReader)
528 {
529   //
530   // Vx of the raw data reading
531   //
532
533   rawReader->Select("TRD"); //[mj]
534
535   AliTRDarrayADC *digits = 0;
536   AliTRDarrayDictionary *track0 = 0;
537   AliTRDarrayDictionary *track1 = 0;
538   AliTRDarrayDictionary *track2 = 0;  
539
540   //AliTRDSignalIndex *indexes = 0;
541   // Create the digits manager
542   AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
543   digitsManager->CreateArrays();
544
545   if (!fTrackletContainer) {
546     const Int_t kTrackletChmb=256;
547     fTrackletContainer = new UInt_t *[2];
548     fTrackletContainer[0] = new UInt_t[kTrackletChmb];
549     fTrackletContainer[1] = new UInt_t[kTrackletChmb];
550     memset(fTrackletContainer[0], 0, kTrackletChmb*sizeof(UInt_t)); //jkl
551     memset(fTrackletContainer[1], 0, kTrackletChmb*sizeof(UInt_t)); //jkl
552   }
553
554   AliTRDrawStreamBase *pinput = AliTRDrawStreamBase::GetRawStream(rawReader);
555   AliTRDrawStreamBase &input = *pinput;
556   input.SetRawVersion( fFee->GetRAWversion() ); //<= ADDED by MinJung
557
558   AliInfo(Form("Stream version: %s", input.IsA()->GetName()));
559
560   // ----- preparing tracklet output -----
561   AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets");
562   if (!trklLoader) {
563     //AliInfo("Could not get the tracklets data loader, adding it now!");
564     trklLoader = new AliDataLoader("TRD.Tracklets.root","tracklets", "tracklets");
565     AliRunLoader::Instance()->GetLoader("TRDLoader")->AddDataLoader(trklLoader);
566   }
567   if (!trklLoader) {
568     return 0x0;
569   }
570   AliTreeLoader *trklTreeLoader = dynamic_cast<AliTreeLoader*> (trklLoader->GetBaseLoader("tracklets-raw"));
571   if (!trklTreeLoader) {
572     trklTreeLoader = new AliTreeLoader("tracklets-raw", trklLoader);
573     trklLoader->AddBaseLoader(trklTreeLoader);
574   }
575   if (!trklTreeLoader) {
576     return 0x0;
577   }
578
579   if (!trklTreeLoader->Tree())
580     trklTreeLoader->MakeTree();
581
582   // Loop through the digits
583   Int_t det    = 0;
584
585   while (det >= 0)
586     {
587       det = input.NextChamber(digitsManager,fTrackletContainer);
588
589       if (*(fTrackletContainer[0]) > 0 || *(fTrackletContainer[1]) > 0) WriteTracklets(det);
590
591       if (det >= 0)
592         {
593           // get...
594           digits = (AliTRDarrayADC *) digitsManager->GetDigits(det);
595           track0 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,0);
596           track1 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,1);
597           track2 = (AliTRDarrayDictionary *) digitsManager->GetDictionary(det,2);
598           // and compress
599           if (digits) digits->Compress();  
600           if (track0) track0->Compress();   
601           if (track1) track1->Compress();     
602           if (track2) track2->Compress();
603         }
604     }
605
606   trklTreeLoader->WriteData("OVERWRITE");
607   trklLoader->UnloadAll();
608
609   if (fTrackletContainer){
610     delete [] fTrackletContainer[0];
611     delete [] fTrackletContainer[1];
612     delete [] fTrackletContainer;
613     fTrackletContainer = NULL;
614   }
615
616   delete pinput;
617   pinput = NULL;
618
619   return digitsManager;
620
621 }
622
623 //_____________________________________________________________________________
624 void AliTRDrawData::WriteIntermediateWords(UInt_t* buf, Int_t& nw, Int_t& of, const Int_t& maxSize, const Int_t& det, const Int_t& side) {
625     // 
626     // write tracklet end marker(0x10001000) 
627     // and half chamber headers(H[0] and H[1])
628     //
629     
630     Int_t        layer = fGeo->GetLayer( det );   // Layer
631     Int_t        stack = fGeo->GetStack( det );   // Stack
632     Int_t         sect = fGeo->GetSector( det );  // Sector (=iDDL)
633     Int_t           rv = fFee->GetRAWversion();
634     const Int_t kNTBin = fDigitsParam->GetNTimeBins(det);
635     Bool_t  trackletOn = fFee->GetTracklet();
636     UInt_t           x = 0;
637
638     // Write end of tracklet marker
639     if (nw < maxSize){
640       buf[nw++] = fgkEndOfTrackletMarker;
641       buf[nw++] = fgkEndOfTrackletMarker;     // the number of tracklet end marker should be more than 2
642     }
643     else {
644       of++;
645     }
646    
647     // Half Chamber header
648     // h[0] (there are 2 HC headers) xmmm mmmm nnnn nnnq qqss sssp ppcc ci01
649     //  , where  x : Raw version speacial number (=1)
650     //               m : Raw version major number (test pattern, ZS, disable tracklet, 0, options)
651     //               n : Raw version minor number
652     //               q : number of addtional header words (default = 1)
653     //               s : SM sector number (ALICE numbering)
654     //               p : plane(layer) number
655     //                   c : chamber(stack) number
656     //                   i : side number (0:A, 1:B)
657     Int_t majorv = 0;   // The major version number 
658     Int_t minorv = 0;   // The minor version number
659     Int_t add    = 1;   // The number of additional header words to follow : now 1, previous 2
660     Int_t tp     = 0;   // test pattern (default=0)
661     Int_t zs     = (rv==3) ? 1 : 0;                     // zero suppression
662     Int_t dt     = (trackletOn) ? 0 : 1;        // disable tracklet 
663     
664     majorv = (tp<<6) | (zs<<5) | (dt<<4) | 1;   // major version
665     
666     x = (1<<31) | (majorv<<24) | (minorv<<17) | (add<<14) | (sect<<9) | (layer<<6) | (stack<<3) | (side<<2) | 1;
667     if (nw < maxSize) buf[nw++] = x; else of++;
668     
669     // h[1]             tttt ttbb bbbb bbbb bbbb bbpp pphh hh01
670     // , where  t : number of time bins
671     //          b : bunch crossing number
672     //          p : pretrigger counter
673     //          h : pretrigger phase
674     Int_t bcCtr   = 99; // bunch crossing counter. Here it is set to 99 always for no reason
675     Int_t ptCtr   = 15; // pretrigger counter. Here it is set to 15 always for no reason
676     Int_t ptPhase = 11; // pretrigger phase. Here it is set to 11 always for no reason
677     //x = (bcCtr<<16) | (ptCtr<<12) | (ptPhase<<8) | ((kNTBin-1)<<2) | 1;       // old format
678     x = ((kNTBin)<<26) | (bcCtr<<10) | (ptCtr<<6) | (ptPhase<<2) | 1;
679     if (nw < maxSize) buf[nw++] = x; else of++;
680 }
681
682 //_____________________________________________________________________________
683 Bool_t AliTRDrawData::WriteTracklets(Int_t det)
684 {
685   //
686   // Write the raw data tracklets into seperate file
687   //
688
689   UInt_t **leaves = new UInt_t *[2];
690   for (Int_t i=0; i<2 ;i++){
691     leaves[i] = new UInt_t[258];
692     leaves[i][0] = det; // det
693     leaves[i][1] = i;   // side
694     memcpy(leaves[i]+2, fTrackletContainer[i], sizeof(UInt_t) * 256);
695   }
696
697   if (!fTrackletTree){
698     AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets");
699     dl->MakeTree();
700     fTrackletTree = dl->Tree();
701   }
702
703   TBranch *trkbranch = fTrackletTree->GetBranch("trkbranch");
704   if (!trkbranch) {
705     trkbranch = fTrackletTree->Branch("trkbranch",leaves[0],"det/i:side/i:tracklets[256]/i");
706   }
707
708   for (Int_t i=0; i<2; i++){
709     if (leaves[i][2]>0) {
710       trkbranch->SetAddress(leaves[i]);
711       fTrackletTree->Fill();
712     }
713   }
714
715   //  AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets"); //jkl: wrong
716   //  dl->WriteData("OVERWRITE"); //jkl: wrong
717   //dl->Unload();
718   delete [] leaves;
719
720   return kTRUE;
721 }