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