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