]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFDDLRawData.cxx
updates to comply with AliTOFGeometryV5 becoming AliTOFGeometry
[u/mrichter/AliRoot.git] / TOF / AliTOFDDLRawData.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.19  2007/06/22 11:37:47  cvetan
19 Fixes in order to write correct raw-data on big-endian platforms (Marco)
20
21 Revision 1.18  2007/05/21 13:26:19  decaro
22 Correction on matching_window control and bug fixing (R.Preghenella)
23
24 Revision 1.17  2007/05/10 09:29:34  hristov
25 Last moment fixes and changes from v4-05-Release (Silvia)
26
27 Revision 1.16  2007/05/03 09:07:22  decaro
28 Double digit in the same TDC channel. Wrong sequence during the raw data writing in unpacked mode: solved
29
30 Revision 1.15  2007/04/23 16:51:39  decaro
31 Digits-to-raw_data conversion: correction for a more real description (A.De Caro, R.Preghenella)
32
33 Revision 1.14  2007/03/28 10:50:33  decaro
34 Rounding off problem in rawData coding/decoding: solved
35
36 Revision 1.13  2007/02/20 15:57:00  decaro
37 Raw data update: to read the TOF raw data defined in UNPACKED mode
38
39 Revision 1.12  2006/08/22 13:29:42  arcelli
40 removal of effective c++ warnings (C.Zampolli)
41
42 Revision 1.11  2006/08/10 14:46:54  decaro
43 TOF raw data format: updated version
44
45 Revision 1.10.1  2006/06/28 A.De Caro
46         Update TOF raw data format
47                according to the final version
48                (see ALICE internal note in preparation
49                 'ALICE TOF raw data format')
50
51 Revision 0.02  2005/7/25 A.De Caro
52         Update number of bits allocated for time-of-flight
53                and 'charge' measurements
54
55 Revision 0.01  2004/6/11 A.De Caro, S.B.Sellitto, R.Silvestri
56         First implementation: global methods RawDataTOF
57                                              GetDigits
58 */
59
60 ////////////////////////////////////////////////////////////////////
61 //                                                                //
62 // This class contains the methods to create the Raw Data files   //
63 // for the TOF detector starting from the Digits.                 //
64 // In this implementation, we defined the structure               //
65 // of the ALICE-TOF raw data (according to the                    //
66 // ALICE technical note, in preparation)                          //
67 // starting from the TOF digit format.                            //
68 //                                                                //
69 ////////////////////////////////////////////////////////////////////
70
71 #include "Riostream.h"
72
73 #include "TBranch.h"
74 #include "TClonesArray.h"
75 #include "TMath.h"
76 #include "TRandom.h"
77
78 #include "AliBitPacking.h"
79 #include "AliDAQ.h"
80 #include "AliLog.h"
81 //#include "AliRawDataHeader.h"
82 #include "AliRawDataHeaderSim.h"
83
84 #include "AliTOFDDLRawData.h"
85 #include "AliTOFDigitMap.h"
86 #include "AliTOFdigit.h"
87 #include "AliTOFGeometry.h"
88 #include "AliTOFRawStream.h"
89 #include "AliFstream.h"
90
91 extern TRandom *gRandom;
92
93 ClassImp(AliTOFDDLRawData)
94
95 //---------------------------------------------------------------------------
96 AliTOFDDLRawData::AliTOFDDLRawData():
97   fVerbose(0),
98   fIndex(-1),
99   fPackedAcquisition(kTRUE),
100   fFakeOrphaneProduction(kFALSE),
101   fMatchingWindow(8192),
102   fTOFdigitMap(new AliTOFDigitMap()),
103   fTOFdigitArray(0x0),
104   fTOFrawStream(new AliTOFRawStream())
105 {
106   //Default constructor
107 }
108 //----------------------------------------------------------------------------
109 AliTOFDDLRawData::AliTOFDDLRawData(const AliTOFDDLRawData &source) :
110   TObject(source),
111   fVerbose(0),
112   fIndex(-1),
113   fPackedAcquisition(kTRUE),
114   fFakeOrphaneProduction(kFALSE),
115   fMatchingWindow(8192),
116   fTOFdigitMap(new AliTOFDigitMap()),
117   fTOFdigitArray(0x0),
118   fTOFrawStream(new AliTOFRawStream())
119  {
120   //Copy Constructor
121   this->fIndex=source.fIndex;
122   this->fVerbose=source.fVerbose;
123   this->fPackedAcquisition=source.fPackedAcquisition;
124   this->fFakeOrphaneProduction=source.fFakeOrphaneProduction;
125   this->fMatchingWindow=source.fMatchingWindow;
126   this->fTOFdigitMap=source.fTOFdigitMap;
127   this->fTOFdigitArray=source.fTOFdigitArray;
128   this->fTOFrawStream=source.fTOFrawStream;
129   return;
130 }
131
132 //----------------------------------------------------------------------------
133 AliTOFDDLRawData& AliTOFDDLRawData::operator=(const AliTOFDDLRawData &source) {
134   //Assigment operator
135   this->fIndex=source.fIndex;
136   this->fVerbose=source.fVerbose;
137   this->fPackedAcquisition=source.fPackedAcquisition;
138   this->fFakeOrphaneProduction=source.fFakeOrphaneProduction;
139   this->fMatchingWindow=source.fMatchingWindow;
140   this->fTOFdigitMap=source.fTOFdigitMap;
141   this->fTOFdigitArray=source.fTOFdigitArray;
142   this->fTOFrawStream=source.fTOFrawStream;
143   return *this;
144 }
145
146 //----------------------------------------------------------------------------
147 AliTOFDDLRawData::~AliTOFDDLRawData()
148 {
149   delete fTOFdigitMap;
150   delete fTOFrawStream;
151 }
152 //----------------------------------------------------------------------------
153 Int_t AliTOFDDLRawData::RawDataTOF(TBranch* branch)
154 {
155   //
156   // This method creates the Raw data files for TOF detector
157   //
158
159   const Int_t kSize = 5000; //max number of digits per DDL file times 2
160
161   UInt_t buf[kSize];
162
163   fIndex = -1;
164
165   fTOFdigitArray = * (TClonesArray**) branch->GetAddress();
166
167   char fileName[15];
168   AliFstream* outfile;         // logical name of the output file 
169
170   //AliRawDataHeader header;
171   AliRawDataHeaderSim header;
172
173   UInt_t sizeRawData = 0;
174
175   branch->GetEvent();
176   
177   GetDigits();
178
179   Int_t jj = -1;
180   Int_t iDDL = -1;
181   Int_t nDDL = -1;
182   Int_t nTRM =  0;
183   Int_t iChain = -1;
184
185   UInt_t nWordsPerTRM = 0;
186
187   //loop over TOF DDL files
188   for (nDDL=0; nDDL<AliDAQ::NumberOfDdls("TOF"); nDDL++) {
189
190     strcpy(fileName,AliDAQ::DdlFileName("TOF",nDDL)); //The name of the output file
191
192     outfile = new AliFstream(fileName);
193     iDDL = fTOFrawStream->GetDDLnumberPerSector(nDDL);
194
195     // write Dummy DATA HEADER
196     UInt_t dataHeaderPosition = outfile->Tellp();
197     outfile->WriteBuffer((char*)(&header),sizeof(header));
198
199     // DRM section: trailer
200     MakeDRMtrailer(buf);
201
202     // LTM section
203     fIndex++;
204     buf[fIndex] = MakeFiller();
205     MakeLTMtrailer(buf);
206     MakeLTMdata(buf);
207     MakeLTMheader(buf);
208
209     // loop on TRM number per DRM
210     for (nTRM=AliTOFGeometry::NTRM(); nTRM>=3; nTRM--) {
211
212       nWordsPerTRM = 0;
213
214       // the slot number 3 of the even DRM (i.e. right) doesn't contain TDC digit data
215       if (iDDL%2==0 && nTRM==3) continue;
216
217       // TRM global trailer
218       MakeTRMtrailer(buf);
219
220       // loop on TRM chain number per TRM
221       for (iChain=AliTOFGeometry::NChain()-1; iChain>=0; iChain--) {
222
223         // TRM chain trailer
224         MakeTRMchainTrailer(iChain, buf);
225         nWordsPerTRM++;
226
227         // TRM TDC digits
228         MakeTDCdigits(nDDL, nTRM, iChain, buf, nWordsPerTRM);
229
230         // TRM chain header
231         MakeTRMchainHeader(nTRM, iChain, buf);
232         nWordsPerTRM++;
233
234       } // end loop on iChain
235
236       // TRM global header
237       MakeTRMheader(nTRM, buf);
238
239       // TRM filler in case where TRM data number is odd
240       if (nWordsPerTRM%2!=0) MakeTRMfiller(buf, nWordsPerTRM);
241
242     } // end loop on nTRM
243        
244     // DRM section: in
245     MakeDRMheader(nDDL, buf);
246
247     ReverseArray(buf, fIndex+1);
248
249     outfile->WriteBuffer((char *)buf,((fIndex+1)*sizeof(UInt_t)));
250
251     for (jj=0; jj<(fIndex+1); jj++) buf[jj]=0;
252     fIndex = -1;
253     
254     //Write REAL DATA HEADER
255     UInt_t currentFilePosition = outfile->Tellp();
256     sizeRawData = currentFilePosition - dataHeaderPosition - sizeof(header);
257     header.fSize = currentFilePosition - dataHeaderPosition;
258     header.SetAttribute(0);  // valid data
259     outfile->Seekp(dataHeaderPosition);
260     outfile->WriteBuffer((char*)(&header),sizeof(header));
261     outfile->Seekp(currentFilePosition);
262
263     delete outfile;
264
265   } //end loop on DDL file number
266
267   return 0;
268
269 }
270
271 //----------------------------------------------------------------------------
272 void AliTOFDDLRawData::GetDigits()
273 {
274   //
275   // Fill the TOF volumes' map with the TOF digit indices
276   //
277
278   Int_t vol[5] = {-1,-1,-1,-1,-1};
279
280   Int_t digit = -1;
281   Int_t ndigits = fTOFdigitArray->GetEntries();
282
283   AliTOFdigit *digs;
284
285   // loop on TOF digits
286   for (digit=0; digit<ndigits; digit++) {
287     digs = (AliTOFdigit*)fTOFdigitArray->UncheckedAt(digit);
288
289     vol[0] = digs->GetSector(); // Sector Number (0-17)
290     vol[1] = digs->GetPlate();  // Plate Number (0-4)
291     vol[2] = digs->GetStrip();  // Strip Number (0-14/18)
292     vol[3] = digs->GetPadx();   // Pad Number in x direction (0-47)
293     vol[4] = digs->GetPadz();   // Pad Number in z direction (0-1)
294
295     fTOFdigitMap->AddDigit(vol, digit);
296
297   } // close loop on digit del TOF
298
299 }
300
301 //----------------------------------------------------------------------------
302 void AliTOFDDLRawData::MakeDRMheader(Int_t nDDL, UInt_t *buf)
303 {
304   //
305   // DRM global header
306   //
307
308   Int_t iDDL = fTOFrawStream->GetDDLnumberPerSector(nDDL);
309
310   Int_t iSector = fTOFrawStream->GetSectorNumber(nDDL);
311
312   UInt_t baseWord=0;
313   UInt_t word;
314
315   // DRM event CRC
316   baseWord=0;
317   word = 1; // 0001 -> DRM data are coming from the VME slot number 1
318   AliBitPacking::PackWord(word,baseWord, 0, 3);
319   word = 0; // event CRC
320   AliBitPacking::PackWord(word,baseWord, 4,19);
321   word = 0; // reserved for future use
322   AliBitPacking::PackWord(word,baseWord,20,27);
323   word = 4; // 0100 -> DRM header ID
324   AliBitPacking::PackWord(word,baseWord,28,31);
325   fIndex++;
326   buf[fIndex]=baseWord;
327
328   // DRM status header 3
329   baseWord=0;
330   word = 1; // 0001 -> DRM data are coming from the VME slot number 1
331   AliBitPacking::PackWord(word,baseWord, 0, 3);
332   word = 0; // TTC event counter
333   AliBitPacking::PackWord(word,baseWord, 4,27);
334   word = 4; // 0100 -> DRM header ID
335   AliBitPacking::PackWord(word,baseWord,28,31);
336   fIndex++;
337   buf[fIndex]=baseWord;
338
339   // DRM status header 2
340   baseWord=0;
341   word = 1; // 0001 -> DRM data are coming from the VME slot number 1
342   AliBitPacking::PackWord(word,baseWord, 0, 3);
343
344   if (iDDL%2==1) {
345     word = 2047; // enable ID: [00000000000;11111111111] for odd
346                  // (i.e. right) crates
347     AliBitPacking::PackWord(word,baseWord, 4,14);
348   } else {
349     word = 2045; // enable ID: [00000000000;11111111101] for even
350                  // (i.e. left) crates
351     AliBitPacking::PackWord(word,baseWord, 4,14);
352   }
353
354   word = 0; //
355   AliBitPacking::PackWord(word,baseWord,15,15);
356   word = 0; // fault ID
357   AliBitPacking::PackWord(word,baseWord,16,27);
358   word = 4; // 0100 -> DRM header ID
359   AliBitPacking::PackWord(word,baseWord,28,31);
360   fIndex++;
361   buf[fIndex]=baseWord;
362   
363   // DRM status header 1
364   baseWord=0;
365   word = 1; // 0001 -> DRM data are coming from the VME slot number 1
366   AliBitPacking::PackWord(word,baseWord, 0, 3);
367
368   if (iDDL%2==1) {
369     word = 2047; // slot ID: [00000000000;11111111111] for odd
370                  // (i.e. right) crates
371     AliBitPacking::PackWord(word,baseWord, 4,14);
372   } else {
373     word = 2045; // slot ID: [00000000000;11111111101] for even
374                  // (i.e. left) crates
375     AliBitPacking::PackWord(word,baseWord, 4,14);
376   }
377       
378   word = 1; // LHC clock status: 1/0
379   AliBitPacking::PackWord(word,baseWord,15,15);
380   word = 0; // reserved for future use
381   AliBitPacking::PackWord(word,baseWord,16,27);
382   word = 4; // 0100 -> DRM header ID
383   AliBitPacking::PackWord(word,baseWord,28,31);
384   fIndex++;
385   buf[fIndex]=baseWord;
386
387   // DRM global header
388   baseWord=0;
389   word = 1; // 0001 -> DRM data are coming from the VME slot number 1
390   AliBitPacking::PackWord(word,baseWord, 0, 3);
391   word = fIndex+1 + 1; // event words
392   AliBitPacking::PackWord(word,baseWord, 4,20);
393   word = iDDL; // crate ID [0;3]
394   AliBitPacking::PackWord(word,baseWord,21,22);
395   word = iSector; // sector ID [0;17]
396   AliBitPacking::PackWord(word,baseWord,23,27);
397   word = 4; // 0100 -> DRM header ID
398   AliBitPacking::PackWord(word,baseWord,28,31);
399   fIndex++;
400   buf[fIndex]=baseWord;
401
402 }
403
404 //----------------------------------------------------------------------------
405 void AliTOFDDLRawData::MakeDRMtrailer(UInt_t *buf)
406 {
407   //
408   // DRM global trailer
409   //
410   
411   UInt_t baseWord;
412   UInt_t word;
413   
414   baseWord=0;
415   word = 1; // 0001 -> DRM data are coming from the VME slot number 1
416   AliBitPacking::PackWord(word,baseWord, 0, 3);
417   word = 0; // local event counter
418   AliBitPacking::PackWord(word,baseWord, 4,15);
419   word = 0; // reserved for future use
420   AliBitPacking::PackWord(word,baseWord,16,27);
421   word = 5; // 0101 -> DRM trailer ID
422   AliBitPacking::PackWord(word,baseWord,28,31);
423   fIndex++;
424   buf[fIndex]=baseWord;
425
426 }
427
428 //----------------------------------------------------------------------------
429 void AliTOFDDLRawData::MakeLTMheader(UInt_t *buf)
430 {
431   //
432   // LTM header
433   //
434
435   UInt_t baseWord;
436   UInt_t word;
437   
438   baseWord=0;
439   word = 2; // 0010 -> LTM data are coming from the VME slot number 2
440   AliBitPacking::PackWord(word,baseWord, 0, 3);
441   word = 35; // event words
442   AliBitPacking::PackWord(word,baseWord, 4,16);
443   word = 0; // crc error
444   AliBitPacking::PackWord(word,baseWord,17,17);
445   word = 0; // fault
446   AliBitPacking::PackWord(word,baseWord,18,23);
447   word = 0;
448   AliBitPacking::PackWord(word,baseWord,24,27);
449   word = 4; // 0100 -> LTM header ID
450   AliBitPacking::PackWord(word,baseWord,28,31);
451   fIndex++;
452   buf[fIndex]=baseWord;
453
454 }
455
456 //----------------------------------------------------------------------------
457 void AliTOFDDLRawData::MakeLTMdata(UInt_t *buf)
458 {
459   //
460   // LTM data
461   //
462
463   UInt_t baseWord;
464   UInt_t word;
465   
466   baseWord=0;
467   word = 100; // Local temperature in LTM5 -> 4 X 25 degree (environment temperature)
468   AliBitPacking::PackWord(word,baseWord, 0, 9);
469   word = 100; // Local temperature in LTM6 -> 4 X 25 degree (environment temperature)
470   AliBitPacking::PackWord(word,baseWord,10,19);
471   word = 100; // Local temperature in LTM7 -> 4 X 25 degree (environment temperature)
472   AliBitPacking::PackWord(word,baseWord,20,29);
473   word = 0;
474   AliBitPacking::PackWord(word,baseWord,30,31);
475
476   fIndex++;
477   buf[fIndex]=baseWord;
478
479   // Local temperature in LTM2, LMT3, LTM4 -> 4 X 25 degree (environment temperature)
480   fIndex++;
481   buf[fIndex]=baseWord;
482
483   // Local temperature in FEAC7, FEAC8, LTM1 -> 4 X 25 degree (environment temperature)
484   fIndex++;
485   buf[fIndex]=baseWord;
486
487   // Local temperature in FEAC4, FEAC5, FEAC6 -> 4 X 25 degree (environment temperature)
488   fIndex++;
489   buf[fIndex]=baseWord;
490
491   // Local temperature in FEAC1, FEAC2, FEAC3 -> 4 X 25 degree (environment temperature)
492   fIndex++;
493   buf[fIndex]=baseWord;
494
495   baseWord=0;
496   word = 0; // GND-FEAC15 -> Voltage drop between GND and FEAC15
497   AliBitPacking::PackWord(word,baseWord, 0, 9);
498   word = 0; // VTH16 -> Thereshould voltage for FEAC16
499   AliBitPacking::PackWord(word,baseWord,10,19);
500   word = 0; // GND-FEAC16 -> Voltage drop between GND and FEAC16
501   AliBitPacking::PackWord(word,baseWord,20,29);
502   word = 0;
503   AliBitPacking::PackWord(word,baseWord,30,31);
504
505   fIndex++;
506   buf[fIndex]=baseWord;
507
508   // VTH14 -> Thereshould voltage for FEAC14
509   // GND-FEAC14 -> Voltage drop between GND and FEAC14
510   // VTH15 -> Thereshould voltage for FEAC15
511   fIndex++;
512   buf[fIndex]=baseWord;
513
514   // GND-FEAC12 -> Voltage drop between GND and FEAC12
515   // VTH13 -> Thereshould voltage for FEAC13
516   // GND-FEAC13 -> Voltage drop between GND and FEAC13
517   fIndex++;
518   buf[fIndex]=baseWord;
519
520   // VTH11 -> Thereshould voltage for FEAC11
521   // GND-FEAC11 -> Voltage drop between GND and FEAC11
522   // VTH12 -> Thereshould voltage for FEAC12
523   fIndex++;
524   buf[fIndex]=baseWord;
525
526   // GND-FEAC9 -> Voltage drop between GND and FEAC9
527   // VTH10 -> Thereshould voltage for FEAC10
528   // GND-FEAC10 -> Voltage drop between GND and FEAC10
529   fIndex++;
530   buf[fIndex]=baseWord;
531
532   // VTH8 -> Thereshould voltage for FEAC8
533   // GND-FEAC8 -> Voltage drop between GND and FEAC8
534   // VTH9 -> Thereshould voltage for FEAC9
535   fIndex++;
536   buf[fIndex]=baseWord;
537
538   // GND-FEAC6 -> Voltage drop between GND and FEAC6
539   // VTH7 -> Thereshould voltage for FEAC7
540   // GND-FEAC7 -> Voltage drop between GND and FEAC7
541   fIndex++;
542   buf[fIndex]=baseWord;
543
544   // VTH5 -> Thereshould voltage for FEAC5
545   // GND-FEAC5 -> Voltage drop between GND and FEAC5
546   // VTH6 -> Thereshould voltage for FEAC6
547   fIndex++;
548   buf[fIndex]=baseWord;
549
550   // GND-FEAC3 -> Voltage drop between GND and FEAC3
551   // VTH4 -> Thereshould voltage for FEAC4
552   // GND-FEAC4 -> Voltage drop between GND and FEAC4
553   fIndex++;
554   buf[fIndex]=baseWord;
555
556   // VTH2 -> Thereshould voltage for FEAC2
557   // GND-FEAC2 -> Voltage drop between GND and FEAC2
558   // VTH3 -> Thereshould voltage for FEAC3
559   fIndex++;
560   buf[fIndex]=baseWord;
561
562   // LV16 -> Low Voltage measured by FEAC16
563   // GND-FEAC1 -> Voltage drop between GND and FEAC1
564   // VTH1 -> Thereshould voltage for FEAC1
565   fIndex++;
566   buf[fIndex]=baseWord;
567
568   // Low Voltage measured by FEAC13, FEAC14, FEAC15
569   fIndex++;
570   buf[fIndex]=baseWord;
571
572   // Low Voltage measured by FEAC10, FEAC11, FEAC12
573   fIndex++;
574   buf[fIndex]=baseWord;
575
576   // Low Voltage measured by FEAC7, FEAC8, FEAC9
577   fIndex++;
578   buf[fIndex]=baseWord;
579
580   // Low Voltage measured by FEAC4, FEAC5, FEAC6
581   fIndex++;
582   buf[fIndex]=baseWord;
583
584   // Low Voltage measured by FEAC1, FEAC2, FEAC3
585   fIndex++;
586   buf[fIndex]=baseWord;
587
588
589   baseWord=0;
590   word = 0; // PDL45 -> Delay Line setting for PDL45
591   AliBitPacking::PackWord(word,baseWord, 0, 7);
592   word = 0; // PDL46 -> Delay Line setting for PDL46
593   AliBitPacking::PackWord(word,baseWord, 8,15);
594   word = 0; // PDL47 -> Delay Line setting for PDL47
595   AliBitPacking::PackWord(word,baseWord,16,23);
596   word = 0; // PDL48 -> Delay Line setting for PDL48
597   AliBitPacking::PackWord(word,baseWord,24,31);
598   fIndex++;
599   buf[fIndex]=baseWord;
600
601   // Delay Line setting for PDL41, PDL42, PDL43, PDL44
602   fIndex++;
603   buf[fIndex]=baseWord;
604
605   // Delay Line setting for PDL37, PDL38, PDL39, PDL40
606   fIndex++;
607   buf[fIndex]=baseWord;
608
609   // Delay Line setting for PDL33, PDL34, PDL35, PDL36
610   fIndex++;
611   buf[fIndex]=baseWord;
612
613   // Delay Line setting for PDL29, PDL30, PDL31, PDL32
614   fIndex++;
615   buf[fIndex]=baseWord;
616
617   // Delay Line setting for PDL25, PDL26, PDL27, PDL28
618   fIndex++;
619   buf[fIndex]=baseWord;
620
621   // Delay Line setting for PDL21, PDL22, PDL23, PDL24
622   fIndex++;
623   buf[fIndex]=baseWord;
624
625   // Delay Line setting for PDL17, PDL18, PDL19, PDL20
626   fIndex++;
627   buf[fIndex]=baseWord;
628
629   // Delay Line setting for PDL13, PDL14, PDL15, PDL16
630   fIndex++;
631   buf[fIndex]=baseWord;
632
633   // Delay Line setting for PDL9, PDL10, PDL11, PDL12
634   fIndex++;
635   buf[fIndex]=baseWord;
636
637   // Delay Line setting for PDL5, PDL6, PDL7, PDL8
638   fIndex++;
639   buf[fIndex]=baseWord;
640
641   // Delay Line setting for PDL1, PDL2, PDL3, PDL4
642   fIndex++;
643   buf[fIndex]=baseWord;
644
645 }
646
647 //----------------------------------------------------------------------------
648 void AliTOFDDLRawData::MakeLTMtrailer(UInt_t *buf)
649 {
650   //
651   // LTM trailer
652   //
653  
654   UInt_t baseWord;
655   UInt_t word;
656   
657   baseWord=0;
658   word = 2; // 0010 -> LTM data are coming from the VME slot number 2
659   AliBitPacking::PackWord(word,baseWord, 0, 3);
660   word = 0; // event crc
661   AliBitPacking::PackWord(word,baseWord, 4,15);
662   word = 0; // event number
663   AliBitPacking::PackWord(word,baseWord,16,27);
664   word = 5; // 0101 -> LTM trailer ID
665   AliBitPacking::PackWord(word,baseWord,28,31);
666   fIndex++;
667   buf[fIndex]=baseWord;
668
669 }
670
671 //----------------------------------------------------------------------------
672 void AliTOFDDLRawData::MakeTRMheader(Int_t nTRM, UInt_t *buf)
673 {
674   //
675   // TRM header for the TRM number nTRM [ 3;12]
676   //
677
678   if (nTRM<3 || nTRM>12) {
679     AliWarning(Form(" TRM number is out of the right range [3;12] (nTRM = %3i",nTRM));
680     return;
681   }
682
683   UInt_t baseWord;
684   UInt_t word;
685
686   baseWord = 0;
687   word = nTRM; // TRM data coming from the VME slot number nTRM
688   AliBitPacking::PackWord(word,baseWord, 0, 3);
689   word = 0; // event words
690   AliBitPacking::PackWord(word,baseWord, 4,16);
691
692   if (fPackedAcquisition)
693     word = 0; // ACQuisition mode: [0;3] see document
694   else
695     word = 3; // ACQuisition mode: [0;3] see document
696   AliBitPacking::PackWord(word,baseWord,17,18);
697   word = 0; // description of a SEU inside LUT tables for INL compensation;
698             // the data are unaffected
699   AliBitPacking::PackWord(word,baseWord,19,19);
700   word = 0; // Must Be Zero (MBZ)
701   AliBitPacking::PackWord(word,baseWord,20,27);
702   word = 4; // 0100 -> TRM header ID
703   AliBitPacking::PackWord(word,baseWord,28,31);
704   fIndex++;
705   buf[fIndex]=baseWord;
706
707 }
708
709 //----------------------------------------------------------------------------
710 void AliTOFDDLRawData::MakeTRMtrailer(UInt_t *buf)
711 {
712   //
713   // TRM trailer
714   //
715
716   UInt_t baseWord;
717   UInt_t word;
718
719   baseWord=0;
720   word = 15; // 1111 -> TRM trailer ID 1
721   AliBitPacking::PackWord(word,baseWord, 0, 3);
722   word = 0; // event CRC
723   AliBitPacking::PackWord(word,baseWord, 4,15);
724   word = 0; // local event counter == DRM local event counter
725   AliBitPacking::PackWord(word,baseWord,16,27);
726   word = 5; // 0101 -> TRM trailer ID 2
727   AliBitPacking::PackWord(word,baseWord,28,31);
728   fIndex++;
729   buf[fIndex]=baseWord;
730
731 }
732
733 //----------------------------------------------------------------------------
734 void AliTOFDDLRawData::MakeTRMchainHeader(Int_t nTRM, Int_t iChain,
735                                           UInt_t *buf)
736 {
737   //
738   // TRM chain header
739   //
740   
741   UInt_t baseWord;
742   UInt_t word;
743
744   if (nTRM<3 || nTRM>12) {
745     AliWarning(Form(" TRM number is out of the right range [3;12] (nTRM = %3i", nTRM));
746     return;
747   }
748   
749   if (iChain<0 || iChain>1) {
750     AliWarning(Form(" Chain number is out of the right range [0;1] (iChain = %3i", iChain));
751     return;
752   }
753
754   baseWord=0;
755   word = nTRM; // TRM data coming from the VME slot ID nTRM
756   AliBitPacking::PackWord(word,baseWord, 0, 3);
757   word = 0; // bunch ID
758   AliBitPacking::PackWord(word,baseWord, 4,15);
759   word = 100; // PB24 temperature -> 4 X 25 degree (environment temperature)
760   AliBitPacking::PackWord(word,baseWord,16,23);
761   word = (Int_t)(5 * gRandom->Rndm()); // PB24 ID [0;4]
762   AliBitPacking::PackWord(word,baseWord,24,26);
763   word = 0; // TS
764   AliBitPacking::PackWord(word,baseWord,27,27);
765   switch (iChain) {
766     case 0:
767       word = 0; // 0000 -> TRM chain 0 ID
768       break;
769     case 1:
770       word = 2; // 0010 -> TRM chain 1 ID
771       break;
772     }
773   AliBitPacking::PackWord(word,baseWord,28,31);
774   fIndex++;
775   buf[fIndex]=baseWord;
776             
777 }
778
779 //----------------------------------------------------------------------------
780 void AliTOFDDLRawData::MakeTRMfiller(UInt_t *buf, UInt_t nWordsPerTRM)
781 {
782   //
783   // TRM filler
784   //
785
786   Int_t jj = -1;
787
788   fIndex++;
789   for (jj=fIndex; jj>fIndex-(Int_t)nWordsPerTRM; jj--) {
790     buf[jj] = buf[jj-1];
791   }
792
793   buf[fIndex-nWordsPerTRM] = MakeFiller();
794
795 }
796   
797 //----------------------------------------------------------------------------
798 UInt_t AliTOFDDLRawData::MakeFiller()
799 {
800   //
801   // Filler word definition: to make even the number of words per TRM/LTM
802   //
803
804   UInt_t baseWord;
805   UInt_t word;
806
807   baseWord=0;
808   word = 0; // 0000 -> filler ID 1
809   AliBitPacking::PackWord(word,baseWord, 0, 3);
810   word = 0; // MBZ
811   AliBitPacking::PackWord(word,baseWord, 4,27);
812   word = 7; // 0111 -> filler ID 2
813   AliBitPacking::PackWord(word,baseWord, 28,31);
814   
815   return baseWord;
816
817 }
818
819 //----------------------------------------------------------------------------
820 void AliTOFDDLRawData::MakeTRMchainTrailer(Int_t iChain, UInt_t *buf)
821 {
822   //
823   // TRM chain trailer
824   //
825
826   if (iChain<0 || iChain>1) {
827     AliWarning(Form(" Chain number is out of the right range [0;1] (iChain = %3i", iChain));
828     return;
829   }
830
831   UInt_t baseWord;
832   UInt_t word;
833   
834   baseWord=0;
835   word = 0; // status
836   AliBitPacking::PackWord(word,baseWord, 0, 3);
837   word = 0; // MBZ
838   AliBitPacking::PackWord(word,baseWord, 4,15);
839   word = 0; // event counter
840   AliBitPacking::PackWord(word,baseWord,16,27);
841   switch (iChain) {
842     case 0:
843       word = 1; // 0001 -> TRM chain 0 trailer ID
844       break;
845     case 1:
846       word = 3; // 0011 -> TRM chain 1 trailer ID
847       break;
848     }
849   AliBitPacking::PackWord(word,baseWord,28,31);
850   fIndex++;
851   buf[fIndex]=baseWord;
852
853 }
854
855 //----------------------------------------------------------------------------
856 void AliTOFDDLRawData::MakeTDCdigits(Int_t nDDL, Int_t nTRM, Int_t iChain,
857                                      UInt_t *buf, UInt_t &nWordsPerTRM)
858 {
859   //
860   // TRM TDC digit
861   //
862
863   const Double_t kOneMoreFilledCell = 1./(AliTOFGeometry::NPadXSector()*AliTOFGeometry::NSectors());
864   Double_t percentFilledCells = Double_t(fTOFdigitMap->GetFilledCellNumber())/(AliTOFGeometry::NPadXSector()*AliTOFGeometry::NSectors());
865
866   if (nDDL<0 || nDDL>71) {
867     AliWarning(Form(" DDL number is out of the right range [0;71] (nDDL = %3i", nDDL));
868     return;
869   }
870   
871   if (nTRM<3 || nTRM>12) {
872     AliWarning(Form(" TRM number is out of the right range [3;12] (nTRM = %3i", nTRM));
873     return;
874   }
875   
876   if (iChain<0 || iChain>1) {
877     AliWarning(Form(" Chain number is out of the right range [0;1] (iChain = %3i", iChain));
878     return;
879   }
880   
881   Int_t psArray[1000];
882   UInt_t localBuffer[1000];
883   Int_t localIndex = -1;
884
885   Int_t iDDL = nDDL%AliTOFGeometry::NDDL();
886
887   Int_t volume[5] = {-1, -1, -1, -1, -1};
888   Int_t indexDigit[3] = {-1, -1, -1};
889
890   Int_t totCharge = -1;
891   Int_t timeOfFlight = -1;
892
893   Int_t trailingSpurious = -1;
894   Int_t leadingSpurious = -1;
895
896   AliTOFdigit *digs;
897
898   UInt_t baseWord=0;
899   UInt_t word=0;
900
901   Int_t jj = -1;
902   Int_t nTDC = -1;
903   Int_t iCH = -1;
904
905   //Int_t numberOfMeasuresPerChannel = 0;
906   //Int_t maxMeasuresPerChannelInTDC = 0;
907
908   Bool_t outOut = HeadOrTail();
909
910   ofstream ftxt;
911
912   if (fVerbose==2) ftxt.open("TOFdigits.txt",ios::app);
913
914   for (jj=0; jj<5; jj++) volume[jj] = -1;
915
916   // loop on TDC number
917   for (nTDC=AliTOFGeometry::NTdc()-1; nTDC>=0; nTDC--) {
918
919     // the DRM odd (i.e. left) slot number 3 doesn't contain TDC digit data
920     // for TDC numbers 3-14
921     if (iDDL%2==1 && nTRM==3 && (Int_t)(nTDC/3.)!=0) continue;
922
923     // loop on TDC channel number
924     for (iCH=AliTOFGeometry::NCh()-1; iCH>=0; iCH--) {
925
926       //numberOfMeasuresPerChannel = 0;
927
928       fTOFrawStream->EquipmentId2VolumeId(nDDL, nTRM, iChain, nTDC, iCH, volume);
929         
930       if (volume[0]==-1 || volume[1]==-1 || volume[2]==-1 ||
931           volume[3]==-1 || volume[4]==-1) continue;
932
933       for (jj=0; jj<3; jj++) indexDigit[jj] = -1;
934
935       fTOFdigitMap->GetDigitIndex(volume, indexDigit);
936
937       if (indexDigit[0]<0) {
938
939         trailingSpurious = Int_t(2097152*gRandom->Rndm());
940         leadingSpurious = Int_t(2097152*gRandom->Rndm());
941
942         if ( fFakeOrphaneProduction &&
943              ( ( fPackedAcquisition && percentFilledCells<0.12 && gRandom->Rndm()<(0.12-percentFilledCells) ) ||
944                (!fPackedAcquisition && percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells) )  )  ) {
945
946           percentFilledCells+=kOneMoreFilledCell;
947
948           Int_t dummyPS = 0;
949
950           if (outOut) {
951             word = trailingSpurious; // trailing edge measurement
952             dummyPS = 2;
953           }
954           else {
955             word = leadingSpurious; // leading edge measurement
956             dummyPS = 1;
957           }
958
959           if (fVerbose==2) {
960             if (nDDL<10) ftxt << "  " << nDDL;
961             else         ftxt << " " << nDDL;
962             if (nTRM<10) ftxt << "  " << nTRM;
963             else         ftxt << " " << nTRM;
964             ftxt << "  " << iChain;
965             if (nTDC<10) ftxt << "  " << nTDC;
966             else         ftxt << " " << nTDC;
967             ftxt << "  " << iCH;
968             if (volume[0]<10) ftxt  << "  ->  " << volume[0];
969             else              ftxt  << "  -> " << volume[0];
970             ftxt << "  " << volume[1];
971             if (volume[2]<10) ftxt << "  " << volume[2];
972             else              ftxt << " " << volume[2];
973             ftxt << "  " << volume[4];
974             if (volume[3]<10) ftxt << "  " << volume[3];
975             else              ftxt << " " << volume[3];
976             ftxt << "   " << -1;
977             if (word<10)                           ftxt << "        " << word;
978             else if (word>=10     && word<100)     ftxt << "       " << word;
979             else if (word>=100    && word<1000)    ftxt << "      " << word;
980             else if (word>=1000   && word<10000)   ftxt << "     " << word;
981             else if (word>=10000  && word<100000)  ftxt << "    " << word;
982             else if (word>=100000 && word<1000000) ftxt << "   " << word;
983             else                                   ftxt << "  " << word;
984             ftxt << "   " << dummyPS << endl;
985           }
986
987           AliBitPacking::PackWord(word,baseWord, 0,20);
988           word = iCH; // TDC channel ID [0;7]
989           AliBitPacking::PackWord(word,baseWord,21,23);
990           word = nTDC; // TDC ID [0;14]
991           AliBitPacking::PackWord(word,baseWord,24,27);
992           word = 0; // error flag
993           AliBitPacking::PackWord(word,baseWord,28,28);
994           word = dummyPS; // Packing Status [0;3]
995           AliBitPacking::PackWord(word,baseWord,29,30);
996           word = 1; // TRM TDC digit ID
997           AliBitPacking::PackWord(word,baseWord,31,31);
998
999           localIndex++;
1000           localBuffer[localIndex]=baseWord;
1001           psArray[localIndex]=dummyPS;
1002
1003           nWordsPerTRM++;
1004           baseWord=0;
1005
1006         } // if ( fFakeOrphaneProduction && ( ( fPackedAcquisition && percentFilledCells<0.12 && gRandom->Rndm()<(0.12-percentFilledCells) ) or ... ) )
1007       } // if (indexDigit[0]<0)
1008
1009       for (jj=0; jj<3;jj++) {
1010
1011         if (indexDigit[jj]<0) continue;
1012
1013         digs = (AliTOFdigit*)fTOFdigitArray->UncheckedAt(indexDigit[jj]);
1014           
1015         if (digs->GetSector()!=volume[0] ||
1016             digs->GetPlate() !=volume[1] ||
1017             digs->GetStrip() !=volume[2] ||
1018             digs->GetPadx()  !=volume[3] ||
1019             digs->GetPadz()  !=volume[4]) AliWarning(" --- ERROR --- ");
1020
1021         timeOfFlight = (Int_t)(digs->GetTdc());
1022
1023         if (timeOfFlight>=fMatchingWindow) continue;
1024
1025         //numberOfMeasuresPerChannel++;
1026
1027         // totCharge = (Int_t)digs->GetAdc(); //Use realistic ToT, for Standard production with no miscalibration/Slewing it == fAdC in digit (see AliTOFDigitizer)
1028         totCharge = (Int_t)(digs->GetToT());
1029         // temporary control
1030         if (totCharge<0) totCharge = 0;//TMath::Abs(totCharge);
1031
1032         if (fPackedAcquisition) {
1033
1034         if (fVerbose==2) {
1035           if (nDDL<10) ftxt << "  " << nDDL;
1036           else         ftxt << " " << nDDL;
1037           if (nTRM<10) ftxt << "  " << nTRM;
1038           else         ftxt << " " << nTRM;
1039           ftxt << "  " << iChain;
1040           if (nTDC<10) ftxt << "  " << nTDC;
1041           else         ftxt << " " << nTDC;
1042           ftxt << "  " << iCH;
1043           if (volume[0]<10) ftxt  << "  ->  " << volume[0];
1044           else              ftxt  << "  -> " << volume[0];
1045           ftxt << "  " << volume[1];
1046           if (volume[2]<10) ftxt << "  " << volume[2];
1047           else              ftxt << " " << volume[2];
1048           ftxt << "  " << volume[4];
1049           if (volume[3]<10) ftxt << "  " << volume[3];
1050           else              ftxt << " " << volume[3];
1051           if (totCharge<10)                        ftxt << "    " << totCharge;
1052           else if (totCharge>=10 && totCharge<100) ftxt << "   " << totCharge;
1053           else                                     ftxt << "  " << totCharge;
1054           if (timeOfFlight<10)                             ftxt << "     " << timeOfFlight << endl;
1055           else if (timeOfFlight>=10  && timeOfFlight<100)  ftxt << "    " << timeOfFlight << endl;
1056           else if (timeOfFlight>=100 && timeOfFlight<1000) ftxt << "   " << timeOfFlight << endl;
1057           else                                             ftxt << "  " << timeOfFlight << endl;
1058         }
1059
1060         word = timeOfFlight%8192; // time-of-fligth measurement
1061         AliBitPacking::PackWord(word,baseWord, 0,12);
1062
1063         if (totCharge>=256) totCharge = 255;
1064         word = totCharge; // time-over-threshould measurement
1065         AliBitPacking::PackWord(word,baseWord,13,20);
1066
1067         word = iCH; // TDC channel ID [0;7]
1068         AliBitPacking::PackWord(word,baseWord,21,23);
1069         word = nTDC; // TDC ID [0;14]
1070         AliBitPacking::PackWord(word,baseWord,24,27);
1071         word = 0; // error flag
1072         AliBitPacking::PackWord(word,baseWord,28,28);
1073         word = 0; // Packing Status [0;3]
1074         AliBitPacking::PackWord(word,baseWord,29,30);
1075         word = 1; // TRM TDC digit ID
1076         AliBitPacking::PackWord(word,baseWord,31,31);
1077
1078         localIndex++;
1079         localBuffer[localIndex]=baseWord;
1080
1081         nWordsPerTRM++;
1082         baseWord=0;
1083
1084         if ( fFakeOrphaneProduction &&
1085              percentFilledCells<0.12 && gRandom->Rndm()<(0.12-percentFilledCells) ) {
1086
1087           percentFilledCells+=kOneMoreFilledCell;
1088
1089           trailingSpurious = Int_t(2097152*gRandom->Rndm());
1090           leadingSpurious = Int_t(2097152*gRandom->Rndm());
1091
1092           Int_t dummyPS = 0;
1093
1094           if (outOut) {
1095             word = trailingSpurious; // trailing edge measurement
1096             dummyPS = 2;
1097           }
1098           else {
1099             word = leadingSpurious; // leading edge measurement
1100             dummyPS = 1;
1101           }
1102
1103           if (fVerbose==2) {
1104             if (nDDL<10) ftxt << "  " << nDDL;
1105             else         ftxt << " " << nDDL;
1106             if (nTRM<10) ftxt << "  " << nTRM;
1107             else         ftxt << " " << nTRM;
1108             ftxt << "  " << iChain;
1109             if (nTDC<10) ftxt << "  " << nTDC;
1110             else         ftxt << " " << nTDC;
1111             ftxt << "  " << iCH;
1112             if (volume[0]<10) ftxt  << "  ->  " << volume[0];
1113             else              ftxt  << "  -> " << volume[0];
1114             ftxt << "  " << volume[1];
1115             if (volume[2]<10) ftxt << "  " << volume[2];
1116             else              ftxt << " " << volume[2];
1117             ftxt << "  " << volume[4];
1118             if (volume[3]<10) ftxt << "  " << volume[3];
1119             else              ftxt << " " << volume[3];
1120             ftxt << "   " << -1;
1121             if (word<10)                           ftxt << "        " << word;
1122             else if (word>=10     && word<100)     ftxt << "       " << word;
1123             else if (word>=100    && word<1000)    ftxt << "      " << word;
1124             else if (word>=1000   && word<10000)   ftxt << "     " << word;
1125             else if (word>=10000  && word<100000)  ftxt << "    " << word;
1126             else if (word>=100000 && word<1000000) ftxt << "   " << word;
1127             else                                   ftxt << "  " << word;
1128             ftxt << "   " << dummyPS << endl;
1129           }
1130
1131           AliBitPacking::PackWord(word,baseWord, 0,20);
1132           word = iCH; // TDC channel ID [0;7]
1133           AliBitPacking::PackWord(word,baseWord,21,23);
1134           word = nTDC; // TDC ID [0;14]
1135           AliBitPacking::PackWord(word,baseWord,24,27);
1136           word = 0; // error flag
1137           AliBitPacking::PackWord(word,baseWord,28,28);
1138           word = dummyPS; // Packing Status [0;3]
1139           AliBitPacking::PackWord(word,baseWord,29,30);
1140           word = 1; // TRM TDC digit ID
1141           AliBitPacking::PackWord(word,baseWord,31,31);
1142
1143           localIndex++;
1144           localBuffer[localIndex]=baseWord;
1145           psArray[localIndex]=dummyPS;
1146
1147           nWordsPerTRM++;
1148           baseWord=0;
1149
1150         } // if ( fFakeOrphaneProduction && percentFilledCells<0.12 && gRandom->Rndm()<(0.12-percentFilledCells) )
1151
1152
1153         } // if (fPackedAcquisition)
1154         else { // if (!fPackedAcquisition)
1155
1156         if ( fFakeOrphaneProduction &&
1157              percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells) && outOut ) {
1158
1159           percentFilledCells+=kOneMoreFilledCell;
1160
1161           trailingSpurious = Int_t(2097152*gRandom->Rndm());
1162           word = trailingSpurious;
1163           Int_t dummyPS = 2;
1164
1165           if (fVerbose==2) {
1166             if (nDDL<10) ftxt << "  " << nDDL;
1167             else         ftxt << " " << nDDL;
1168             if (nTRM<10) ftxt << "  " << nTRM;
1169             else         ftxt << " " << nTRM;
1170             ftxt << "  " << iChain;
1171             if (nTDC<10) ftxt << "  " << nTDC;
1172             else         ftxt << " " << nTDC;
1173             ftxt << "  " << iCH;
1174             if (volume[0]<10) ftxt  << "  ->  " << volume[0];
1175             else              ftxt  << "  -> " << volume[0];
1176             ftxt << "  " << volume[1];
1177             if (volume[2]<10) ftxt << "  " << volume[2];
1178             else              ftxt << " " << volume[2];
1179             ftxt << "  " << volume[4];
1180             if (volume[3]<10) ftxt << "  " << volume[3];
1181             else              ftxt << " " << volume[3];
1182             ftxt << "   " << -1;
1183             if (word<10)                           ftxt << "        " << word;
1184             else if (word>=10     && word<100)     ftxt << "       " << word;
1185             else if (word>=100    && word<1000)    ftxt << "      " << word;
1186             else if (word>=1000   && word<10000)   ftxt << "     " << word;
1187             else if (word>=10000  && word<100000)  ftxt << "    " << word;
1188             else if (word>=100000 && word<1000000) ftxt << "   " << word;
1189             else                                   ftxt << "  " << word;
1190             ftxt << "   " << dummyPS << endl;
1191           }
1192
1193           AliBitPacking::PackWord(word,baseWord, 0,20);
1194           word = iCH; // TDC channel ID [0;7]
1195           AliBitPacking::PackWord(word,baseWord,21,23);
1196           word = nTDC; // TDC ID [0;14]
1197           AliBitPacking::PackWord(word,baseWord,24,27);
1198           word = 0; // error flag
1199           AliBitPacking::PackWord(word,baseWord,28,28);
1200           word = dummyPS; // Packing Status [0;3]
1201           AliBitPacking::PackWord(word,baseWord,29,30);
1202           word = 1; // TRM TDC digit ID
1203           AliBitPacking::PackWord(word,baseWord,31,31);
1204
1205           localIndex++;
1206           localBuffer[localIndex]=baseWord;
1207           psArray[localIndex]=dummyPS;
1208
1209           nWordsPerTRM++;
1210           baseWord=0;
1211
1212         } // if ( fFakeOrphaneProduction && percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells)  && outOut )
1213
1214
1215         word = (timeOfFlight + Int_t(totCharge*AliTOFGeometry::ToTBinWidth()/AliTOFGeometry::TdcBinWidth()))%2097152; // trailing edge measurement
1216
1217         if (fVerbose==2) {
1218           if (nDDL<10) ftxt << "  " << nDDL;
1219           else         ftxt << " " << nDDL;
1220           if (nTRM<10) ftxt << "  " << nTRM;
1221           else         ftxt << " " << nTRM;
1222           ftxt << "  " << iChain;
1223           if (nTDC<10) ftxt << "  " << nTDC;
1224           else         ftxt << " " << nTDC;
1225           ftxt << "  " << iCH;
1226           if (volume[0]<10) ftxt  << "  ->  " << volume[0];
1227           else              ftxt  << "  -> " << volume[0];
1228           ftxt << "  " << volume[1];
1229           if (volume[2]<10) ftxt << "  " << volume[2];
1230           else              ftxt << " " << volume[2];
1231           ftxt << "  " << volume[4];
1232           if (volume[3]<10) ftxt << "  " << volume[3];
1233           else              ftxt << " " << volume[3];
1234           ftxt << "   " << -1;
1235           if (word<10)                           ftxt << "        " << word;
1236           else if (word>=10     && word<100)     ftxt << "       " << word;
1237           else if (word>=100    && word<1000)    ftxt << "      " << word;
1238           else if (word>=1000   && word<10000)   ftxt << "     " << word;
1239           else if (word>=10000  && word<100000)  ftxt << "    " << word;
1240           else if (word>=100000 && word<1000000) ftxt << "   " << word;
1241           else                                   ftxt << "  " << word;
1242           ftxt << "   " << 2 << endl;
1243         }
1244
1245         AliBitPacking::PackWord(word,baseWord, 0,20);
1246
1247         word = iCH; // TDC channel ID [0;7]
1248         AliBitPacking::PackWord(word,baseWord,21,23);
1249         word = nTDC; // TDC ID [0;14]
1250         AliBitPacking::PackWord(word,baseWord,24,27);
1251         word = 0; // error flag
1252         AliBitPacking::PackWord(word,baseWord,28,28);
1253         word = 2; // Packing Status [0;3]
1254         AliBitPacking::PackWord(word,baseWord,29,30);
1255         word = 1; // TRM TDC digit ID
1256         AliBitPacking::PackWord(word,baseWord,31,31);
1257
1258         localIndex++;
1259         localBuffer[localIndex]=baseWord;
1260         psArray[localIndex]=2;
1261
1262         nWordsPerTRM++;
1263         baseWord=0;
1264
1265         word = timeOfFlight%2097152; // leading edge measurement
1266
1267         if (fVerbose==2) {
1268           if (nDDL<10) ftxt << "  " << nDDL;
1269           else         ftxt << " " << nDDL;
1270           if (nTRM<10) ftxt << "  " << nTRM;
1271           else         ftxt << " " << nTRM;
1272           ftxt << "  " << iChain;
1273           if (nTDC<10) ftxt << "  " << nTDC;
1274           else         ftxt << " " << nTDC;
1275           ftxt << "  " << iCH;
1276           if (volume[0]<10) ftxt  << "  ->  " << volume[0];
1277           else              ftxt  << "  -> " << volume[0];
1278           ftxt << "  " << volume[1];
1279           if (volume[2]<10) ftxt << "  " << volume[2];
1280           else              ftxt << " " << volume[2];
1281           ftxt << "  " << volume[4];
1282           if (volume[3]<10) ftxt << "  " << volume[3];
1283           else              ftxt << " " << volume[3];
1284           ftxt << "   " << -1;
1285           if (word<10)                           ftxt << "        " << word;
1286           else if (word>=10     && word<100)     ftxt << "       " << word;
1287           else if (word>=100    && word<1000)    ftxt << "      " << word;
1288           else if (word>=1000   && word<10000)   ftxt << "     " << word;
1289           else if (word>=10000  && word<100000)  ftxt << "    " << word;
1290           else if (word>=100000 && word<1000000) ftxt << "   " << word;
1291           else                                   ftxt << "  " << word;
1292           ftxt << "   " << 1 << endl;
1293         }
1294
1295         AliBitPacking::PackWord(word,baseWord, 0,20);
1296
1297         word = iCH; // TDC channel ID [0;7]
1298         AliBitPacking::PackWord(word,baseWord,21,23);
1299         word = nTDC; // TDC ID [0;14]
1300         AliBitPacking::PackWord(word,baseWord,24,27);
1301         word = 0; // error flag
1302         AliBitPacking::PackWord(word,baseWord,28,28);
1303         word = 1; // Packing Status [0;3]
1304         AliBitPacking::PackWord(word,baseWord,29,30);
1305         word = 1; // TRM TDC digit ID
1306         AliBitPacking::PackWord(word,baseWord,31,31);
1307
1308         localIndex++;
1309         localBuffer[localIndex]=baseWord;
1310         psArray[localIndex]=1;
1311
1312         nWordsPerTRM++;
1313         baseWord=0;
1314
1315
1316         if ( fFakeOrphaneProduction &&
1317              percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells) && !outOut ) {
1318
1319           percentFilledCells+=kOneMoreFilledCell;
1320
1321           leadingSpurious = Int_t(2097152*gRandom->Rndm());
1322           word = leadingSpurious;
1323           Int_t dummyPS = 1;
1324
1325           if (fVerbose==2) {
1326             if (nDDL<10) ftxt << "  " << nDDL;
1327             else         ftxt << " " << nDDL;
1328             if (nTRM<10) ftxt << "  " << nTRM;
1329             else         ftxt << " " << nTRM;
1330             ftxt << "  " << iChain;
1331             if (nTDC<10) ftxt << "  " << nTDC;
1332             else         ftxt << " " << nTDC;
1333             ftxt << "  " << iCH;
1334             if (volume[0]<10) ftxt  << "  ->  " << volume[0];
1335             else              ftxt  << "  -> " << volume[0];
1336             ftxt << "  " << volume[1];
1337             if (volume[2]<10) ftxt << "  " << volume[2];
1338             else              ftxt << " " << volume[2];
1339             ftxt << "  " << volume[4];
1340             if (volume[3]<10) ftxt << "  " << volume[3];
1341             else              ftxt << " " << volume[3];
1342             ftxt << "   " << -1;
1343             if (word<10)                           ftxt << "        " << word;
1344             else if (word>=10     && word<100)     ftxt << "       " << word;
1345             else if (word>=100    && word<1000)    ftxt << "      " << word;
1346             else if (word>=1000   && word<10000)   ftxt << "     " << word;
1347             else if (word>=10000  && word<100000)  ftxt << "    " << word;
1348             else if (word>=100000 && word<1000000) ftxt << "   " << word;
1349             else                                   ftxt << "  " << word;
1350             ftxt << "   " << dummyPS << endl;
1351           }
1352
1353           AliBitPacking::PackWord(word,baseWord, 0,20);
1354           word = iCH; // TDC channel ID [0;7]
1355           AliBitPacking::PackWord(word,baseWord,21,23);
1356           word = nTDC; // TDC ID [0;14]
1357           AliBitPacking::PackWord(word,baseWord,24,27);
1358           word = 0; // error flag
1359           AliBitPacking::PackWord(word,baseWord,28,28);
1360           word = dummyPS; // Packing Status [0;3]
1361           AliBitPacking::PackWord(word,baseWord,29,30);
1362           word = 1; // TRM TDC digit ID
1363           AliBitPacking::PackWord(word,baseWord,31,31);
1364
1365           localIndex++;
1366           localBuffer[localIndex]=baseWord;
1367           psArray[localIndex]=dummyPS;
1368
1369           nWordsPerTRM++;
1370           baseWord=0;
1371
1372         } // if ( fFakeOrphaneProduction && percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells) && !outOut )
1373
1374
1375         } // if (!fPackedAcquisition)
1376
1377       } //end loop on digits in the same volume
1378
1379       //if (numberOfMeasuresPerChannel>maxMeasuresPerChannelInTDC)
1380       //maxMeasuresPerChannelInTDC = numberOfMeasuresPerChannel;
1381
1382     } // end loop on TDC channel number
1383
1384     //AliInfo(Form(" TDC number %2i:  numberOfMeasuresPerChannel = %2i  ---  maxMeasuresPerChannelInTDC = %2i ", nTDC, numberOfMeasuresPerChannel, maxMeasuresPerChannelInTDC));
1385
1386     if (localIndex==-1) continue;
1387
1388     if (fPackedAcquisition) {
1389
1390       for (Int_t jj=0; jj<=localIndex; jj++) {
1391         fIndex++;
1392         buf[fIndex] = localBuffer[jj];
1393         localBuffer[jj] = 0;
1394         psArray[jj] = -1;
1395       }
1396
1397     }
1398     else {
1399       /*
1400       if (maxMeasuresPerChannelInTDC = 1) {
1401
1402         for (Int_t jj=0; jj<=localIndex; jj++) {
1403           if (psArray[jj]==2) {
1404             fIndex++;
1405             buf[fIndex] = localBuffer[jj];
1406             localBuffer[jj] = 0;
1407             psArray[jj] = -1;
1408           }
1409         }
1410         for (Int_t jj=0; jj<=localIndex; jj++) {
1411           if (psArray[jj]==1) {
1412             fIndex++;
1413             buf[fIndex] = localBuffer[jj];
1414             localBuffer[jj] = 0;
1415             psArray[jj] = -1;
1416           }
1417         }
1418
1419       } // if (maxMeasuresPerChannelInTDC = 1)
1420       else if (maxMeasuresPerChannelInTDC>1) {
1421
1422         AliInfo(Form(" In the TOF DDL %2i, TRM %2i, TDC %2i, chain %1i, the maximum number of t.o.f. good measurements per channel is %2i",
1423                      nDDL, nTRM, iChain, nTDC, iCH, maxMeasuresPerChannelInTDC));
1424       */
1425         for (Int_t jj=0; jj<=localIndex; jj++) {
1426             fIndex++;
1427             buf[fIndex] = localBuffer[jj];
1428             localBuffer[jj] = 0;
1429             psArray[jj] = -1;
1430         }
1431
1432         //} // else if (maxMeasuresPerChannelInTDC>1)
1433
1434     } // else (!fPackedAcquisition)
1435
1436     localIndex = -1;
1437
1438     //maxMeasuresPerChannelInTDC = 0;
1439
1440   } // end loop on TDC number
1441
1442
1443   if (fVerbose==2) ftxt.close();
1444
1445 }
1446
1447 //----------------------------------------------------------------------------
1448 void AliTOFDDLRawData::ReverseArray(UInt_t a[], Int_t n) const
1449 {
1450   //
1451   // Reverses the n elements of array a
1452   //
1453
1454   Int_t ii, temp;
1455
1456   for (ii=0; ii<n/2; ii++) {
1457     temp      = a[ii];
1458     a[ii]     = a[n-ii-1];
1459     a[n-ii-1] = temp;
1460   }
1461
1462   return;
1463
1464 }
1465
1466 //----------------------------------------------------------------------------
1467 Bool_t AliTOFDDLRawData::HeadOrTail() const
1468 {
1469   //
1470   // Returns the result of a 'pitch and toss'
1471   //
1472
1473   Double_t dummy = gRandom->Rndm();
1474
1475   if (dummy<0.5) return kFALSE;
1476   else return kTRUE;
1477
1478 }