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