]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSRawStreamSDD.cxx
SDD Quality Assurance protoype macro - SDD RawStreamer code changed in order to be...
[u/mrichter/AliRoot.git] / ITS / AliITSRawStreamSDD.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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(0),
64 fEventId(0),
65 fChannel(0),
66 fJitter(0),
67 fNCarlos(kModulesPerDDL),
68 fDDL(0){
69 // create an object to read ITS SDD raw digits
70   
71   for(Int_t i=0;i<2;i++){
72     for(Int_t ic=0;ic<kModulesPerDDL;ic++){
73       fChannelData[ic][i]=0;
74       fLastBit[ic][i]=0;
75       fChannelCode[ic][i]=0;
76       fReadCode[ic][i]=kFALSE;
77       fReadBits[ic][i]=0;
78       fTimeBin[ic][i]=0;
79       fAnode[ic][i]=0;     
80     }
81     fLowThreshold[i]=0;
82   }
83   for(Int_t i=0;i<kFifoWords;i++) fNfifo[i]=0;
84   for(Int_t i=0;i<kDDLsNumber;i++) fSkip[i]=0;
85   fRawReader->Reset();
86   fRawReader->Select("ITSSDD");
87
88   //  fNCarlos = 8; //select this for test data
89   for(Short_t i=0; i<kCarlosWords; i++) iCarlosWord[i]=805306368+i;
90   for(Short_t i=0; i<kFifoWords; i++) iFifoWord[i]=805306384+i;
91 }
92
93 UInt_t AliITSRawStreamSDD::ReadBits()
94 {
95 // read bits from the given channel
96   UInt_t result = (fChannelData[fCarlosId][fChannel] & ((1<<fReadBits[fCarlosId][fChannel]) - 1));
97   fChannelData[fCarlosId][fChannel] >>= fReadBits[fCarlosId][fChannel]; 
98   fLastBit[fCarlosId][fChannel] -= fReadBits[fCarlosId][fChannel];
99   return result;
100 }
101
102 Int_t AliITSRawStreamSDD::DecompAmbra(Int_t value) const
103 {
104   // AMBRA decompression (from 8 to 10 bit)
105   
106   if ((value & 0x80) == 0) {
107     return value & 0x7f;
108   } else if ((value & 0x40) == 0) {
109     return 0x081 + ((value & 0x3f) << 1);
110   } else if ((value & 0x20) == 0) {
111     return 0x104 + ((value & 0x1f) << 3);
112   } else {
113     return 0x208 + ((value & 0x1f) << 4);
114   }
115   
116 }
117
118 Bool_t AliITSRawStreamSDD::Next()
119 {
120 // read the next raw digit
121 // returns kFALSE if there is no digit left
122
123   fPrevModuleID = fModuleID;
124   fDDL=fRawReader->GetDDLID();
125   Int_t ddln = fRawReader->GetDDLID();
126   if(ddln <0) ddln=0;
127   Bool_t kSkip = ResetSkip(ddln);
128   if(!kSkip) return kSkip;
129
130   for(Int_t i=0;i<kModulesPerDDL;i++) iCountFoot[i]=0; 
131
132   while (kTRUE) {
133     if ((fChannel < 0) || (fLastBit[fCarlosId][fChannel] < fReadBits[fCarlosId][fChannel])) {
134
135       if (!fRawReader->ReadNextInt(fData)){ 
136         return kFALSE;  // read next word 
137       }
138       ddln = fRawReader->GetDDLID();
139
140       if(ddln!=fDDL) Reset();
141       if(ddln < 0 || ddln > (kDDLsNumber-1)) ddln  = 0;
142       kSkip = ResetSkip(ddln);
143       if(!kSkip) return kFALSE;  // read next word
144
145       fChannel = -1;
146       if((fData>=iCarlosWord[0]&&fData<=iCarlosWord[11])||(fData>=iFifoWord[0]&&fData<=iFifoWord[3])){
147         for(Short_t ik=0;ik<kCarlosWords;ik++) {
148           if(fData==iCarlosWord[ik]) {
149             fCarlosId = ik;
150             Int_t iFifoIdx = ik/3;
151             if(fNCarlos == 8) {
152               if(ik==2) iFifoIdx = 1;
153               if(ik==4 || ik==5)  iFifoIdx = 2;
154               if(ik==6 || ik==7)  iFifoIdx = 3 ;
155             }       
156             fNfifo[iFifoIdx] = fCarlosId;
157           }
158         }
159         for(Short_t ik=0;ik<kFifoWords;ik++) {
160           if (fData==iFifoWord[ik]) { 
161             fCarlosId = fNfifo[ik];         
162           }
163         }
164         /*
165         if(fData==iCarlosWord[0]){
166           fCarlosId = 0;
167           fNfifo[0] = fCarlosId;
168         }
169         else if (fData==iCarlosWord[1]){
170           fCarlosId = 1;
171           fNfifo[0] = fCarlosId;
172         }
173         else if (fData==iCarlosWord[2]){
174           fCarlosId = 2;
175           if(fNCarlos == 8) fNfifo[1] = fCarlosId;
176           else fNfifo[0] = fCarlosId;
177         }
178         else if (fData==iCarlosWord[3]){
179           fCarlosId = 3;
180           fNfifo[1] = fCarlosId;
181         }
182         else if (fData==iCarlosWord[4]){
183           fCarlosId = 4;
184           if(fNCarlos == 8) fNfifo[2] = fCarlosId;
185           else fNfifo[1] = fCarlosId;
186         }
187         else if (fData==iCarlosWord[5]){
188           fCarlosId = 5;
189           if(fNCarlos == 8) fNfifo[2] = fCarlosId;
190           else fNfifo[1] = fCarlosId;
191         }
192         else if (fData==iCarlosWord[6]){
193           fCarlosId = 6;
194           if(fNCarlos == 8) fNfifo[3] = fCarlosId;
195           else fNfifo[2] = fCarlosId;
196         }
197         else if (fData==iCarlosWord[7]){
198           fCarlosId = 7;
199           if(fNCarlos == 8) fNfifo[3] = fCarlosId;
200           else fNfifo[2] = fCarlosId;
201         }
202         else if (fData==iCarlosWord[8]){
203           fCarlosId = 8;
204           fNfifo[2] = fCarlosId;
205         }
206         else if (fData==iCarlosWord[9]){
207           fCarlosId = 9;
208           fNfifo[3] = fCarlosId;
209         }
210         else if (fData==iCarlosWord[10]){
211           fCarlosId = 10;
212           fNfifo[3] = fCarlosId;
213         }
214         else if (fData==iCarlosWord[11]){
215           fCarlosId = 11;
216           fNfifo[3] = fCarlosId;
217         }
218         */
219
220       }
221       if(fData==1059004191) continue;
222       if (fNCarlos == 8 && (fCarlosId == 8 || fCarlosId == 9 || 
223                             fCarlosId ==10 || fCarlosId == 11))continue;
224       
225       fModuleID = fgkDDLModuleMap[fRawReader->GetDDLID()][fCarlosId];
226       
227       if ((fData >> 28) == 0x02) {           // header
228         fEventId = (fData >> 3) & 0x07FF;
229       } else if ((fData >> 28) == 0x03) {    // footer
230         iCountFoot[fCarlosId]++; // stop before the last word (last word=jitter)
231         Bool_t exit=kTRUE;
232         for(Int_t ic=0;ic<fNCarlos;ic++){
233           if(iCountFoot[ic]<3) exit=kFALSE;
234         }
235         if(exit){
236           return kFALSE;
237         }        
238       } else if ((fData >> 29) == 0x00) {    // error
239
240         if ((fData & 0x00000163) != 0) {         
241           fRawReader->AddMajorErrorLog(kDataError,Form("Error code = %8.8x",fData));     
242           AliWarning(Form("error codes = %8.8x",fData));
243           return kFALSE;          
244         }
245       } else if ((fData >> 30) == 0x01) {    // JTAG word
246         // ignored
247       } else if ((fData >> 30) == 0x02) {    // channel 0 data
248         
249         fChannel = 0;
250       } else if ((fData >> 30) == 0x03) {    // channel 1 data
251         fChannel = 1;
252       } else {                               // unknown data format
253         fRawReader->AddMajorErrorLog(kDataFormatErr,Form("Invalid data %8.8x",fData));
254         AliWarning(Form("invalid data: %8.8x\n", fData));
255         return kFALSE;
256       }
257       
258       if (fChannel >= 0) {   // add read word to the data
259         fChannelData[fCarlosId][fChannel] += 
260           (ULong64_t(fData & 0x3FFFFFFF) << fLastBit[fCarlosId][fChannel]);
261         fLastBit[fCarlosId][fChannel] += 30;
262       }
263
264     } else {  // decode data
265
266       if (fReadCode[fCarlosId][fChannel]) {// read the next code word
267         fChannelCode[fCarlosId][fChannel] = ReadBits();
268         fReadCode[fCarlosId][fChannel] = kFALSE;
269         fReadBits[fCarlosId][fChannel] = fgkCodeLength[fChannelCode[fCarlosId][fChannel]];
270
271       } else {                      // read the next data word
272         UInt_t data = ReadBits();
273         fReadCode[fCarlosId][fChannel] = kTRUE;
274         fReadBits[fCarlosId][fChannel] = 3;
275         if (fChannelCode[fCarlosId][fChannel] == 0) {         // set the time bin         
276           fTimeBin[fCarlosId][fChannel] = data;
277         } else if (fChannelCode[fCarlosId][fChannel] == 1) {  // next anode
278           fTimeBin[fCarlosId][fChannel] = 0;
279           fAnode[fCarlosId][fChannel]++;
280
281         }else {                                   // ADC signal data
282           fSignal = DecompAmbra(data + (1 << fChannelCode[fCarlosId][fChannel]) + fLowThreshold[fChannel]);
283           fCoord1 = fAnode[fCarlosId][fChannel];
284           fCoord2 = fTimeBin[fCarlosId][fChannel];
285           fTimeBin[fCarlosId][fChannel]++;
286           return kTRUE;
287         }
288       }
289     }
290   }
291   return kFALSE;  
292
293 }
294
295 void AliITSRawStreamSDD::Reset(){
296
297   //reset data member for a new ddl
298   for(Int_t i=0;i<2;i++){
299     for(Int_t ic=0;ic<kModulesPerDDL;ic++){
300       fChannelData[ic][i]=0;
301       fLastBit[ic][i]=0;
302       fChannelCode[ic][i]=0;
303       fReadCode[ic][i]=kFALSE;
304       fReadBits[ic][i]=0;
305       fTimeBin[ic][i]=0;
306       fAnode[ic][i]=0;     
307     }
308     fLowThreshold[i]=0;
309   }
310   for(Int_t i=0;i<kModulesPerDDL;i++) iCountFoot[i]=0;
311   
312   fChannel=-1;
313   fDDL=fRawReader->GetDDLID();
314 }
315
316 Bool_t AliITSRawStreamSDD::ResetSkip(Int_t ddln){
317   while (fSkip[ddln] < 9) {
318     if (!fRawReader->ReadNextInt(fData))return kFALSE;
319     if ((fData >> 30) == 0x01) continue;  // JTAG word
320     fSkip[ddln]++;
321   }
322   return kTRUE;
323 }
324