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