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