1 /**************************************************************************
2 * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.13 2007/02/20 15:57:00 decaro
19 Raw data update: to read the TOF raw data defined in UNPACKED mode
21 Revision 1.12 2006/08/22 13:29:42 arcelli
22 removal of effective c++ warnings (C.Zampolli)
24 Revision 1.11 2006/08/10 14:46:54 decaro
25 TOF raw data format: updated version
27 Revision 1.10.1 2006/06/28 A.De Caro
28 Update TOF raw data format
29 according to the final version
30 (see ALICE internal note in preparation
31 'ALICE TOF raw data format')
33 Revision 0.02 2005/7/25 A.De Caro
34 Update number of bits allocated for time-of-flight
35 and 'charge' measurements
37 Revision 0.01 2004/6/11 A.De Caro, S.B.Sellitto, R.Silvestri
38 First implementation: global methods RawDataTOF
42 ////////////////////////////////////////////////////////////////////
44 // This class contains the methods to create the Raw Data files //
45 // for the TOF detector starting from the Digits. //
46 // In this implementation, we defined the structure //
47 // of the ALICE-TOF raw data (according to the //
48 // ALICE technical note, in preparation) //
49 // starting from the TOF digit format. //
51 ////////////////////////////////////////////////////////////////////
53 #include "Riostream.h"
56 #include "TClonesArray.h"
60 #include "AliBitPacking.h"
63 //#include "AliRawDataHeader.h"
64 #include "AliRawDataHeaderSim.h"
66 #include "AliTOFDDLRawData.h"
67 #include "AliTOFDigitMap.h"
68 #include "AliTOFdigit.h"
69 #include "AliTOFGeometry.h"
70 #include "AliTOFRawStream.h"
72 extern TRandom *gRandom;
74 ClassImp(AliTOFDDLRawData)
76 //---------------------------------------------------------------------------
77 AliTOFDDLRawData::AliTOFDDLRawData():
80 fPackedAcquisition(kTRUE),
82 fTOFdigitMap(new AliTOFDigitMap()),
84 fTOFrawStream(new AliTOFRawStream())
89 //----------------------------------------------------------------------------
90 AliTOFDDLRawData::AliTOFDDLRawData(AliTOFGeometry *tofGeom):
93 fPackedAcquisition(kTRUE),
94 fTOFgeometry(tofGeom),
95 fTOFdigitMap(new AliTOFDigitMap()),
97 fTOFrawStream(new AliTOFRawStream())
103 //----------------------------------------------------------------------------
104 AliTOFDDLRawData::AliTOFDDLRawData(const AliTOFDDLRawData &source) :
108 fPackedAcquisition(kTRUE),
110 fTOFdigitMap(new AliTOFDigitMap()),
112 fTOFrawStream(new AliTOFRawStream())
115 this->fIndex=source.fIndex;
116 this->fVerbose=source.fVerbose;
117 this->fPackedAcquisition=source.fPackedAcquisition;
118 this->fTOFgeometry=source.fTOFgeometry;
119 this->fTOFdigitMap=source.fTOFdigitMap;
120 this->fTOFdigitArray=source.fTOFdigitArray;
121 this->fTOFrawStream=source.fTOFrawStream;
125 //----------------------------------------------------------------------------
126 AliTOFDDLRawData& AliTOFDDLRawData::operator=(const AliTOFDDLRawData &source) {
128 this->fIndex=source.fIndex;
129 this->fVerbose=source.fVerbose;
130 this->fPackedAcquisition=source.fPackedAcquisition;
131 this->fTOFgeometry=source.fTOFgeometry;
132 this->fTOFdigitMap=source.fTOFdigitMap;
133 this->fTOFdigitArray=source.fTOFdigitArray;
134 this->fTOFrawStream=source.fTOFrawStream;
138 //----------------------------------------------------------------------------
139 Int_t AliTOFDDLRawData::RawDataTOF(TBranch* branch)
142 // This method creates the Raw data files for TOF detector
145 const Int_t kSize = 5000; //max number of digits per DDL file times 2
151 fTOFdigitArray = * (TClonesArray**) branch->GetAddress();
154 ofstream outfile; // logical name of the output file
156 //AliRawDataHeader header;
157 AliRawDataHeaderSim header;
159 UInt_t sizeRawData = 0;
171 UInt_t nWordsPerTRM = 0;
173 //loop over TOF DDL files
174 for (nDDL=0; nDDL<AliDAQ::NumberOfDdls("TOF"); nDDL++) {
176 strcpy(fileName,AliDAQ::DdlFileName("TOF",nDDL)); //The name of the output file
178 outfile.open(fileName,ios::binary);
180 outfile.open(fileName);
183 iDDL = fTOFrawStream->GetDDLnumberPerSector(nDDL);
185 // write Dummy DATA HEADER
186 UInt_t dataHeaderPosition = outfile.tellp();
187 outfile.write((char*)(&header),sizeof(header));
189 // DRM section: trailer
194 buf[fIndex] = MakeFiller();
199 // loop on TRM number per DRM
200 for (nTRM=AliTOFGeometry::NTRM(); nTRM>=3; nTRM--) {
204 // the slot number 3 of the even DRM (i.e. right) doesn't contain TDC digit data
205 if (iDDL%2==0 && nTRM==3) continue;
207 // TRM global trailer
210 // loop on TRM chain number per TRM
211 for (iChain=AliTOFGeometry::NChain()-1; iChain>=0; iChain--) {
214 MakeTRMchainTrailer(iChain, buf);
218 MakeTDCdigits(nDDL, nTRM, iChain, buf, nWordsPerTRM);
221 MakeTRMchainHeader(nTRM, iChain, buf);
224 } // end loop on iChain
227 MakeTRMheader(nTRM, buf);
229 // TRM filler in case where TRM data number is odd
230 if (nWordsPerTRM%2!=0) MakeTRMfiller(buf, nWordsPerTRM);
232 } // end loop on nTRM
235 MakeDRMheader(nDDL, buf);
237 ReverseArray(buf, fIndex+1);
239 outfile.write((char *)buf,((fIndex+1)*sizeof(UInt_t)));
241 for (jj=0; jj<(fIndex+1); jj++) buf[jj]=0;
244 //Write REAL DATA HEADER
245 UInt_t currentFilePosition = outfile.tellp();
246 sizeRawData = currentFilePosition - dataHeaderPosition - sizeof(header);
247 header.fSize = currentFilePosition - dataHeaderPosition;
248 header.SetAttribute(0); // valid data
249 outfile.seekp(dataHeaderPosition);
250 outfile.write((char*)(&header),sizeof(header));
251 outfile.seekp(currentFilePosition);
255 } //end loop on DDL file number
261 //----------------------------------------------------------------------------
262 void AliTOFDDLRawData::GetDigits()
265 // Fill the TOF volumes' map with the TOF digit indices
268 Int_t vol[5] = {-1,-1,-1,-1,-1};
271 Int_t ndigits = fTOFdigitArray->GetEntries();
275 // loop on TOF digits
276 for (digit=0; digit<ndigits; digit++) {
277 digs = (AliTOFdigit*)fTOFdigitArray->UncheckedAt(digit);
279 vol[0] = digs->GetSector(); // Sector Number (0-17)
280 vol[1] = digs->GetPlate(); // Plate Number (0-4)
281 vol[2] = digs->GetStrip(); // Strip Number (0-14/18)
282 vol[3] = digs->GetPadx(); // Pad Number in x direction (0-47)
283 vol[4] = digs->GetPadz(); // Pad Number in z direction (0-1)
285 fTOFdigitMap->AddDigit(vol, digit);
287 } // close loop on digit del TOF
291 //----------------------------------------------------------------------------
292 void AliTOFDDLRawData::MakeDRMheader(Int_t nDDL, UInt_t *buf)
298 Int_t iDDL = fTOFrawStream->GetDDLnumberPerSector(nDDL);
300 Int_t iSector = fTOFrawStream->GetSectorNumber(nDDL);
307 word = 1; // 0001 -> DRM data are coming from the VME slot number 1
308 AliBitPacking::PackWord(word,baseWord, 0, 3);
309 word = 0; // event CRC
310 AliBitPacking::PackWord(word,baseWord, 4,19);
311 word = 0; // reserved for future use
312 AliBitPacking::PackWord(word,baseWord,20,27);
313 word = 4; // 0100 -> DRM header ID
314 AliBitPacking::PackWord(word,baseWord,28,31);
316 buf[fIndex]=baseWord;
318 // DRM status header 3
320 word = 1; // 0001 -> DRM data are coming from the VME slot number 1
321 AliBitPacking::PackWord(word,baseWord, 0, 3);
322 word = 0; // TTC event counter
323 AliBitPacking::PackWord(word,baseWord, 4,27);
324 word = 4; // 0100 -> DRM header ID
325 AliBitPacking::PackWord(word,baseWord,28,31);
327 buf[fIndex]=baseWord;
329 // DRM status header 2
331 word = 1; // 0001 -> DRM data are coming from the VME slot number 1
332 AliBitPacking::PackWord(word,baseWord, 0, 3);
335 word = 2047; // enable ID: [00000000000;11111111111] for odd
336 // (i.e. right) crates
337 AliBitPacking::PackWord(word,baseWord, 4,14);
339 word = 2045; // enable ID: [00000000000;11111111101] for even
340 // (i.e. left) crates
341 AliBitPacking::PackWord(word,baseWord, 4,14);
345 AliBitPacking::PackWord(word,baseWord,15,15);
346 word = 0; // fault ID
347 AliBitPacking::PackWord(word,baseWord,16,27);
348 word = 4; // 0100 -> DRM header ID
349 AliBitPacking::PackWord(word,baseWord,28,31);
351 buf[fIndex]=baseWord;
353 // DRM status header 1
355 word = 1; // 0001 -> DRM data are coming from the VME slot number 1
356 AliBitPacking::PackWord(word,baseWord, 0, 3);
359 word = 2047; // slot ID: [00000000000;11111111111] for odd
360 // (i.e. right) crates
361 AliBitPacking::PackWord(word,baseWord, 4,14);
363 word = 2045; // slot ID: [00000000000;11111111101] for even
364 // (i.e. left) crates
365 AliBitPacking::PackWord(word,baseWord, 4,14);
368 word = 1; // LHC clock status: 1/0
369 AliBitPacking::PackWord(word,baseWord,15,15);
370 word = 0; // reserved for future use
371 AliBitPacking::PackWord(word,baseWord,16,27);
372 word = 4; // 0100 -> DRM header ID
373 AliBitPacking::PackWord(word,baseWord,28,31);
375 buf[fIndex]=baseWord;
379 word = 1; // 0001 -> DRM data are coming from the VME slot number 1
380 AliBitPacking::PackWord(word,baseWord, 0, 3);
381 word = fIndex+1 + 1; // event words
382 AliBitPacking::PackWord(word,baseWord, 4,20);
383 word = iDDL; // crate ID [0;3]
384 AliBitPacking::PackWord(word,baseWord,21,22);
385 word = iSector; // sector ID [0;17]
386 AliBitPacking::PackWord(word,baseWord,23,27);
387 word = 4; // 0100 -> DRM header ID
388 AliBitPacking::PackWord(word,baseWord,28,31);
390 buf[fIndex]=baseWord;
394 //----------------------------------------------------------------------------
395 void AliTOFDDLRawData::MakeDRMtrailer(UInt_t *buf)
398 // DRM global trailer
405 word = 1; // 0001 -> DRM data are coming from the VME slot number 1
406 AliBitPacking::PackWord(word,baseWord, 0, 3);
407 word = 0; // local event counter
408 AliBitPacking::PackWord(word,baseWord, 4,15);
409 word = 0; // reserved for future use
410 AliBitPacking::PackWord(word,baseWord,16,27);
411 word = 5; // 0101 -> DRM trailer ID
412 AliBitPacking::PackWord(word,baseWord,28,31);
414 buf[fIndex]=baseWord;
418 //----------------------------------------------------------------------------
419 void AliTOFDDLRawData::MakeLTMheader(UInt_t *buf)
429 word = 2; // 0010 -> LTM data are coming from the VME slot number 2
430 AliBitPacking::PackWord(word,baseWord, 0, 3);
431 word = 35; // event words
432 AliBitPacking::PackWord(word,baseWord, 4,16);
433 word = 0; // crc error
434 AliBitPacking::PackWord(word,baseWord,17,17);
436 AliBitPacking::PackWord(word,baseWord,18,23);
438 AliBitPacking::PackWord(word,baseWord,24,27);
439 word = 4; // 0100 -> LTM header ID
440 AliBitPacking::PackWord(word,baseWord,28,31);
442 buf[fIndex]=baseWord;
446 //----------------------------------------------------------------------------
447 void AliTOFDDLRawData::MakeLTMdata(UInt_t *buf)
457 word = 100; // Local temperature in LTM5 -> 4 X 25 degree (environment temperature)
458 AliBitPacking::PackWord(word,baseWord, 0, 9);
459 word = 100; // Local temperature in LTM6 -> 4 X 25 degree (environment temperature)
460 AliBitPacking::PackWord(word,baseWord,10,19);
461 word = 100; // Local temperature in LTM7 -> 4 X 25 degree (environment temperature)
462 AliBitPacking::PackWord(word,baseWord,20,29);
464 AliBitPacking::PackWord(word,baseWord,30,31);
467 buf[fIndex]=baseWord;
469 // Local temperature in LTM2, LMT3, LTM4 -> 4 X 25 degree (environment temperature)
471 buf[fIndex]=baseWord;
473 // Local temperature in FEAC7, FEAC8, LTM1 -> 4 X 25 degree (environment temperature)
475 buf[fIndex]=baseWord;
477 // Local temperature in FEAC4, FEAC5, FEAC6 -> 4 X 25 degree (environment temperature)
479 buf[fIndex]=baseWord;
481 // Local temperature in FEAC1, FEAC2, FEAC3 -> 4 X 25 degree (environment temperature)
483 buf[fIndex]=baseWord;
486 word = 0; // GND-FEAC15 -> Voltage drop between GND and FEAC15
487 AliBitPacking::PackWord(word,baseWord, 0, 9);
488 word = 0; // VTH16 -> Thereshould voltage for FEAC16
489 AliBitPacking::PackWord(word,baseWord,10,19);
490 word = 0; // GND-FEAC16 -> Voltage drop between GND and FEAC16
491 AliBitPacking::PackWord(word,baseWord,20,29);
493 AliBitPacking::PackWord(word,baseWord,30,31);
496 buf[fIndex]=baseWord;
498 // VTH14 -> Thereshould voltage for FEAC14
499 // GND-FEAC14 -> Voltage drop between GND and FEAC14
500 // VTH15 -> Thereshould voltage for FEAC15
502 buf[fIndex]=baseWord;
504 // GND-FEAC12 -> Voltage drop between GND and FEAC12
505 // VTH13 -> Thereshould voltage for FEAC13
506 // GND-FEAC13 -> Voltage drop between GND and FEAC13
508 buf[fIndex]=baseWord;
510 // VTH11 -> Thereshould voltage for FEAC11
511 // GND-FEAC11 -> Voltage drop between GND and FEAC11
512 // VTH12 -> Thereshould voltage for FEAC12
514 buf[fIndex]=baseWord;
516 // GND-FEAC9 -> Voltage drop between GND and FEAC9
517 // VTH10 -> Thereshould voltage for FEAC10
518 // GND-FEAC10 -> Voltage drop between GND and FEAC10
520 buf[fIndex]=baseWord;
522 // VTH8 -> Thereshould voltage for FEAC8
523 // GND-FEAC8 -> Voltage drop between GND and FEAC8
524 // VTH9 -> Thereshould voltage for FEAC9
526 buf[fIndex]=baseWord;
528 // GND-FEAC6 -> Voltage drop between GND and FEAC6
529 // VTH7 -> Thereshould voltage for FEAC7
530 // GND-FEAC7 -> Voltage drop between GND and FEAC7
532 buf[fIndex]=baseWord;
534 // VTH5 -> Thereshould voltage for FEAC5
535 // GND-FEAC5 -> Voltage drop between GND and FEAC5
536 // VTH6 -> Thereshould voltage for FEAC6
538 buf[fIndex]=baseWord;
540 // GND-FEAC3 -> Voltage drop between GND and FEAC3
541 // VTH4 -> Thereshould voltage for FEAC4
542 // GND-FEAC4 -> Voltage drop between GND and FEAC4
544 buf[fIndex]=baseWord;
546 // VTH2 -> Thereshould voltage for FEAC2
547 // GND-FEAC2 -> Voltage drop between GND and FEAC2
548 // VTH3 -> Thereshould voltage for FEAC3
550 buf[fIndex]=baseWord;
552 // LV16 -> Low Voltage measured by FEAC16
553 // GND-FEAC1 -> Voltage drop between GND and FEAC1
554 // VTH1 -> Thereshould voltage for FEAC1
556 buf[fIndex]=baseWord;
558 // Low Voltage measured by FEAC13, FEAC14, FEAC15
560 buf[fIndex]=baseWord;
562 // Low Voltage measured by FEAC10, FEAC11, FEAC12
564 buf[fIndex]=baseWord;
566 // Low Voltage measured by FEAC7, FEAC8, FEAC9
568 buf[fIndex]=baseWord;
570 // Low Voltage measured by FEAC4, FEAC5, FEAC6
572 buf[fIndex]=baseWord;
574 // Low Voltage measured by FEAC1, FEAC2, FEAC3
576 buf[fIndex]=baseWord;
580 word = 0; // PDL45 -> Delay Line setting for PDL45
581 AliBitPacking::PackWord(word,baseWord, 0, 7);
582 word = 0; // PDL46 -> Delay Line setting for PDL46
583 AliBitPacking::PackWord(word,baseWord, 8,15);
584 word = 0; // PDL47 -> Delay Line setting for PDL47
585 AliBitPacking::PackWord(word,baseWord,16,23);
586 word = 0; // PDL48 -> Delay Line setting for PDL48
587 AliBitPacking::PackWord(word,baseWord,24,31);
589 buf[fIndex]=baseWord;
591 // Delay Line setting for PDL41, PDL42, PDL43, PDL44
593 buf[fIndex]=baseWord;
595 // Delay Line setting for PDL37, PDL38, PDL39, PDL40
597 buf[fIndex]=baseWord;
599 // Delay Line setting for PDL33, PDL34, PDL35, PDL36
601 buf[fIndex]=baseWord;
603 // Delay Line setting for PDL29, PDL30, PDL31, PDL32
605 buf[fIndex]=baseWord;
607 // Delay Line setting for PDL25, PDL26, PDL27, PDL28
609 buf[fIndex]=baseWord;
611 // Delay Line setting for PDL21, PDL22, PDL23, PDL24
613 buf[fIndex]=baseWord;
615 // Delay Line setting for PDL17, PDL18, PDL19, PDL20
617 buf[fIndex]=baseWord;
619 // Delay Line setting for PDL13, PDL14, PDL15, PDL16
621 buf[fIndex]=baseWord;
623 // Delay Line setting for PDL9, PDL10, PDL11, PDL12
625 buf[fIndex]=baseWord;
627 // Delay Line setting for PDL5, PDL6, PDL7, PDL8
629 buf[fIndex]=baseWord;
631 // Delay Line setting for PDL1, PDL2, PDL3, PDL4
633 buf[fIndex]=baseWord;
637 //----------------------------------------------------------------------------
638 void AliTOFDDLRawData::MakeLTMtrailer(UInt_t *buf)
648 word = 2; // 0010 -> LTM data are coming from the VME slot number 2
649 AliBitPacking::PackWord(word,baseWord, 0, 3);
650 word = 0; // event crc
651 AliBitPacking::PackWord(word,baseWord, 4,15);
652 word = 0; // event number
653 AliBitPacking::PackWord(word,baseWord,16,27);
654 word = 5; // 0101 -> LTM trailer ID
655 AliBitPacking::PackWord(word,baseWord,28,31);
657 buf[fIndex]=baseWord;
661 //----------------------------------------------------------------------------
662 void AliTOFDDLRawData::MakeTRMheader(Int_t nTRM, UInt_t *buf)
665 // TRM header for the TRM number nTRM [ 3;12]
668 if (nTRM<3 || nTRM>12) {
669 AliWarning(Form(" TRM number is out of the right range [3;12] (nTRM = %3i",nTRM));
677 word = nTRM; // TRM data coming from the VME slot number nTRM
678 AliBitPacking::PackWord(word,baseWord, 0, 3);
679 word = 0; // event words
680 AliBitPacking::PackWord(word,baseWord, 4,16);
682 if (fPackedAcquisition)
683 word = 0; // ACQuisition mode: [0;3] see document
685 word = 3; // ACQuisition mode: [0;3] see document
686 AliBitPacking::PackWord(word,baseWord,17,18);
687 word = 0; // description of a SEU inside LUT tables for INL compensation;
688 // the data are unaffected
689 AliBitPacking::PackWord(word,baseWord,19,19);
690 word = 0; // Must Be Zero (MBZ)
691 AliBitPacking::PackWord(word,baseWord,20,27);
692 word = 4; // 0100 -> TRM header ID
693 AliBitPacking::PackWord(word,baseWord,28,31);
695 buf[fIndex]=baseWord;
699 //----------------------------------------------------------------------------
700 void AliTOFDDLRawData::MakeTRMtrailer(UInt_t *buf)
710 word = 15; // 1111 -> TRM trailer ID 1
711 AliBitPacking::PackWord(word,baseWord, 0, 3);
712 word = 0; // event CRC
713 AliBitPacking::PackWord(word,baseWord, 4,15);
714 word = 0; // local event counter == DRM local event counter
715 AliBitPacking::PackWord(word,baseWord,16,27);
716 word = 5; // 0101 -> TRM trailer ID 2
717 AliBitPacking::PackWord(word,baseWord,28,31);
719 buf[fIndex]=baseWord;
723 //----------------------------------------------------------------------------
724 void AliTOFDDLRawData::MakeTRMchainHeader(Int_t nTRM, Int_t iChain,
734 if (nTRM<3 || nTRM>12) {
735 AliWarning(Form(" TRM number is out of the right range [3;12] (nTRM = %3i", nTRM));
739 if (iChain<0 || iChain>1) {
740 AliWarning(Form(" Chain number is out of the right range [0;1] (iChain = %3i", iChain));
745 word = nTRM; // TRM data coming from the VME slot ID nTRM
746 AliBitPacking::PackWord(word,baseWord, 0, 3);
747 word = 0; // bunch ID
748 AliBitPacking::PackWord(word,baseWord, 4,15);
749 word = 100; // PB24 temperature -> 4 X 25 degree (environment temperature)
750 AliBitPacking::PackWord(word,baseWord,16,23);
751 word = (Int_t)(5 * gRandom->Rndm()); // PB24 ID [0;4]
752 AliBitPacking::PackWord(word,baseWord,24,26);
754 AliBitPacking::PackWord(word,baseWord,27,27);
757 word = 0; // 0000 -> TRM chain 0 ID
760 word = 2; // 0010 -> TRM chain 1 ID
763 AliBitPacking::PackWord(word,baseWord,28,31);
765 buf[fIndex]=baseWord;
769 //----------------------------------------------------------------------------
770 void AliTOFDDLRawData::MakeTRMfiller(UInt_t *buf, UInt_t nWordsPerTRM)
779 for (jj=fIndex; jj>fIndex-(Int_t)nWordsPerTRM; jj--) {
783 buf[fIndex-nWordsPerTRM] = MakeFiller();
787 //----------------------------------------------------------------------------
788 UInt_t AliTOFDDLRawData::MakeFiller()
791 // Filler word definition: to make even the number of words per TRM/LTM
798 word = 0; // 0000 -> filler ID 1
799 AliBitPacking::PackWord(word,baseWord, 0, 3);
801 AliBitPacking::PackWord(word,baseWord, 4,27);
802 word = 7; // 0111 -> filler ID 2
803 AliBitPacking::PackWord(word,baseWord, 28,31);
809 //----------------------------------------------------------------------------
810 void AliTOFDDLRawData::MakeTRMchainTrailer(Int_t iChain, UInt_t *buf)
816 if (iChain<0 || iChain>1) {
817 AliWarning(Form(" Chain number is out of the right range [0;1] (iChain = %3i", iChain));
826 AliBitPacking::PackWord(word,baseWord, 0, 3);
828 AliBitPacking::PackWord(word,baseWord, 4,15);
829 word = 0; // event counter
830 AliBitPacking::PackWord(word,baseWord,16,27);
833 word = 1; // 0001 -> TRM chain 0 trailer ID
836 word = 3; // 0011 -> TRM chain 1 trailer ID
839 AliBitPacking::PackWord(word,baseWord,28,31);
841 buf[fIndex]=baseWord;
845 //----------------------------------------------------------------------------
846 void AliTOFDDLRawData::MakeTDCdigits(Int_t nDDL, Int_t nTRM, Int_t iChain,
847 UInt_t *buf, UInt_t &nWordsPerTRM)
853 const Double_t kOneMoreFilledCell = 1./(fTOFgeometry->NPadXSector()*fTOFgeometry->NSectors());
854 Double_t percentFilledCells = Double_t(fTOFdigitMap->GetFilledCellNumber())/(fTOFgeometry->NPadXSector()*fTOFgeometry->NSectors());
856 if (nDDL<0 || nDDL>71) {
857 AliWarning(Form(" DDL number is out of the right range [0;71] (nDDL = %3i", nDDL));
861 if (nTRM<3 || nTRM>12) {
862 AliWarning(Form(" TRM number is out of the right range [3;12] (nTRM = %3i", nTRM));
866 if (iChain<0 || iChain>1) {
867 AliWarning(Form(" Chain number is out of the right range [0;1] (iChain = %3i", iChain));
872 UInt_t localBuffer[1000];
873 Int_t localIndex = -1;
875 Int_t iDDL = nDDL%AliTOFGeometry::NDDL();
877 Int_t volume[5] = {-1, -1, -1, -1, -1};
878 Int_t indexDigit[3] = {-1, -1, -1};
880 Int_t totCharge = -1;
881 Int_t timeOfFlight = -1;
883 Int_t trailingSpurious = -1;
884 Int_t leadingSpurious = -1;
897 if (fVerbose==2) ftxt.open("TOFdigits.txt",ios::app);
899 for (jj=0; jj<5; jj++) volume[jj] = -1;
901 // loop on TDC number
902 for (nTDC=AliTOFGeometry::NTdc()-1; nTDC>=0; nTDC--) {
904 // the DRM odd (i.e. left) slot number 3 doesn't contain TDC digit data
905 // for TDC numbers 3-14
906 if (iDDL%2==1 && nTRM==3 && (Int_t)(nTDC/3.)!=0) continue;
908 // loop on TDC channel number
909 for (iCH=AliTOFGeometry::NCh()-1; iCH>=0; iCH--) {
911 fTOFrawStream->EquipmentId2VolumeId(nDDL, nTRM, iChain, nTDC, iCH, volume);
913 if (volume[0]==-1 || volume[1]==-1 || volume[2]==-1 ||
914 volume[3]==-1 || volume[4]==-1) continue;
916 for (jj=0; jj<3; jj++) indexDigit[jj] = -1;
918 fTOFdigitMap->GetDigitIndex(volume, indexDigit);
920 if (indexDigit[0]<0) {
922 trailingSpurious = Int_t(8192*gRandom->Rndm()) + Int_t(Int_t(256*gRandom->Rndm())*AliTOFGeometry::ToTBinWidth()/AliTOFGeometry::TdcBinWidth());
923 leadingSpurious = Int_t(8192*gRandom->Rndm());
925 if ( ( fPackedAcquisition && percentFilledCells<0.12 && gRandom->Rndm()<(0.12-percentFilledCells)) ||
926 (!fPackedAcquisition && percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells)) ) {
928 percentFilledCells+=kOneMoreFilledCell;
933 word = trailingSpurious; // trailing edge measurement
937 word = leadingSpurious; // leading edge measurement
942 if (nDDL<10) ftxt << " " << nDDL;
943 else ftxt << " " << nDDL;
944 if (nTRM<10) ftxt << " " << nTRM;
945 else ftxt << " " << nTRM;
946 ftxt << " " << iChain;
947 if (nTDC<10) ftxt << " " << nTDC;
948 else ftxt << " " << nTDC;
950 if (volume[0]<10) ftxt << " -> " << volume[0];
951 else ftxt << " -> " << volume[0];
952 ftxt << " " << volume[1];
953 if (volume[2]<10) ftxt << " " << volume[2];
954 else ftxt << " " << volume[2];
955 ftxt << " " << volume[4];
956 if (volume[3]<10) ftxt << " " << volume[3];
957 else ftxt << " " << volume[3];
959 if (word<10) ftxt << " " << word;// << endl;
960 else if (word>=10 && word<100) ftxt << " " << word;// << endl;
961 else if (word>=100 && word<1000) ftxt << " " << word;// << endl;
962 else ftxt << " " << word;// << endl;
963 ftxt << " " << dummyPS << endl;
966 AliBitPacking::PackWord(word,baseWord, 0,20);
967 word = iCH; // TDC channel ID [0;7]
968 AliBitPacking::PackWord(word,baseWord,21,23);
969 word = nTDC; // TDC ID [0;14]
970 AliBitPacking::PackWord(word,baseWord,24,27);
971 word = 0; // error flag
972 AliBitPacking::PackWord(word,baseWord,28,28);
973 word = dummyPS; // Packing Status [0;3]
974 AliBitPacking::PackWord(word,baseWord,29,30);
975 word = 1; // TRM TDC digit ID
976 AliBitPacking::PackWord(word,baseWord,31,31);
979 localBuffer[localIndex]=baseWord;
980 psArray[localIndex]=dummyPS;
985 } // if (fPackedAcquisition && percentFilledCells<0.12 && gRandom->Rndm()<(0.12-percentFilledCells)) or ...
987 } // if (indexDigit[0]<0)
989 for (jj=0; jj<3;jj++) {
991 if (indexDigit[jj]<0) continue;
992 digs = (AliTOFdigit*)fTOFdigitArray->UncheckedAt(indexDigit[jj]);
994 if (digs->GetSector()!=volume[0] ||
995 digs->GetPlate() !=volume[1] ||
996 digs->GetStrip() !=volume[2] ||
997 digs->GetPadx() !=volume[3] ||
998 digs->GetPadz() !=volume[4]) AliWarning(" --- ERROR --- ");
1000 //timeOfFlight = (Int_t)(digs->GetTdc())%8192;
1001 timeOfFlight = (Int_t)(digs->GetTdc());
1002 if (timeOfFlight>=8192) timeOfFlight = 0;
1004 totCharge = (Int_t)(digs->GetToT());//digs->GetAdc();
1005 // temporary control
1006 if (totCharge<0) totCharge = 0;//TMath::Abs(totCharge);
1007 if (totCharge>=256) totCharge = 255;
1009 if (fPackedAcquisition) {
1012 if (nDDL<10) ftxt << " " << nDDL;
1013 else ftxt << " " << nDDL;
1014 if (nTRM<10) ftxt << " " << nTRM;
1015 else ftxt << " " << nTRM;
1016 ftxt << " " << iChain;
1017 if (nTDC<10) ftxt << " " << nTDC;
1018 else ftxt << " " << nTDC;
1020 if (volume[0]<10) ftxt << " -> " << volume[0];
1021 else ftxt << " -> " << volume[0];
1022 ftxt << " " << volume[1];
1023 if (volume[2]<10) ftxt << " " << volume[2];
1024 else ftxt << " " << volume[2];
1025 ftxt << " " << volume[4];
1026 if (volume[3]<10) ftxt << " " << volume[3];
1027 else ftxt << " " << volume[3];
1028 if (totCharge<10) ftxt << " " << totCharge;
1029 else if (totCharge>=10 && totCharge<100) ftxt << " " << totCharge;
1030 else ftxt << " " << totCharge;
1031 if (timeOfFlight<10) ftxt << " " << timeOfFlight << endl;
1032 else if (timeOfFlight>=10 && timeOfFlight<100) ftxt << " " << timeOfFlight << endl;
1033 else if (timeOfFlight>=100 && timeOfFlight<1000) ftxt << " " << timeOfFlight << endl;
1034 else ftxt << " " << timeOfFlight << endl;
1037 word = timeOfFlight; // time-of-fligth measurement
1038 AliBitPacking::PackWord(word,baseWord, 0,12);
1040 word = totCharge; // time-over-threshould measurement
1041 AliBitPacking::PackWord(word,baseWord,13,20);
1043 word = iCH; // TDC channel ID [0;7]
1044 AliBitPacking::PackWord(word,baseWord,21,23);
1045 word = nTDC; // TDC ID [0;14]
1046 AliBitPacking::PackWord(word,baseWord,24,27);
1047 word = 0; // error flag
1048 AliBitPacking::PackWord(word,baseWord,28,28);
1049 word = 0; // Packing Status [0;3]
1050 AliBitPacking::PackWord(word,baseWord,29,30);
1051 word = 1; // TRM TDC digit ID
1052 AliBitPacking::PackWord(word,baseWord,31,31);
1055 localBuffer[localIndex]=baseWord;
1060 if (percentFilledCells<0.12 && gRandom->Rndm()<(0.12-percentFilledCells)) {
1062 percentFilledCells+=kOneMoreFilledCell;
1064 trailingSpurious = Int_t(8192*gRandom->Rndm()) + Int_t(Int_t(256*gRandom->Rndm())*AliTOFGeometry::ToTBinWidth()/AliTOFGeometry::TdcBinWidth());
1065 leadingSpurious = Int_t(8192*gRandom->Rndm());
1070 word = trailingSpurious; // trailing edge measurement
1074 word = leadingSpurious; // leading edge measurement
1079 if (nDDL<10) ftxt << " " << nDDL;
1080 else ftxt << " " << nDDL;
1081 if (nTRM<10) ftxt << " " << nTRM;
1082 else ftxt << " " << nTRM;
1083 ftxt << " " << iChain;
1084 if (nTDC<10) ftxt << " " << nTDC;
1085 else ftxt << " " << nTDC;
1087 if (volume[0]<10) ftxt << " -> " << volume[0];
1088 else ftxt << " -> " << volume[0];
1089 ftxt << " " << volume[1];
1090 if (volume[2]<10) ftxt << " " << volume[2];
1091 else ftxt << " " << volume[2];
1092 ftxt << " " << volume[4];
1093 if (volume[3]<10) ftxt << " " << volume[3];
1094 else ftxt << " " << volume[3];
1096 if (word<10) ftxt << " " << word;
1097 else if (word>=10 && word<100) ftxt << " " << word;
1098 else if (word>=100 && word<1000) ftxt << " " << word;
1099 else ftxt << " " << word;
1100 ftxt << " " << dummyPS << endl;
1103 AliBitPacking::PackWord(word,baseWord, 0,20);
1104 word = iCH; // TDC channel ID [0;7]
1105 AliBitPacking::PackWord(word,baseWord,21,23);
1106 word = nTDC; // TDC ID [0;14]
1107 AliBitPacking::PackWord(word,baseWord,24,27);
1108 word = 0; // error flag
1109 AliBitPacking::PackWord(word,baseWord,28,28);
1110 word = dummyPS; // Packing Status [0;3]
1111 AliBitPacking::PackWord(word,baseWord,29,30);
1112 word = 1; // TRM TDC digit ID
1113 AliBitPacking::PackWord(word,baseWord,31,31);
1116 localBuffer[localIndex]=baseWord;
1117 psArray[localIndex]=dummyPS;
1122 } // if (percentFilledCells<0.12 && gRandom->Rndm()<(0.12-percentFilledCells))
1125 } // if (fPackedAcquisition)
1126 else { // if (!fPackedAcquisition)
1128 if (percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells) && HeadOrTail()) {
1130 percentFilledCells+=kOneMoreFilledCell;
1132 trailingSpurious = Int_t(8192*gRandom->Rndm()) + Int_t(Int_t(256*gRandom->Rndm())*AliTOFGeometry::ToTBinWidth()/AliTOFGeometry::TdcBinWidth());
1133 word = trailingSpurious;
1137 if (nDDL<10) ftxt << " " << nDDL;
1138 else ftxt << " " << nDDL;
1139 if (nTRM<10) ftxt << " " << nTRM;
1140 else ftxt << " " << nTRM;
1141 ftxt << " " << iChain;
1142 if (nTDC<10) ftxt << " " << nTDC;
1143 else ftxt << " " << nTDC;
1145 if (volume[0]<10) ftxt << " -> " << volume[0];
1146 else ftxt << " -> " << volume[0];
1147 ftxt << " " << volume[1];
1148 if (volume[2]<10) ftxt << " " << volume[2];
1149 else ftxt << " " << volume[2];
1150 ftxt << " " << volume[4];
1151 if (volume[3]<10) ftxt << " " << volume[3];
1152 else ftxt << " " << volume[3];
1154 if (word<10) ftxt << " " << word;
1155 else if (word>=10 && word<100) ftxt << " " << word;
1156 else if (word>=100 && word<1000) ftxt << " " << word;
1157 else ftxt << " " << word;
1158 ftxt << " " << dummyPS << endl;
1161 AliBitPacking::PackWord(word,baseWord, 0,20);
1162 word = iCH; // TDC channel ID [0;7]
1163 AliBitPacking::PackWord(word,baseWord,21,23);
1164 word = nTDC; // TDC ID [0;14]
1165 AliBitPacking::PackWord(word,baseWord,24,27);
1166 word = 0; // error flag
1167 AliBitPacking::PackWord(word,baseWord,28,28);
1168 word = dummyPS; // Packing Status [0;3]
1169 AliBitPacking::PackWord(word,baseWord,29,30);
1170 word = 1; // TRM TDC digit ID
1171 AliBitPacking::PackWord(word,baseWord,31,31);
1174 localBuffer[localIndex]=baseWord;
1175 psArray[localIndex]=dummyPS;
1180 } // if (percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells))
1183 word = timeOfFlight + Int_t(totCharge*AliTOFGeometry::ToTBinWidth()/AliTOFGeometry::TdcBinWidth()); // trailing edge measurement
1186 if (nDDL<10) ftxt << " " << nDDL;
1187 else ftxt << " " << nDDL;
1188 if (nTRM<10) ftxt << " " << nTRM;
1189 else ftxt << " " << nTRM;
1190 ftxt << " " << iChain;
1191 if (nTDC<10) ftxt << " " << nTDC;
1192 else ftxt << " " << nTDC;
1194 if (volume[0]<10) ftxt << " -> " << volume[0];
1195 else ftxt << " -> " << volume[0];
1196 ftxt << " " << volume[1];
1197 if (volume[2]<10) ftxt << " " << volume[2];
1198 else ftxt << " " << volume[2];
1199 ftxt << " " << volume[4];
1200 if (volume[3]<10) ftxt << " " << volume[3];
1201 else ftxt << " " << volume[3];
1203 if (word<10) ftxt << " " << word;
1204 else if (word>=10 && word<100) ftxt << " " << word;
1205 else if (word>=100 && word<1000) ftxt << " " << word;
1206 else ftxt << " " << word;
1207 ftxt << " " << 2 << endl;
1210 AliBitPacking::PackWord(word,baseWord, 0,20);
1212 word = iCH; // TDC channel ID [0;7]
1213 AliBitPacking::PackWord(word,baseWord,21,23);
1214 word = nTDC; // TDC ID [0;14]
1215 AliBitPacking::PackWord(word,baseWord,24,27);
1216 word = 0; // error flag
1217 AliBitPacking::PackWord(word,baseWord,28,28);
1218 word = 2; // Packing Status [0;3]
1219 AliBitPacking::PackWord(word,baseWord,29,30);
1220 word = 1; // TRM TDC digit ID
1221 AliBitPacking::PackWord(word,baseWord,31,31);
1224 localBuffer[localIndex]=baseWord;
1225 psArray[localIndex]=2;
1232 word = timeOfFlight; // leading edge measurement
1235 if (nDDL<10) ftxt << " " << nDDL;
1236 else ftxt << " " << nDDL;
1237 if (nTRM<10) ftxt << " " << nTRM;
1238 else ftxt << " " << nTRM;
1239 ftxt << " " << iChain;
1240 if (nTDC<10) ftxt << " " << nTDC;
1241 else ftxt << " " << nTDC;
1243 if (volume[0]<10) ftxt << " -> " << volume[0];
1244 else ftxt << " -> " << volume[0];
1245 ftxt << " " << volume[1];
1246 if (volume[2]<10) ftxt << " " << volume[2];
1247 else ftxt << " " << volume[2];
1248 ftxt << " " << volume[4];
1249 if (volume[3]<10) ftxt << " " << volume[3];
1250 else ftxt << " " << volume[3];
1252 if (word<10) ftxt << " " << word;
1253 else if (word>=10 && word<100) ftxt << " " << word;
1254 else ftxt << " " << word;
1255 ftxt << " " << 1 << endl;
1258 AliBitPacking::PackWord(word,baseWord, 0,20);
1260 word = iCH; // TDC channel ID [0;7]
1261 AliBitPacking::PackWord(word,baseWord,21,23);
1262 word = nTDC; // TDC ID [0;14]
1263 AliBitPacking::PackWord(word,baseWord,24,27);
1264 word = 0; // error flag
1265 AliBitPacking::PackWord(word,baseWord,28,28);
1266 word = 1; // Packing Status [0;3]
1267 AliBitPacking::PackWord(word,baseWord,29,30);
1268 word = 1; // TRM TDC digit ID
1269 AliBitPacking::PackWord(word,baseWord,31,31);
1272 localBuffer[localIndex]=baseWord;
1273 psArray[localIndex]=1;
1279 if (percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells) && !HeadOrTail()) {
1281 percentFilledCells+=kOneMoreFilledCell;
1283 leadingSpurious = Int_t(8192*gRandom->Rndm());
1284 word = leadingSpurious;
1288 if (nDDL<10) ftxt << " " << nDDL;
1289 else ftxt << " " << nDDL;
1290 if (nTRM<10) ftxt << " " << nTRM;
1291 else ftxt << " " << nTRM;
1292 ftxt << " " << iChain;
1293 if (nTDC<10) ftxt << " " << nTDC;
1294 else ftxt << " " << nTDC;
1296 if (volume[0]<10) ftxt << " -> " << volume[0];
1297 else ftxt << " -> " << volume[0];
1298 ftxt << " " << volume[1];
1299 if (volume[2]<10) ftxt << " " << volume[2];
1300 else ftxt << " " << volume[2];
1301 ftxt << " " << volume[4];
1302 if (volume[3]<10) ftxt << " " << volume[3];
1303 else ftxt << " " << volume[3];
1305 if (word<10) ftxt << " " << word;// << endl;
1306 else if (word>=10 && word<100) ftxt << " " << word;// << endl;
1307 else if (word>=100 && word<1000) ftxt << " " << word;// << endl;
1308 else ftxt << " " << word;// << endl;
1309 ftxt << " " << dummyPS << endl;
1312 AliBitPacking::PackWord(word,baseWord, 0,20);
1313 word = iCH; // TDC channel ID [0;7]
1314 AliBitPacking::PackWord(word,baseWord,21,23);
1315 word = nTDC; // TDC ID [0;14]
1316 AliBitPacking::PackWord(word,baseWord,24,27);
1317 word = 0; // error flag
1318 AliBitPacking::PackWord(word,baseWord,28,28);
1319 word = dummyPS; // Packing Status [0;3]
1320 AliBitPacking::PackWord(word,baseWord,29,30);
1321 word = 1; // TRM TDC digit ID
1322 AliBitPacking::PackWord(word,baseWord,31,31);
1325 localBuffer[localIndex]=baseWord;
1326 psArray[localIndex]=dummyPS;
1331 } // if (percentFilledCells<0.24 && gRandom->Rndm()<(0.24-percentFilledCells))
1334 } // if (!fPackedAcquisition)
1336 } //end loop on digits in the same volume
1338 } // end loop on TDC channel number
1340 } // end loop on TDC number
1342 if (fPackedAcquisition) {
1344 for (Int_t jj=0; jj<=localIndex; jj++) {
1346 buf[fIndex] = localBuffer[jj];
1352 for (Int_t jj=0; jj<=localIndex; jj++) {
1353 if (psArray[jj]==2) {
1355 buf[fIndex] = localBuffer[jj];
1358 for (Int_t jj=0; jj<=localIndex; jj++) {
1359 if (psArray[jj]==1) {
1361 buf[fIndex] = localBuffer[jj];
1367 if (fVerbose==2) ftxt.close();
1371 //----------------------------------------------------------------------------
1372 void AliTOFDDLRawData::ReverseArray(UInt_t a[], Int_t n) const
1375 // Reverses the n elements of array a
1380 for (ii=0; ii<n/2; ii++) {
1390 //----------------------------------------------------------------------------
1391 Bool_t AliTOFDDLRawData::HeadOrTail() const
1394 // Returns the result of a 'pitch and toss'
1397 Double_t dummy = gRandom->Rndm();
1399 if (dummy<0.5) return kFALSE;