]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUSuze02.cxx
New vertexer with IB tracklets
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUSuze02.cxx
1 #include "TMath.h"
2 #include "AliITSUSuze02.h"
3
4 //*******************************************************************
5 //
6 //  Simulation of the SUZE02 readout
7 //  Origin: Serhiy.Senuykov@cern.ch
8 //  
9 //*******************************************************************
10
11 ClassImp(AliITSUSuze02)
12
13 AliITSUSuze02::AliITSUSuze02(Int_t Nrows, Int_t Ncols):
14   fNRowsChip(Nrows),
15   fNColsChip(Ncols),
16   fChip(new TMatrixF(Nrows,Ncols)), 
17   fTestColumnSize(2),
18   fTestRowSize(2),
19   fNWindowsPer32colsMax(0),
20   fNWindowsPerHalfFSBBMax(0),
21   fNWindowsPerFSBBMax(0),
22   fNDigitsEncoded(0), 
23   fNEncodedWindows(0),
24   fNDigitsLost(0),
25   fNLostWindows(0),
26   fDataSizePerChip(0),
27   fNWindowsPer32colsMin(0),
28   fNWindowsPerHalfFSBBMin(0),
29   fNWindowsPerFSBBMin(0)
30 {
31   if (Ncols%(kNumberOfFSBB*kNumberOfHalfFSBB) != 0) {
32     printf("Number of columns should be multiple of %d. SUZE matrix wasn't created\n",kNumberOfFSBB*kNumberOfHalfFSBB);
33   }
34 }  
35   
36 AliITSUSuze02::AliITSUSuze02(const AliITSUSuze02& suze): 
37   fNRowsChip(suze.fNRowsChip),
38   fNColsChip(suze.fNColsChip),
39   fChip(new TMatrixF(*suze.fChip)), 
40   fTestColumnSize(suze.fTestColumnSize),
41   fTestRowSize(suze.fTestRowSize),
42   fNWindowsPer32colsMax(suze.fNWindowsPer32colsMax),
43   fNWindowsPerHalfFSBBMax(suze.fNWindowsPerHalfFSBBMax),
44   fNWindowsPerFSBBMax(suze.fNWindowsPerFSBBMax),
45   fNDigitsEncoded(suze.fNDigitsEncoded), 
46   fNEncodedWindows(suze.fNEncodedWindows),
47   fNDigitsLost(suze.fNDigitsLost),
48   fNLostWindows(suze.fNLostWindows),
49   fDataSizePerChip(suze.fDataSizePerChip),
50   fNWindowsPer32colsMin(suze.fNWindowsPer32colsMin),
51   fNWindowsPerHalfFSBBMin(suze.fNWindowsPerHalfFSBBMin),
52   fNWindowsPerFSBBMin(suze.fNWindowsPerFSBBMin)
53 {
54 }
55
56 AliITSUSuze02 &AliITSUSuze02::operator=(const AliITSUSuze02& suze) {
57   if (&suze == this) return *this;
58
59   fNRowsChip = suze.fNRowsChip;
60   fNColsChip = suze.fNColsChip;
61   fChip = new TMatrixF(*suze.fChip);  
62   fTestColumnSize = suze.fTestColumnSize;
63   fTestRowSize = suze.fTestRowSize;
64   fNWindowsPer32colsMax = suze.fNWindowsPer32colsMax;
65   fNWindowsPerHalfFSBBMax = suze.fNWindowsPerHalfFSBBMax;
66   fNWindowsPerFSBBMax = suze.fNWindowsPerFSBBMax;
67   fNDigitsEncoded = suze.fNDigitsEncoded;
68   fNEncodedWindows = suze.fNEncodedWindows;
69   fNDigitsLost = suze.fNDigitsLost;
70   fNLostWindows = suze.fNLostWindows;
71   fDataSizePerChip = suze.fDataSizePerChip;
72   fNWindowsPer32colsMin = suze.fNWindowsPer32colsMin;
73   fNWindowsPerHalfFSBBMin = suze.fNWindowsPerHalfFSBBMin;
74   fNWindowsPerFSBBMin = suze.fNWindowsPerFSBBMin;
75
76   return *this;
77 }
78
79 AliITSUSuze02::~AliITSUSuze02() {
80   if(fChip) delete fChip;
81 }
82
83 void AliITSUSuze02::SetEncodingWindowSize(Int_t Wrows, Int_t Wcols){
84   fTestColumnSize=Wrows;
85   fTestRowSize=Wcols;
86 }
87
88 void AliITSUSuze02::SetQuotas(Int_t q32, Int_t qHalfFSBB, Int_t qFSBB){
89   fNWindowsPer32colsMax=q32;
90   fNWindowsPerHalfFSBBMax=qHalfFSBB;
91   fNWindowsPerFSBBMax=qFSBB;
92 }
93
94 void AliITSUSuze02::AddDigit(Int_t row, Int_t col){
95   (*fChip)(row,col)++;
96 }
97
98 //void AliITSUSuze02::Process(Bool_t Verbose){  
99 void AliITSUSuze02::Process(TH1F* OverflowCodes, TH1F* NDigitsPerEncodingWindowDist, Bool_t Verbose) {
100
101   //cout<<"Processing"<<endl;
102   //fChip->Print(); 
103   
104   Int_t NRowsFSBB=fNRowsChip;
105   Int_t NColsFSBB=fNColsChip/kNumberOfFSBB;
106   TMatrixF FSBB(NRowsFSBB,NColsFSBB);
107   
108   Int_t NRowsSuperLine=fTestColumnSize;
109   Int_t NColsSuperLine=NColsFSBB;   
110   TMatrixF SuperLineUp(NRowsSuperLine,NColsSuperLine);
111   TMatrixF SuperLineDown(NRowsSuperLine,NColsSuperLine);
112   
113   Int_t NRowsSuperLineX2=NRowsSuperLine*2; 
114   Int_t NColsSuperLineX2=NColsSuperLine; //SuperLineX2 and SuperLine size in columns is equal to FSBB.
115   TMatrixF SuperLineX2(NRowsSuperLineX2,NColsSuperLineX2); 
116   
117   TMatrixF TestRow(1,fTestRowSize);
118   TMatrixF TestColumn(fTestColumnSize,1);      
119   Int_t TestRowSum=0;
120   TMatrixF EncodingWindow(fTestColumnSize,fTestRowSize);
121   
122   Int_t EncodingWindowStartRow=0;
123   Int_t EncodingWindowStopRow=0;
124   Int_t EncodingWindowStartCol=0;
125   Int_t EncodingWindowStopCol=0;
126   
127   Int_t nMasks=fTestRowSize-1;
128   Int_t MaskSize=NRowsSuperLineX2-1;
129   TMatrixF Masks(MaskSize,nMasks);
130
131    //parameters for internal data size calculation
132   Int_t DataSizePerWindowInRow=8+fTestColumnSize*fTestRowSize+TMath::Ceil(TMath::Log2(fTestColumnSize));
133   Int_t DataSizePerSuperLineX2Header=30;   
134   
135   Int_t DataSizePerSuperLineX2=0;
136   
137   Bool_t Overflow32=kFALSE;
138   Bool_t OverflowHalfFSBB=kFALSE;
139   Bool_t OverflowFSBB=kFALSE;
140   
141   Int_t nWindowsPer32cols=fNWindowsPer32colsMax;
142   Int_t nWindowsPerHalfFSBB=fNWindowsPerHalfFSBBMax;
143   Int_t nWindowsPerFSBB=fNWindowsPerFSBBMax; 
144   
145   fNWindowsPer32colsMin=fNWindowsPer32colsMax;
146   fNWindowsPerHalfFSBBMin=fNWindowsPerHalfFSBBMax;
147   fNWindowsPerFSBBMin=fNWindowsPerFSBBMax;  
148   
149   fNEncodedWindows=0;
150   fNDigitsEncoded=0;
151   fNLostWindows=0;
152   fNDigitsLost=0;
153   fDataSizePerChip=0;    
154   
155   for(Int_t FSBBindex=0; FSBBindex<kNumberOfFSBB; FSBBindex++){
156     FSBB=fChip->GetSub(0,NRowsFSBB-1,FSBBindex*NColsFSBB,(FSBBindex+1)*NColsFSBB-1);
157     SuperLineDown=FSBB.GetSub(0,NRowsSuperLine-1,0,NColsSuperLine-1);
158     for(Int_t SuperLineX2StartRow=0; SuperLineX2StartRow<NRowsFSBB; SuperLineX2StartRow+=NRowsSuperLine){
159       if(nWindowsPerFSBB<fNWindowsPerFSBBMin) {fNWindowsPerFSBBMin=nWindowsPerFSBB;} //saving the lowest number of remaining windows
160       nWindowsPerFSBB=fNWindowsPerFSBBMax; //Reset number of available encoding windows for the new double SuperLine
161       OverflowFSBB=kFALSE;
162       SuperLineUp=SuperLineDown;
163       if(SuperLineX2StartRow+NRowsSuperLineX2<=NRowsFSBB){
164         SuperLineDown=FSBB.GetSub(SuperLineX2StartRow+NRowsSuperLine,SuperLineX2StartRow+NRowsSuperLineX2-1,0,NColsSuperLine-1);
165       }
166       else if(SuperLineX2StartRow+NRowsSuperLine<NRowsFSBB){
167         SuperLineDown.Zero();
168         SuperLineDown.SetSub(0,0,FSBB.GetSub(SuperLineX2StartRow+NRowsSuperLine,NRowsFSBB-1,0,NColsSuperLine-1));
169       } 
170       else{
171         SuperLineDown.Zero();
172       }
173       if(SuperLineUp.Sum()>0){
174         DataSizePerSuperLineX2=0; 
175         SuperLineX2.SetSub(0,0,SuperLineUp);
176         SuperLineX2.SetSub(NRowsSuperLine,0,SuperLineDown);
177         for(Int_t HalfFSBBindex=0; HalfFSBBindex<kNumberOfHalfFSBB; HalfFSBBindex++){
178           if(nWindowsPerHalfFSBB<fNWindowsPerHalfFSBBMin) {fNWindowsPerHalfFSBBMin=nWindowsPerHalfFSBB;}  
179           nWindowsPerHalfFSBB=fNWindowsPerHalfFSBBMax; //reset counter per HalfFSBB
180           OverflowHalfFSBB=kFALSE;
181           for(Int_t i=0; i<MaskSize; i++){ //reset masks to 1111111
182             for(Int_t j=0; j<nMasks; j++){
183               Masks(i,j)=1;
184             }
185           }
186           
187           for(Int_t TestRowStartCol=HalfFSBBindex*NColsSuperLineX2/kNumberOfHalfFSBB; TestRowStartCol<(HalfFSBBindex+1)*NColsSuperLineX2/kNumberOfHalfFSBB; TestRowStartCol++){
188             if(TestRowStartCol%32==0){
189                     if(nWindowsPer32cols<fNWindowsPer32colsMin) fNWindowsPer32colsMin=nWindowsPer32cols;
190               nWindowsPer32cols=fNWindowsPer32colsMax; //reset nWindowsPer32cols counter every 32 columns
191               Overflow32=kFALSE;
192             }
193             //apply masks
194             for(Int_t RowIndex=0; RowIndex<MaskSize; RowIndex++){
195               for(Int_t MaskIndex=0; MaskIndex<nMasks; MaskIndex++){
196                 if(Masks(RowIndex,MaskIndex)==0){
197                   // cout<<"Mask has zero bit at pos:"<<RowIndex<<":"<<MaskIndex<<endl;
198                   // cout<<"will clean the pixels at row:"<<RowIndex<<" from pos.:"<<TestRowStartCol<<" till "<< TestRowStartCol+(nMasks-MaskIndex)<<endl;
199                   for(Int_t ColIndex=TestRowStartCol; ColIndex<TestRowStartCol+(nMasks-MaskIndex); ColIndex++){
200                     if(ColIndex<(HalfFSBBindex+1)*NColsSuperLineX2/kNumberOfHalfFSBB){
201                       SuperLineX2(RowIndex,ColIndex)=0;
202                       // cout<<"Mask has zero bit. Cleaning at pos:"<<RowIndex<<":"<<ColIndex<<endl;
203                       if(RowIndex>=fTestColumnSize){
204                         SuperLineDown(RowIndex-fTestColumnSize,ColIndex)=0;
205                       }
206                     }
207                   }
208                   break;
209                 }
210               }
211             }
212             //shift masks
213             for(Int_t RowIndex=0; RowIndex<MaskSize; RowIndex++){
214               for(Int_t MaskIndex=nMasks-1; MaskIndex>0; MaskIndex--){
215                 Masks(RowIndex,MaskIndex)=Masks(RowIndex,MaskIndex-1);
216               }
217               Masks(RowIndex,0)=1;
218             }  
219             
220             for(Int_t TestRowStartRow=0; TestRowStartRow<fTestColumnSize; TestRowStartRow++){
221               TestRowSum=0;
222               for(Int_t TestRowIndex=0; TestRowIndex<fTestRowSize; TestRowIndex++){
223                 if(TestRowStartCol+TestRowIndex<(HalfFSBBindex+1)*NColsSuperLineX2/kNumberOfHalfFSBB){
224                   TestRowSum+=SuperLineX2(TestRowStartRow,TestRowStartCol+TestRowIndex);
225                 }
226                 if(TestRowSum>0){
227                   //cout<<"TestR at col n."<<TestRowStartCol<<" and row n."<<TestRowStartRow<<" has a hit"<<endl;
228                   break;
229                 }
230               }
231               if(TestRowSum>0){
232                 TestColumn=SuperLineX2.GetSub(TestRowStartRow,TestRowStartRow+fTestColumnSize-1,TestRowStartCol,TestRowStartCol);   
233                 if(TestColumn.Sum()>0){
234                   EncodingWindowStartRow=TestRowStartRow;
235                   EncodingWindowStopRow=EncodingWindowStartRow+fTestColumnSize-1;
236                   EncodingWindowStartCol=TestRowStartCol; 
237                   
238                   if(TestRowStartCol+fTestRowSize>(HalfFSBBindex+1)*NColsSuperLineX2/kNumberOfHalfFSBB){
239                     EncodingWindowStopCol=((HalfFSBBindex+1)*NColsSuperLineX2/kNumberOfHalfFSBB)-1;
240                   }
241                   else{
242                     EncodingWindowStopCol=EncodingWindowStartCol+fTestRowSize-1;
243                   }
244
245                   EncodingWindow.Zero();
246                   EncodingWindow.SetSub(0,0,SuperLineX2.GetSub(EncodingWindowStartRow,EncodingWindowStopRow,EncodingWindowStartCol,EncodingWindowStopCol));
247                   
248                   if(nWindowsPer32cols && nWindowsPerHalfFSBB && nWindowsPerFSBB){
249                     //cout<<"Will encode window starting at "<<TestRowStartRow<<":"<<TestRowStartCol<<endl;
250                     fNDigitsEncoded+=EncodingWindow.Sum();
251                     fNEncodedWindows++;
252                     OverflowCodes->Fill(0);
253                                 NDigitsPerEncodingWindowDist->Fill(EncodingWindow.Sum());
254                     nWindowsPerFSBB--;
255                     nWindowsPerHalfFSBB--;
256                     nWindowsPer32cols--;
257                     DataSizePerSuperLineX2+=DataSizePerWindowInRow;
258                     //cout<<"Windows left:"<<nWindowsPerFSBB<<":"<<nWindowsPerHalfFSBB<<":"<<nWindowsPer32cols<<endl;
259                   }
260                   else{
261                     fNDigitsLost+=EncodingWindow.Sum();
262                     fNLostWindows++;
263                     //cout<<"------ No encoding at col.:"<<TestRowStartCol<<" at SuperLineX2 that starts at row:"<<SuperLineX2StartRow<<endl;
264                     if(!nWindowsPer32cols) Overflow32=kTRUE;
265                     if(!nWindowsPerHalfFSBB) OverflowHalfFSBB=kTRUE;
266                     if(!nWindowsPerFSBB) OverflowFSBB=kTRUE;
267                     OverflowCodes->Fill(Overflow32+2*OverflowHalfFSBB+4*OverflowFSBB);
268                   } 
269                   for(Int_t k=TestRowStartRow;k<TestRowStartRow+fTestColumnSize; k++){
270                     Masks(k,0)=0;
271                     //cout<<"Setting Mask to 0 at "<<k<<endl;
272                     //Setting to 0 the first column of the encoding window. This part if probably present in the real SUZE.
273                     SuperLineX2(k,TestRowStartCol)=0;
274                     if(k>=fTestColumnSize){
275                       SuperLineDown(k-fTestColumnSize,TestRowStartCol)=0;
276                     }
277                   }
278                 }  
279                 break;
280               }
281             }
282           }
283         }
284         fDataSizePerChip+=(DataSizePerSuperLineX2+DataSizePerSuperLineX2Header);
285       }
286     }
287   }
288   if(Verbose){
289     cout<<fNDigitsEncoded<<" digits encoded in "<<fNEncodedWindows<<" windows"<<endl;
290     cout<<fNDigitsLost<<" digits lost in "<<fNLostWindows<<" windows"<<endl;
291   }
292 }
293
294 void AliITSUSuze02::GetResults(){
295   
296 }
297 /*
298 void AliITSUSuze02::InitHistos(){
299   fOverflowCodes = new TH1F("OverflowCodes","Overflow codes",8,0,8);
300   if(fNRowsChip*fNColsChip){
301     fNDigitsPerEncodingWindowDist = new TH1F("nDigitsPerEncodingWindowPerChip","nDigitsPerEncodingWindowPerChip",fTestColumnSize*fTestRowSize,1,fTestColumnSize*fTestRowSize+1);
302   }
303   else{
304     printf("Run AliITSUSuze02::SetEncodingWindowSize first\n");
305   }
306 }
307 */
308 void AliITSUSuze02::ResetChip(){
309   fChip->Zero();
310 }