]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0RawReader.cxx
Dependencies on MC truth in reconstruction limited to rec. of MC data for comparison...
[u/mrichter/AliRoot.git] / T0 / AliT0RawReader.cxx
1 #include "AliT0RawReader.h"
2 #include "AliT0Parameters.h"
3 #include "AliBitPacking.h"
4 #include "TBits.h"
5
6 #include <Riostream.h>
7 #include "TMath.h"
8 #include "TH1F.h"
9 #include "TArrayI.h"
10 #include "AliLog.h"
11  
12 ClassImp(AliT0RawReader)
13   
14   AliT0RawReader::AliT0RawReader (AliRawReader *rawReader)
15     :  TTask("T0RawReader","read raw T0 data"),
16        fRawReader(rawReader),
17        fData(NULL),
18        fPosition(0)
19        
20 {
21   //
22 // create an object to read T0raw digits
23   AliDebug(1,"Start ");
24  
25   fRawReader->Reset();
26   fRawReader->Select("T0");
27  
28   cout<<"  AliT0RawReader::AliT0RawReaderfRawReader->Select "<<endl;
29 }
30  AliT0RawReader::~AliT0RawReader ()
31 {
32   // 
33 }
34 /*
35 AliT0RawReader::AliT0RawReader(const AliT0RawReader& o): TTask(o),
36      fRawReader(rawReader),
37        fData(NULL),
38        fPosition(0)
39 {
40   //
41 }
42 */
43
44
45 Bool_t  AliT0RawReader::Next()
46 {
47 // read the next raw digit
48 // returns kFALSE if there is no digit left
49 //"LookUpTable":
50 // Amplitude LED TRM=0; chain=0; TDC 0 -5        channel 0,2,4,6
51 // Time CFD      TRM=0; chain=0; TDC 6 - 11      channel 0,2,4,6
52 // mean time     TRM=0; chain=0; TDC 12          channel 0
53 // T0A           TRM=0; chain=0; TDC 12          channel 2
54 // T0C           TRM=0; chain=0; TDC 12          channel 4
55 // vertex        TRM=0; chain=0; TDC 12          channel 6
56 // mult QTC0        TRM=0; chain=0; TDC 13          channel 0
57 // mult QTC1        TRM=0; chain=0; TDC 13          channel 2
58
59 // Charge QTC0   TRM=1; chain=0; TDC 0 -5        channel 0,2,4,6
60 // Charge QTC1   TRM=1; chain=0; TDC 6 - 11      channel 0,2,4,6
61 // T0A trigger          TRM=1; chain=0; TDC 12          channel 0
62 // T0C trigger          TRM=1; chain=0; TDC 12          channel 2
63 // vertex trigger       TRM=1; chain=0; TDC 12          channel 4
64 // trigger central      TRM=1; chain=0; TDC 13          channel 0
65 // tigger semicenral    TRM=1; chain=0; TDC 13          channel 2
66 //
67 // allData array collect data from all channels in one :
68 // allData[0] - allData[23] 24 CFD channels
69 // allData[24] -   allData[47] 24 LED channels
70 //  allData[48]  mean (T0) signal  
71 // allData[49]   time difference (vertex)
72
73
74  //  if (fDigits == 0x0) fDigits = new AliT0digit(); 
75   // fTree->Branch("T0","AliT0digit",&fDigits,405,1);
76  
77   UInt_t word;
78   Int_t time=0,  itdc=0, ichannel=0; 
79   Int_t numberOfWordsInTRM=0, iTRM=0;
80   Int_t tdcTime, koef,hit, meanTime, timeDiff ;
81   Int_t numTRM=2; // number of TRMs in game For test =1 !!!!!
82
83
84
85   AliT0Parameters* param = AliT0Parameters::Instance();   //-->Zhenya
86   param->Init();
87  
88  for ( Int_t k=0; k<110; k++) {
89     for ( Int_t jj=0; jj<5; jj++) {
90       fAllData[k][jj]=0;
91     }
92   }
93     do {
94     if (!fRawReader->ReadNextData(fData)) return kFALSE;
95   } while (fRawReader->GetDataSize() == 0);
96   
97   //  fPosition = GetPosition();
98   fPosition = 0;
99
100   //DRM header
101   for (Int_t i=0; i<4; i++) {
102     word = GetNextWord();
103   }
104
105   for (Int_t ntrm=0; ntrm< numTRM; ntrm++)
106     {
107       //TRMheader  
108       word = GetNextWord();
109       numberOfWordsInTRM=AliBitPacking::UnpackWord(word,4,16);
110       iTRM=AliBitPacking::UnpackWord(word,0,3);
111       
112       //chain header
113       Int_t ichain=0;
114       word = GetNextWord();
115       
116       for (Int_t i=0; i<numberOfWordsInTRM; i++) {
117         word = GetNextWord();
118         tdcTime =  AliBitPacking::UnpackWord(word,31,31);   
119
120         if ( tdcTime == 1)
121           {
122             itdc=AliBitPacking::UnpackWord(word,24,27);
123             ichannel=AliBitPacking::UnpackWord(word,21,23);
124             time=AliBitPacking::UnpackWord(word,0,20);
125             //  koef = itdc*4 + ichannel/2;
126             koef = param->GetChannel(iTRM,itdc,ichain,ichannel);
127             //   cout<<" RawReader::Next ::"<<iTRM<<"  "<<itdc<<" "<<ichain<<" "<<ichannel<<" "<<  koef<<" "<<time<<endl;
128             if(fAllData[koef][0] == 0)  fAllData[koef][0]=time;  // yield only 1st particle
129             
130           }
131       }
132       word = GetNextWord(); //chain trailer
133       word = GetNextWord(); //TRM trailer
134     } //TRM loop
135    return kTRUE;
136 }
137 //_____________________________________________________________________________
138 Int_t AliT0RawReader::GetPosition()
139 {
140   // Sets the position in the
141   // input stream
142   if (((fRawReader->GetDataSize() * 8) % 32) != 0)
143     AliFatal(Form("Incorrect raw data size ! %d words are found !",fRawReader->GetDataSize()));
144   return (fRawReader->GetDataSize() * 8) / 32;
145 }
146 //_____________________________________________________________________________
147 UInt_t AliT0RawReader::GetNextWord()
148 {
149   // Read the next 32 bit word in backward direction
150   // The input stream access is given by fData and fPosition
151
152
153   //   fPosition--;
154   Int_t iBit = fPosition * 32;
155   Int_t iByte = iBit / 8;
156
157   UInt_t word = 0;
158   word  = fData[iByte+3]<<24;
159   word |= fData[iByte+2]<<16;
160   word |= fData[iByte+1]<<8;
161   word |= fData[iByte];
162    fPosition++;
163
164   return word;
165
166 }
167