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