]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSRawStreamSDD.cxx
adding default constructor that does not need a name
[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 fEndWords(0),
70 fResetSkip(0)
71 {
72 // create an object to read ITS SDD raw digits
73
74   Reset();
75   for(Int_t i=0;i<kFifoWords;i++) fNfifo[i]=0;
76   for(Int_t i=0;i<kDDLsNumber;i++) fSkip[i]=0;
77   fRawReader->Reset();
78   fRawReader->Select("ITSSDD");
79
80   for(Short_t i=0; i<kCarlosWords; i++) fICarlosWord[i]=0x30000000 + i; // 805306368+i;
81   for(Short_t i=0; i<kFifoWords; i++) fIFifoWord[i]=0x30000010 + i;  // 805306384+i;
82 }
83
84 UInt_t AliITSRawStreamSDD::ReadBits()
85 {
86 // read bits from the given channel
87   UInt_t result = (fChannelData[fCarlosId][fChannel] & ((1<<fReadBits[fCarlosId][fChannel]) - 1));
88   fChannelData[fCarlosId][fChannel] >>= fReadBits[fCarlosId][fChannel]; 
89   fLastBit[fCarlosId][fChannel] -= fReadBits[fCarlosId][fChannel];
90   return result;
91 }
92
93 Int_t AliITSRawStreamSDD::DecompAmbra(Int_t value) const
94 {
95   // AMBRA decompression (from 8 to 10 bit)
96   
97   if ((value & 0x80) == 0) {
98     return value & 0x7f;
99   } else if ((value & 0x40) == 0) {
100     return 0x081 + ((value & 0x3f) << 1);
101   } else if ((value & 0x20) == 0) {
102     return 0x104 + ((value & 0x1f) << 3);
103   } else {
104     return 0x208 + ((value & 0x1f) << 4);
105   }
106   
107 }
108
109 Bool_t AliITSRawStreamSDD::Next()
110 {
111 // read the next raw digit
112 // returns kFALSE if there is no digit left
113 // returns kTRUE and fCompletedModule=kFALSE when a digit is found
114 // returns kTRUE and fCompletedModule=kTRUE  when a module is completed (=3x3FFFFFFF footer words)
115
116   fPrevModuleID = fModuleID;
117   fDDL=fRawReader->GetDDLID();
118   //cout << "fDDL: " << fDDL;
119   Int_t ddln = fRawReader->GetDDLID();
120   //cout << ", ddln: " << ddln << endl;
121   if(ddln <0) ddln=0;
122   fCompletedModule=kFALSE;
123
124   while (kTRUE) {
125     if(fResetSkip==0){
126       Bool_t kSkip = ResetSkip(ddln);
127       fResetSkip=1;
128       if(!kSkip) return kSkip;
129     }
130   
131     if ((fChannel < 0) || (fLastBit[fCarlosId][fChannel] < fReadBits[fCarlosId][fChannel])) {
132       if (!fRawReader->ReadNextInt(fData)) {
133         //cout << "read word in Next and skip, fData: ";
134         //printf("%x\n",fData);
135         return kFALSE;  // read next word
136       }
137       //cout << "read word in Next, fData: ";
138       //printf("%x\n",fData);
139
140       ddln = fRawReader->GetDDLID();
141       if(ddln!=fDDL) { 
142         Reset();
143         fChannel=-1;
144         fDDL=fRawReader->GetDDLID();
145       }
146       if(ddln < 0 || ddln > (kDDLsNumber-1)) ddln  = 0;
147       //cout << "in the while loop, fDDL: " << fDDL << ", ddln: " << ddln << endl;
148
149       fChannel = -1;
150       if(fData>=fICarlosWord[0]&&fData<=fICarlosWord[11]) { // Carlos Word
151         if(fEndWords==12) continue; // out of event
152         else if(fEndWords<12){
153           fCarlosId = fData-fICarlosWord[0];
154           Int_t iFifoIdx = fCarlosId/3;
155           if(fNCarlos == 8) {
156             if(fCarlosId==2) iFifoIdx = 1;
157             if(fCarlosId==4 || fCarlosId==5)  iFifoIdx = 2;
158             if(fCarlosId==6 || fCarlosId==7)  iFifoIdx = 3 ;
159           }         
160           fNfifo[iFifoIdx] = fCarlosId;
161           //cout <<  "set fCarlosId to " << fCarlosId << " and fNfifo[" << iFifoIdx << "] to " << fNfifo[iFifoIdx] << endl;
162         }
163       } else if (fData>=fIFifoWord[0]&&fData<=fIFifoWord[3]){
164         if(fEndWords==12) continue; // out of event
165         else if(fEndWords<12){
166           //cout << "fData-fIFifoWord[0]: " << fData-fIFifoWord[0] << endl;
167           fCarlosId = fNfifo[fData-fIFifoWord[0]];          
168           //cout << "fCarlosId set to " << fCarlosId << " from FIFO Word" << endl;
169         }
170       }
171       if((fData >> 4) == 0xFF00000){   // jitter word (to be improved !!!!)
172         for(Int_t i=0;i<kDDLsNumber;i++){fSkip[i]=0;}
173         fResetSkip=0;
174         fEndWords=0;
175         continue;
176       }
177
178       if(fData==0x3F1F1F1F){
179         fEndWords++;
180         if(fEndWords<=12) continue;
181       }
182       if (fNCarlos == 8 && (fCarlosId == 8 || fCarlosId == 9 || 
183                             fCarlosId ==10 || fCarlosId == 11)) continue;  // old data, fNCarlos = 8;
184             
185       //cout << "set module from DDLID " << fRawReader->GetDDLID() << " and fCarlosId: " << fCarlosId << endl;
186       fModuleID = fgkDDLModuleMap[fRawReader->GetDDLID()][fCarlosId];     
187       //cout<<"AliITSRawStreamSDD: fModuleID: "<<fModuleID<<endl;
188       
189       if ((fData >> 28) == 0x02) {           // header
190         fEventId = (fData >> 3) & 0x07FF;
191         //cout<<"AliITSRawStreamSDD:fEventID: "<<fEventId<<endl;
192       } else if ((fData >> 28)== 0x03) {    // footer
193         if(fData==0x3FFFFFFF){
194           fICountFoot[fCarlosId]++; // stop before the last word (last word=jitter)
195           if(fICountFoot[fCarlosId]==3){
196             fCompletedModule=kTRUE;
197             return kTRUE;
198           }
199         }               
200       } else if ((fData >> 29) == 0x00) {    // error
201         if ((fData & 0x00000163) != 0) {         
202           fRawReader->AddMajorErrorLog(kDataError,Form("Error code = %8.8x",fData));     
203           AliWarning(Form("error codes = %8.8x",fData));
204           return kFALSE;          
205         }
206       } else if ((fData >> 30) == 0x01) {    // JTAG word
207         // ignored
208       } else if ((fData >> 30) == 0x02) {    // channel 0 data
209         fChannel = 0;
210         //cout << "fChannel set to " << fChannel << endl;
211       } else if ((fData >> 30) == 0x03) {    // channel 1 data
212         fChannel = 1;
213         //cout << "fChannel set to " << fChannel << endl;
214       } else {                               // unknown data format
215         fRawReader->AddMajorErrorLog(kDataFormatErr,Form("Invalid data %8.8x",fData));
216         AliWarning(Form("invalid data: %8.8x\n", fData));
217         return kFALSE;
218       }
219       
220       if (fChannel >= 0) {   // add read word to the data
221         //cout << "add read word to the data" << endl;
222         fChannelData[fCarlosId][fChannel] += 
223           (ULong64_t(fData & 0x3FFFFFFF) << fLastBit[fCarlosId][fChannel]);
224         fLastBit[fCarlosId][fChannel] += 30;
225       }
226     } else {  // decode data
227       if (fReadCode[fCarlosId][fChannel]) {// read the next code word
228         fChannelCode[fCarlosId][fChannel] = ReadBits();
229         fReadCode[fCarlosId][fChannel] = kFALSE;
230         fReadBits[fCarlosId][fChannel] = fgkCodeLength[fChannelCode[fCarlosId][fChannel]];
231       } else {                      // read the next data word
232         UInt_t data = ReadBits();
233         fReadCode[fCarlosId][fChannel] = kTRUE;
234         fReadBits[fCarlosId][fChannel] = 3;
235         if (fChannelCode[fCarlosId][fChannel] == 0) {         // set the time bin         
236           fTimeBin[fCarlosId][fChannel] = data;
237         } else if (fChannelCode[fCarlosId][fChannel] == 1) {  // next anode
238           fTimeBin[fCarlosId][fChannel] = 0;
239           fAnode[fCarlosId][fChannel]++;
240         } else {                                   // ADC signal data
241           fSignal = DecompAmbra(data + (1 << fChannelCode[fCarlosId][fChannel]) + fLowThreshold[fChannel]);
242           fCoord1 = fAnode[fCarlosId][fChannel];
243           fCoord2 = fTimeBin[fCarlosId][fChannel];
244           fTimeBin[fCarlosId][fChannel]++;
245           //cout << "data read, Module, Anode, Time, Charge = " << fModuleID << "," << fCoord1 << "," << fCoord2 << "," << fSignal << endl;
246           return kTRUE;
247         }
248       }
249     }
250   }
251   return kFALSE;  
252
253 }
254
255 void AliITSRawStreamSDD::Reset(){
256
257   //reset data member for a new ddl
258   for(Int_t i=0;i<2;i++){
259     for(Int_t ic=0;ic<kModulesPerDDL;ic++){
260       fChannelData[ic][i]=0;
261       fLastBit[ic][i]=0;
262       fChannelCode[ic][i]=0;
263       fReadCode[ic][i]=kFALSE;
264       fReadBits[ic][i]=0;
265       fTimeBin[ic][i]=0;
266       fAnode[ic][i]=0;     
267     }
268     fLowThreshold[i]=0;
269   }
270   for(Int_t i=0;i<kModulesPerDDL;i++) fICountFoot[i]=0;
271 }
272
273 Bool_t AliITSRawStreamSDD::ResetSkip(Int_t ddln){
274   Bool_t startCount=kFALSE;
275   while (fSkip[ddln] < 9) {
276     if (!fRawReader->ReadNextInt(fData)) { 
277       //cout << "read word in ResetSkip, fData: ";
278       //printf("%x\n",fData);
279       //cout << "return kFALSE in resetskip" << endl; 
280       return kFALSE;
281     }
282     if(fData==0xFFFFFFFF) startCount=kTRUE;
283     //cout << "read word in ResetSkip, fData: ";
284     //printf("%x\n",fData);
285     if ((fData >> 30) == 0x01) continue;  // JTAG word
286     if(startCount) fSkip[ddln]++;
287   }
288   return kTRUE;
289 }
290