]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0RawReader.cxx
Modify assignment operator
[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        fParam(NULL),
20        fIsOnline(kFALSE)
21 {
22   //
23 // create an object to read T0raw digits
24   AliDebug(1,"Start ");
25  
26   fRawReader->Reset();
27   fRawReader->Select("T0");
28  // AliT0Parameters* fParam = AliT0Parameters::Instance();  
29  
30  // fIsOnline = isOnline; 
31   // fIsOnline=kTRUE;
32  /* if (fIsOnline)
33     fParam->InitIfOnline();
34   else 
35     fParam->Init();
36  */
37 }
38  AliT0RawReader::~AliT0RawReader ()
39 {
40   // 
41 }
42 /*
43 AliT0RawReader::AliT0RawReader(const AliT0RawReader& o): TTask(o),
44      fRawReader(rawReader),
45        fData(NULL),
46        fPosition(0)
47 {
48   //
49 }
50 */
51
52
53 Bool_t  AliT0RawReader::Next()
54 {
55 // read the next raw digit
56 // returns kFALSE if there is no digit left
57 //"LookUpTable":
58 // Amplitude LED TRM=0; chain=0; TDC 0 -5        channel 0,2,4,6
59 // Time CFD      TRM=0; chain=0; TDC 6 - 11      channel 0,2,4,6
60 // mean time     TRM=0; chain=0; TDC 12          channel 0
61 // T0A           TRM=0; chain=0; TDC 12          channel 2
62 // T0C           TRM=0; chain=0; TDC 12          channel 4
63 // vertex        TRM=0; chain=0; TDC 12          channel 6
64 // mult QTC0        TRM=0; chain=0; TDC 13          channel 0
65 // mult QTC1        TRM=0; chain=0; TDC 13          channel 2
66
67 // Charge QTC0   TRM=1; chain=0; TDC 0 -5        channel 0,2,4,6
68 // Charge QTC1   TRM=1; chain=0; TDC 6 - 11      channel 0,2,4,6
69 // T0A trigger          TRM=1; chain=0; TDC 12          channel 0
70 // T0C trigger          TRM=1; chain=0; TDC 12          channel 2
71 // vertex trigger       TRM=1; chain=0; TDC 12          channel 4
72 // trigger central      TRM=1; chain=0; TDC 13          channel 0
73 // tigger semicenral    TRM=1; chain=0; TDC 13          channel 2
74 //
75 // allData array collect data from all channels in one :
76 // allData[0] - allData[23] 24 CFD channels
77 // allData[24] -   allData[47] 24 LED channels
78 //  allData[48]  mean (T0) signal  
79 // allData[49]   time difference (vertex)
80
81  
82   UInt_t word;
83   Int_t time=0,  itdc=0, ichannel=0, uu; 
84   Int_t numberOfWordsInTRM=0, iTRM=0;
85   Int_t tdcTime, koef,hit=0;
86   Int_t koefhits[110];
87   //  Int_t fDRM_GLOBAL_HEADER = 0x40000001;
88   //  Int_t fDRM_GLOBAL_TRAILER  = 0x50000001;
89   
90   //  Int_t  TRM_GLOBAL_HEADER  = 0x40000000;
91     Int_t  TRM_CHAIN_0_HEADER =  0x00000000;
92     //  Int_t  TRM_CHAIN_1_HEADER =  0x20000000;
93   Int_t  TRM_CHAIN_0_TRAILER =  0x10000000;
94   //  Int_t  TRM_CHAIN_1_TRAILER =  0x30000000;
95   //  Int_t  TRM_GLOBAL_TRAILER =  0x5000000f;
96   
97   Int_t  FILLER =  0x70000000;
98   Bool_t correct=kTRUE;
99   Int_t header;
100   
101     AliT0Parameters* fParam = AliT0Parameters::Instance();   
102    if (fIsOnline)
103    fParam->InitIfOnline();
104         else
105    fParam->Init();
106
107    Int_t fNTRM = fParam->GetNumberOfTRMs();
108
109    for ( Int_t k=0; k<110; k++) {
110     koefhits[k]=0;
111     for ( Int_t jj=0; jj<5; jj++) {
112       fAllData[k][jj]=0;
113      }
114    }
115     do {
116       if (!fRawReader->ReadNextData(fData)) return kFALSE;
117     } while (fRawReader->GetDataSize() == 0);
118     
119     //  fPosition = GetPosition();
120     fPosition = 0;
121     cout.setf( ios_base::hex, ios_base::basefield );
122     
123     //DRM header
124     for (Int_t i=0; i<6; i++) {
125       word = GetNextWord();
126       //      cout<<" DRM header "<<i<<" "<<word<<endl;
127       header = AliBitPacking::UnpackWord(word,28,31);
128       if( header !=4 )
129         //      uu = word&fDRM_GLOBAL_HEADER;
130         //      if(uu != fDRM_GLOBAL_HEADER ) 
131         {
132           AliError(Form(" !!!! wrong  DRM header  %x!!!!", word));
133           break;
134         }
135     }
136     
137     for (Int_t ntrm=0; ntrm< fNTRM; ntrm++)
138       {
139         //TRMheader  
140         word = GetNextWord();
141         //      cout<<" TRM header "<<word<<" "<<AliBitPacking::UnpackWord(word,28,31)<<endl;
142       header = AliBitPacking::UnpackWord(word,28,31);
143         
144       //        uu = word&TRM_GLOBAL_HEADER;
145         //      if(uu != TRM_GLOBAL_HEADER &&
146           if ( header != 4 )
147         {
148           AliError(Form(" !!!! wrong TRM header  %x!!!!", word));
149           break;
150         }
151         numberOfWordsInTRM=AliBitPacking::UnpackWord(word,4,16);
152         iTRM=AliBitPacking::UnpackWord(word,0,3);
153         //      cout<<" iTRM "<<iTRM<<" Nworrds in TRM "<< numberOfWordsInTRM  <<endl;   
154
155         for( Int_t ichain=0; ichain<2; ichain++)
156           {
157             //chain header
158             word = GetNextWord();
159             uu = word & TRM_CHAIN_0_HEADER;
160             if(uu != TRM_CHAIN_0_HEADER) 
161               {
162                 AliError(Form(" !!!! wrong CHAIN  0  header %x!!!!", word));
163                 break;
164               }
165             word = GetNextWord();
166             //      cout<<" ???data?? "<<word<<endl;
167             tdcTime =  AliBitPacking::UnpackWord(word,31,31);   
168             for (; tdcTime==1; tdcTime) 
169               {
170                 itdc=AliBitPacking::UnpackWord(word,24,27);
171                 ichannel=AliBitPacking::UnpackWord(word,21,23);
172                 time=AliBitPacking::UnpackWord(word,0,20);
173                 
174                 koef = fParam->GetChannel(iTRM,itdc,ichain,ichannel);
175                 //              cout<<"koef "<<koef<<" iTRM "<<iTRM<<" itdc "<<itdc<<" iChain "<<ichain<<" ichannel  "<<ichannel<<" "<<time<<endl;
176                 if (koef ==-1 ){
177                   AliError(Form("Incorrect lookup table ! ")); 
178                   correct=kFALSE;
179                 }
180                 if(correct){
181                   hit=koefhits[koef];
182                   //if(fAllData[koef][hit] == 0)  
183                   fAllData[koef][hit]=time; 
184                   koefhits[koef]++;
185                 }
186                   word = GetNextWord();
187                   //            cout<<" next word "<<word<<endl;
188                   tdcTime =  AliBitPacking::UnpackWord(word,31,31);   
189                 
190             }
191             
192
193             uu = word&TRM_CHAIN_0_TRAILER;
194             if(uu != TRM_CHAIN_0_TRAILER )
195               {
196                 AliError(Form(" !!!! wrong CHAIN 0 trailer %x !!!!", word));
197                 break;
198               }
199           }
200             
201         word = GetNextWord(); //TRM trailer
202         //      cout<<" TRM traier "<<word<<endl;
203         header = AliBitPacking::UnpackWord(word,28,31);
204         if( header !=5 )
205           {
206             AliError(Form(" !!!! wrong TRM GLOBAL trailer  %x!!!!", word));
207             break;
208           }
209       } //TRM loop
210     word = GetNextWord(); //
211        header = AliBitPacking::UnpackWord(word,28,31);
212     uu = word&FILLER;
213     if (word == FILLER )  word = GetNextWord(); 
214     //          cout<<" word after FILLER "<<word<<endl;
215     if( header !=5 )
216       {
217         AliError(Form(" !!!! wrong DRM GLOBAL trailer  %x!!!!", word));
218           //      break;
219       }
220     cout.setf( ios_base::dec, ios_base::basefield );
221     
222     return kTRUE;
223 }
224 //_____________________________________________________________________________
225 Int_t AliT0RawReader::GetPosition()
226 {
227   // Sets the position in the
228   // input stream
229   if (((fRawReader->GetDataSize() * 8) % 32) != 0)
230     AliFatal(Form("Incorrect raw data size ! %d words are found !",fRawReader->GetDataSize()));
231   return (fRawReader->GetDataSize() * 8) / 32;
232 }
233 //_____________________________________________________________________________
234 UInt_t AliT0RawReader::GetNextWord()
235 {
236   // Read the next 32 bit word in backward direction
237   // The input stream access is given by fData and fPosition
238
239
240   //   fPosition--;
241   Int_t iBit = fPosition * 32;
242   Int_t iByte = iBit / 8;
243
244   UInt_t word = 0;
245   word  = fData[iByte+3]<<24;
246   word |= fData[iByte+2]<<16;
247   word |= fData[iByte+1]<<8;
248   word |= fData[iByte];
249    fPosition++;
250
251   return word;
252
253 }
254