]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSRawStreamSDD.cxx
Update of PID 2 to the present recpoint normalisation
[u/mrichter/AliRoot.git] / ITS / AliITSRawStreamSDD.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, 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 /* $Id$*/
17
18 ///////////////////////////////////////////////////////////////////////////////
19 ///
20 /// This class provides access to ITS SDD digits in raw data.
21 ///
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include "AliITSRawStreamSDD.h"
25 #include "AliRawReader.h"
26 #include "AliLog.h"
27 #include "Riostream.h"
28
29 ClassImp(AliITSRawStreamSDD)
30
31 const Int_t AliITSRawStreamSDD::fgkDDLModuleMap[kDDLsNumber][kModulesPerDDL] = {
32  
33   {240,241,242,246,247,248,252,253,254,258,259,260},
34   {264,265,266,270,271,272,276,277,278,282,283,284},
35   {288,289,290,294,295,296,300,301,302,306,307,308},
36   {312,313,314,318,319,320,-1,-1,-1,-1,-1,-1},
37   {243,244,245,249,250,251,255,256,257,261,262,263},
38   {267,268,269,273,274,275,279,280,281,285,286,287},
39   {291,292,293,297,298,299,303,304,305,309,310,311},
40   {315,316,317,321,322,323,-1,-1,-1,-1,-1,-1},
41   {324,325,326,327,332,333,334,335,340,341,342,343},
42   {348,349,350,351,356,357,358,359,364,365,366,367},
43   {372,373,374,375,380,381,382,383,388,389,390,391},
44   {396,397,398,399,404,405,406,407,412,413,414,415},
45   {420,421,422,423,428,429,430,431,436,437,438,439},
46   {444,445,446,447,452,453,454,455,460,461,462,463},
47   {468,469,470,471,476,477,478,479,484,485,486,487},
48   {492,493,494,495,-1,-1,-1,-1,-1,-1,-1,-1},
49   {328,329,330,331,336,337,338,339,344,345,346,347},
50   {352,353,354,355,360,361,362,363,368,369,370,371},
51   {376,377,378,379,384,385,386,387,392,393,394,395},
52   {400,401,402,403,408,409,410,411,416,417,418,419},
53   {424,425,426,427,432,433,434,435,440,441,442,443},
54   {448,449,450,451,456,457,458,459,464,465,466,467},
55   {472,473,474,475,480,481,482,483,488,489,490,491},
56   {496,497,498,499,-1,-1,-1,-1,-1,-1,-1,-1}};
57   
58 const UInt_t AliITSRawStreamSDD::fgkCodeLength[8] =  {8, 18, 2, 3, 4, 5, 6, 7};
59
60 AliITSRawStreamSDD::AliITSRawStreamSDD(AliRawReader* rawReader) :
61   AliITSRawStream(rawReader),
62 fData(0),
63 fCarlosId(-1),
64 fEventId(0),
65 fChannel(0),
66 fJitter(0),
67 fNCarlos(kModulesPerDDL),
68 fDDL(0),
69 fIdcd(0),
70 fEndWords(0),
71 fResetSkip(0)
72 {
73 // create an object to read ITS SDD raw digits
74
75   Reset();
76   for(Int_t i=0;i<kFifoWords;i++) fNfifo[i]=0;
77   for(Int_t i=0;i<kDDLsNumber;i++) fSkip[i]=0;
78   fRawReader->Reset();
79   fRawReader->Select("ITSSDD");
80
81   for(Short_t i=0; i<kCarlosWords; i++) fICarlosWord[i]=0x30000000 + i; // 805306368+i;
82   for(Short_t i=0; i<kFifoWords; i++) fIFifoWord[i]=0x30000010 + i;  // 805306384+i;
83 }
84
85 UInt_t AliITSRawStreamSDD::ReadBits()
86 {
87 // read bits from the given channel
88   UInt_t result = (fChannelData[fCarlosId][fChannel] & ((1<<fReadBits[fCarlosId][fChannel]) - 1));
89   fChannelData[fCarlosId][fChannel] >>= fReadBits[fCarlosId][fChannel]; 
90   fLastBit[fCarlosId][fChannel] -= fReadBits[fCarlosId][fChannel];
91   return result;
92 }
93
94 Int_t AliITSRawStreamSDD::DecompAmbra(Int_t value) const
95 {
96   // AMBRA decompression (from 8 to 10 bit)
97   
98   if ((value & 0x80) == 0) {
99     return value & 0x7f;
100   } else if ((value & 0x40) == 0) {
101     return 0x081 + ((value & 0x3f) << 1);
102   } else if ((value & 0x20) == 0) {
103     return 0x104 + ((value & 0x1f) << 3);
104   } else {
105     return 0x208 + ((value & 0x1f) << 4);
106   }
107   
108 }
109
110 Bool_t AliITSRawStreamSDD::Next()
111 {
112 // read the next raw digit
113 // returns kFALSE if there is no digit left
114   fPrevModuleID = fModuleID;
115   fDDL=fRawReader->GetDDLID();
116   //cout << "fDDL: " << fDDL;
117   Int_t ddln = fRawReader->GetDDLID();
118   //cout << ", ddln: " << ddln << endl;
119   if(ddln <0) ddln=0;
120   if(fResetSkip==0){
121     Bool_t kSkip = ResetSkip(ddln);
122     fResetSkip=1;
123     if(!kSkip) return kSkip;
124   }
125
126
127   for(Int_t i=0;i<kModulesPerDDL;i++) fICountFoot[i]=0; 
128
129   while (kTRUE) {
130   
131     if ((fChannel < 0) || (fLastBit[fCarlosId][fChannel] < fReadBits[fCarlosId][fChannel])) {
132       //cout << "fCarlosId: " << fCarlosId << ", fChannel: " << fChannel << ", fLastBit[fCarlosId][fChannel]: " << fLastBit[fCarlosId][fChannel] << ", fReadBits[fCarlosId][fChannel]: " << fReadBits[fCarlosId][fChannel] << endl;
133       if (!fRawReader->ReadNextInt(fData)) {
134         //cout << "read word in Next and skip, fData: ";
135         //printf("%x\n",fData);
136         return kFALSE;  // read next word
137       }
138       //cout << "read word in Next, fData: ";
139       //printf("%x\n",fData);
140
141       ddln = fRawReader->GetDDLID();
142       if(ddln!=fDDL) { 
143         Reset();
144         fChannel=-1;
145         fDDL=fRawReader->GetDDLID();
146       }
147       if(ddln < 0 || ddln > (kDDLsNumber-1)) ddln  = 0;
148       //cout << "in the while loop, fDDL: " << fDDL << ", ddln: " << ddln << endl;
149 /*
150       kSkip = ResetSkip(ddln);
151           cout << "kSkip: " << kSkip << endl;
152       if(!kSkip) return kFALSE;  // read next word
153 */
154       fChannel = -1;
155       if(fData>=fICarlosWord[0]&&fData<=fICarlosWord[11]) { // Carlos Word
156         if(fEndWords==12) continue; // out of event
157         else if(fEndWords<12){
158           fCarlosId = fData-fICarlosWord[0];
159           //cout << "set fCarlosId to " << fCarlosId;
160           Int_t iFifoIdx = fCarlosId/3;
161           if(fNCarlos == 8) {
162             if(fCarlosId==2) iFifoIdx = 1;
163             if(fCarlosId==4 || fCarlosId==5)  iFifoIdx = 2;
164             if(fCarlosId==6 || fCarlosId==7)  iFifoIdx = 3 ;
165           }         
166           fNfifo[iFifoIdx] = fCarlosId;
167           //cout << " and fNfifo[" << iFifoIdx << "] to " << fNfifo[iFifoIdx] << endl;
168         }
169       } else if (fData>=fIFifoWord[0]&&fData<=fIFifoWord[3]){
170         //cout << "fIdcd: " << fIdcd << endl;
171         fIdcd=0;
172         if(fEndWords==12) continue; // out of event
173         else if(fEndWords<12){
174           //cout << "fData-fIFifoWord[0]: " << fData-fIFifoWord[0] << endl;
175           fCarlosId = fNfifo[fData-fIFifoWord[0]];          
176           //cout << "fCarlosId set to " << fCarlosId << " from FIFO Word" << endl;
177         }
178       }
179           
180       if((fData >> 4) == 0xFF00000){   // jitter word
181         for(Int_t i=0;i<kDDLsNumber;i++){fSkip[i]=0;}
182         fResetSkip=0;
183         fEndWords=0;
184         return kTRUE;
185       }
186       if(fData==0x3F1F1F1F){
187         fEndWords++;
188         if(fEndWords<=12) continue;
189       }
190       if (fNCarlos == 8 && (fCarlosId == 8 || fCarlosId == 9 || 
191                             fCarlosId ==10 || fCarlosId == 11)) continue;  // old data, fNCarlos = 8;
192             
193       //cout << "set module from DDLID " << fRawReader->GetDDLID() << " and fCarlosId: " << fCarlosId << endl;
194       fModuleID = fgkDDLModuleMap[fRawReader->GetDDLID()-1][fCarlosId];
195       //cout<<"AliITSRawStreamSDD: fModuleID: "<<fModuleID<<endl;
196       
197       if ((fData >> 28) == 0x02) {           // header
198         fEventId = (fData >> 3) & 0x07FF;
199         //cout<<"AliITSRawStreamSDD:fEventID: "<<fEventId<<endl;
200       } else if ((fData >> 28)== 0x03) {    // footer
201         if((fData>=fIFifoWord[0]&&fData<=fIFifoWord[3])||(fData>=fICarlosWord[0]&&fData<=fICarlosWord[11]) )            {
202           // Carlos and fifo words should not be counted as footers
203         }else{   // footer word
204           fICountFoot[fCarlosId]++; // stop before the last word (last word=jitter)
205           Bool_t exit=kTRUE;
206           for(Int_t ic=0;ic<fNCarlos;ic++){
207             if(fICountFoot[ic]<3) exit=kFALSE;
208           }
209           if(exit) return kFALSE;
210         }               
211       } else if ((fData >> 29) == 0x00) {    // error
212         if ((fData & 0x00000163) != 0) {         
213           fRawReader->AddMajorErrorLog(kDataError,Form("Error code = %8.8x",fData));     
214           AliWarning(Form("error codes = %8.8x",fData));
215           return kFALSE;          
216         }
217       } else if ((fData >> 30) == 0x01) {    // JTAG word
218         // ignored
219       } else if ((fData >> 30) == 0x02) {    // channel 0 data
220         fChannel = 0;
221         //cout << "fChannel set to " << fChannel << endl;
222       } else if ((fData >> 30) == 0x03) {    // channel 1 data
223         fChannel = 1;
224         //cout << "fChannel set to " << fChannel << endl;
225       } else {                               // unknown data format
226         fRawReader->AddMajorErrorLog(kDataFormatErr,Form("Invalid data %8.8x",fData));
227         AliWarning(Form("invalid data: %8.8x\n", fData));
228         return kFALSE;
229       }
230       
231       if (fChannel >= 0) {   // add read word to the data
232         //cout << "add read word to the data" << endl;
233         fChannelData[fCarlosId][fChannel] += 
234           (ULong64_t(fData & 0x3FFFFFFF) << fLastBit[fCarlosId][fChannel]);
235         fLastBit[fCarlosId][fChannel] += 30;
236       }
237     } else {  // decode data
238       fIdcd++;
239       if (fReadCode[fCarlosId][fChannel]) {// read the next code word
240         fChannelCode[fCarlosId][fChannel] = ReadBits();
241         fReadCode[fCarlosId][fChannel] = kFALSE;
242         fReadBits[fCarlosId][fChannel] = fgkCodeLength[fChannelCode[fCarlosId][fChannel]];
243       } else {                      // read the next data word
244         UInt_t data = ReadBits();
245         fReadCode[fCarlosId][fChannel] = kTRUE;
246         fReadBits[fCarlosId][fChannel] = 3;
247         if (fChannelCode[fCarlosId][fChannel] == 0) {         // set the time bin         
248           fTimeBin[fCarlosId][fChannel] = data;
249         } else if (fChannelCode[fCarlosId][fChannel] == 1) {  // next anode
250           fTimeBin[fCarlosId][fChannel] = 0;
251           fAnode[fCarlosId][fChannel]++;
252         } else {                                   // ADC signal data
253           fSignal = DecompAmbra(data + (1 << fChannelCode[fCarlosId][fChannel]) + fLowThreshold[fChannel]);
254           fCoord1 = fAnode[fCarlosId][fChannel];
255           fCoord2 = fTimeBin[fCarlosId][fChannel];
256           fTimeBin[fCarlosId][fChannel]++;
257           //cout << "data read, Module, Anode, Time, Charge = " << fModuleID << "," << fCoord1 << "," << fCoord2 << "," << fSignal << endl;
258           return kTRUE;
259         }
260       }
261     }
262   }
263   return kFALSE;  
264
265 }
266
267 void AliITSRawStreamSDD::Reset(){
268
269   //reset data member for a new ddl
270   for(Int_t i=0;i<2;i++){
271     for(Int_t ic=0;ic<kModulesPerDDL;ic++){
272       fChannelData[ic][i]=0;
273       fLastBit[ic][i]=0;
274       fChannelCode[ic][i]=0;
275       fReadCode[ic][i]=kFALSE;
276       fReadBits[ic][i]=0;
277       fTimeBin[ic][i]=0;
278       fAnode[ic][i]=0;     
279     }
280     fLowThreshold[i]=0;
281   }
282   for(Int_t i=0;i<kModulesPerDDL;i++) fICountFoot[i]=0;
283 }
284
285 Bool_t AliITSRawStreamSDD::ResetSkip(Int_t ddln){
286   while (fSkip[ddln] < 9) {
287     if (!fRawReader->ReadNextInt(fData)) { 
288       //cout << "read word in ResetSkip, fData: ";
289       //printf("%x\n",fData);
290       //cout << "return kFALSE in resetskip" << endl; 
291       return kFALSE;
292     }
293     //cout << "read word in ResetSkip, fData: ";
294     //printf("%x\n",fData);
295     if ((fData >> 30) == 0x01) continue;  // JTAG word
296     fSkip[ddln]++;
297   }
298   return kTRUE;
299 }
300