ATLAS Offline Software
Loading...
Searching...
No Matches
AFP_NoisyPixelTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7int AFP_NoisyPixelTool::Identify(std::shared_ptr<const TH2F> input, std::vector<TH2F>& output) const
8{
9 if(output.size()!=1)
10 {
11 return 0;
12 }
13
14 output[0].Reset();
15 output.resize(4, output[0]);//container now has 4 copies of the 0'th element
16 auto name = [&output](const char * n)->TString{
17 return Form(n,output[0].GetName());
18 };
19 auto title = [&output](const char * t)->TString{
20 return Form(t,output[0].GetTitle());
21 };
22 output[1].SetNameTitle(name("leffpixels_found_%s"), title("low efficiency pixels, found, %s"));
23 output[2].SetNameTitle(name("noisypixels_eff_%s"), title("noisy pixels, efficiency, %s"));
24 output[3].SetNameTitle(name("leffpixels_eff_%s"), title("low efficiency pixels, efficiency, %s"));
25 output[0].SetNameTitle(name("noisypixels_found_%s"), title("noisy pixels, found, %s"));
26
27 TH2F& noisypixels_found_output = output[0];
28 TH2F& leffpixels_found_output = output[1];
29 TH2F& noisypixels_eff_output = output[2];
30 TH2F& leffpixels_eff_output = output[3];
31
32
33 if(input->GetMaximum()<0.5) return 0;
34
35 std::vector<TH2F> vec_found_leff{},vec_found_noisy{},vec_eff_leff{},vec_eff_noisy{};
36
37 for(const auto& method : m_methods)
38 {
39 TH2F tmp_eff_noisy, tmp_eff_leff;
40 TH2F tmp_found_leff, tmp_found_noisy;
41
42 if(method!="FIT")
43 {
44 // count inactive pixels around each pixel
45 TH2F tmp_inact_pix_around=countInactivePixelsAround(input);
46
47 // find low eff. and noisy pixels
48 std::tie(tmp_found_leff,tmp_found_noisy,tmp_eff_leff, tmp_eff_noisy) = findLEffAndNoisyPixels(input, tmp_inact_pix_around, method);
49
50 // reconsider some pixels
51 filterLEffPixelsAroundNoisy(input, tmp_inact_pix_around, tmp_found_leff,tmp_found_noisy,tmp_eff_leff,tmp_eff_noisy, method);
52 }
53 else
54 {
55 // method is "FIT"
56
57 tmp_eff_noisy=TH2F(*input);
58 tmp_eff_noisy.SetName("tmp_eff_noisy");
59 tmp_eff_noisy.Reset();
60 tmp_eff_leff=TH2F(*input);
61 tmp_eff_leff.SetName("tmp_eff_leff");
62 tmp_eff_leff.Reset();
63 tmp_found_noisy=TH2F(*input);
64 tmp_found_noisy.SetName("tmp_found_noisy");
65 tmp_found_noisy.Reset();
66 tmp_found_leff=TH2F(*input);
67 tmp_found_leff.SetName("tmp_found_leff");
68 tmp_found_leff.Reset();
69
70 auto p=makeFits(input);
71
72 const double Threshold_LEFF=5.;
73 const double Threshold_NOISY=(m_threshNoisy-1)*100;
74
75 for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
76 {
77 for(int row_ID=1; row_ID<=input->GetNbinsX(); row_ID++)
78 {
79 double ratio_leff=0, ratio_noisy=0;
80 double gauss_helper_12 = (col_ID-p[1][row_ID])/p[2][row_ID];
81 double gauss_helper_45 = (col_ID-p[4][row_ID])/p[5][row_ID];
82 double fit = p[0][row_ID]*exp(-0.5*gauss_helper_12*gauss_helper_12)
83 + p[3][row_ID]*exp(-0.5*gauss_helper_45*gauss_helper_45) + p[6][row_ID];
84 double sigma=std::sqrt(fit);
85 if(sigma<1) sigma=1;
86
87 double bin_content=input->GetBinContent(row_ID,col_ID);
88
89 if(fit>=bin_content) ratio_leff = std::abs(fit-bin_content)/sigma;
90 else ratio_noisy = std::abs(fit-bin_content)/sigma;
91
92 if(bin_content!=0)
93 {
94 if( ratio_leff > Threshold_LEFF) tmp_found_leff.SetBinContent(row_ID,col_ID,1);
95 if( ratio_noisy > Threshold_NOISY) tmp_found_noisy.SetBinContent(row_ID,col_ID,1);
96 }
97 tmp_eff_leff.SetBinContent(row_ID,col_ID,100.*ratio_leff);
98 tmp_eff_noisy.SetBinContent(row_ID,col_ID,100.*ratio_noisy);
99 }
100 }
101 }
102
103 // save output for this method
104 vec_found_leff.push_back(tmp_found_leff);
105 vec_found_noisy.push_back(tmp_found_noisy);
106 vec_eff_leff.push_back(tmp_eff_leff);
107 vec_eff_noisy.push_back(tmp_eff_noisy);
108 }
109
110 // sum found pixels for all methods, set efficiency to the maximum of all methods
111 for(unsigned int m=0; m<m_methods.size(); m++)
112 {
113 for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
114 {
115 for(int row_ID=1; row_ID<=input->GetNbinsX(); row_ID++)
116 {
117 noisypixels_found_output.Fill(row_ID,col_ID, vec_found_noisy[m].GetBinContent(row_ID,col_ID));
118 leffpixels_found_output.Fill(row_ID,col_ID, vec_found_leff[m].GetBinContent(row_ID,col_ID));
119
120 if(m==0)
121 {
122 noisypixels_eff_output.SetBinContent(row_ID,col_ID, vec_eff_noisy[m].GetBinContent(row_ID,col_ID));
123 leffpixels_eff_output.SetBinContent(row_ID,col_ID, vec_eff_leff[m].GetBinContent(row_ID,col_ID));
124 }
125 else
126 {
127 if(noisypixels_eff_output.GetBinContent(row_ID,col_ID)<vec_eff_noisy[m].GetBinContent(row_ID,col_ID))
128 {
129 noisypixels_eff_output.SetBinContent(row_ID,col_ID, vec_eff_noisy[m].GetBinContent(row_ID,col_ID));
130 }
131 if(leffpixels_eff_output.GetBinContent(row_ID,col_ID)>vec_eff_leff[m].GetBinContent(row_ID,col_ID))
132 {
133 leffpixels_eff_output.SetBinContent(row_ID,col_ID, vec_eff_leff[m].GetBinContent(row_ID,col_ID));
134 }
135 }
136 }
137 }
138 }
139
140 // report only pixels that are identified by all methods
141 for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
142 {
143 for(int row_ID=1; row_ID<=input->GetNbinsX(); row_ID++)
144 {
145 if(noisypixels_found_output.GetBinContent(row_ID,col_ID)==m_methods.size()) noisypixels_found_output.SetBinContent(row_ID,col_ID,1);
146 else noisypixels_found_output.SetBinContent(row_ID,col_ID,0);
147
148 if(leffpixels_found_output.GetBinContent(row_ID,col_ID)==m_methods.size()) leffpixels_found_output.SetBinContent(row_ID,col_ID,1);
149 else leffpixels_found_output.SetBinContent(row_ID,col_ID,0);
150 }
151 }
152
153 return 0;
154}
155
156
157std::vector<std::pair<int,int>> AFP_NoisyPixelTool::getLegitPixels(std::shared_ptr<const TH2F> input, const int col_ID, const int row_ID, const std::string& method="8_PIX") const
158{
159 // method="8_PIX" means to investigate all pixels around
160
161 std::vector<std::pair<int,int>> legit_pixels={};
162
163 for(int r=-1; r<2; r++)
164 {
165 for(int c=-1; c<2; c++)
166 {
167 if( (row_ID+r)>=1 && (col_ID+c)>=1 && (col_ID+c)<=input->GetNbinsY() && (row_ID+r)<=input->GetNbinsX() && (r!=0 || c!=0))
168 {
169 if( (method=="2_ROW" && c==0) || (method=="2_COL" && r==0)
170 || (method=="4_PIX" && (c==0 || r==0)) || (method=="8_PIX") )
171 {
172 legit_pixels.push_back(std::pair<int,int>(row_ID+r,col_ID+c));
173 }
174 }
175 }
176 }
177
178 return legit_pixels;
179}
180
181
182TH2F AFP_NoisyPixelTool::countInactivePixelsAround(std::shared_ptr<const TH2F> input) const
183{
184 TH2F tmp_inact_pix_around(*input);
185 tmp_inact_pix_around.SetName("tmp_inact_pix_around");
186 tmp_inact_pix_around.Reset();
187
188 for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
189 {
190 for(int row_ID=1; row_ID<=input->GetNbinsX(); row_ID++)
191 {
192 int inactive_pixels_around=0;
193 auto legit_pixels=getLegitPixels(input, col_ID, row_ID);
194
195 for(auto legpix : legit_pixels)
196 {
197 if(input->GetBinContent(legpix.first, legpix.second)<1) ++inactive_pixels_around;
198 }
199
200 tmp_inact_pix_around.SetBinContent(row_ID,col_ID, inactive_pixels_around+0.01);
201 }
202 }
203
204 return tmp_inact_pix_around;
205}
206
207
208double AFP_NoisyPixelTool::getNeighbours(std::shared_ptr<const TH2F> input, int row_ID, int col_ID) const
209{
210 if( row_ID!=input->GetNbinsX() && row_ID!=1 && col_ID!=1 && col_ID!=input->GetNbinsY()) return 8.;
211 else if( (row_ID==input->GetNbinsX() || row_ID==1) && (col_ID==1 || col_ID==input->GetNbinsY())) return 3.;
212 else return 5;
213}
214
215
216std::tuple<TH2F,TH2F,TH2F,TH2F> AFP_NoisyPixelTool::findLEffAndNoisyPixels(std::shared_ptr<const TH2F> input, const TH2F& tmp_inact_pix_around, const std::string& method) const
217{
218 TH2F tmp_found_leff(*input);
219 tmp_found_leff.SetName("tmp_found_leff");
220 tmp_found_leff.Reset();
221
222 TH2F tmp_found_noisy(*input);
223 tmp_found_noisy.SetName("tmp_found_noisy");
224 tmp_found_noisy.Reset();
225
226 TH2F tmp_eff_leff(*input);
227 tmp_eff_leff.SetName("tmp_eff_leff");
228 tmp_eff_leff.Reset();
229
230 TH2F tmp_eff_noisy(*input);
231 tmp_eff_noisy.SetName("tmp_eff_noisy");
232 tmp_eff_noisy.Reset();
233
234 for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
235 {
236 for(int row_ID=1; row_ID<=input->GetNbinsX(); row_ID++)
237 {
238 if(input->GetBinContent(row_ID, col_ID)!=0 && tmp_inact_pix_around.GetBinContent(row_ID, col_ID)<0.1)
239 {
240 double npixels=getNeighbours(input, row_ID, col_ID);
241
242 double sum=0.;
243 auto legit_pixels=getLegitPixels(input, col_ID, row_ID, method);
244 for(auto legpix : legit_pixels)
245 {
246 sum+=input->GetBinContent(legpix.first,legpix.second);
247 }
248
249 double av_noisy=0;
250 double av_leff=0;
251
252 // TODO: something about edges?
253 if(method=="2_ROW" || method=="2_COL")
254 {
255 av_noisy=sum/(std::abs(npixels-4)*0.5);
256 av_leff=sum/2.;
257 }
258 else if(method=="4_PIX")
259 {
260 av_noisy=sum/(std::ceil(npixels*0.5));
261 av_leff=sum/4.;
262 }
263 else if(method=="8_PIX")
264 {
265 av_noisy=sum/(npixels);
266 av_leff=sum/8.;
267 }
268
269
270 double ratio_leff = input->GetBinContent(row_ID,col_ID)/av_leff;
271 double ratio_noisy = input->GetBinContent(row_ID,col_ID)/av_noisy;
272 tmp_eff_leff.SetBinContent(row_ID,col_ID,100.*ratio_leff);
273 tmp_eff_noisy.SetBinContent(row_ID,col_ID,100.*ratio_noisy);
274
275 if(ratio_leff < m_threshLEff) tmp_found_leff.SetBinContent(row_ID,col_ID,1.);
276 if(ratio_noisy > m_threshNoisy) tmp_found_noisy.SetBinContent(row_ID,col_ID,1.);
277 }
278 }
279 }
280
281 return {tmp_found_leff,tmp_found_noisy,tmp_eff_leff,tmp_eff_noisy};
282}
283
284
285void AFP_NoisyPixelTool::filterLEffPixelsAroundNoisy(std::shared_ptr<const TH2F> input, const TH2F& tmp_inact_pix_around, TH2F& tmp_found_leff, TH2F& tmp_found_noisy, TH2F& tmp_eff_leff, TH2F& tmp_eff_noisy, const std::string& method) const
286{
287 for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
288 {
289 for(int row_ID=1; row_ID<=input->GetNbinsX(); row_ID++)
290 {
291 if(input->GetBinContent(row_ID,col_ID)==0) continue;
292
293 int inactive_pixels_around=tmp_inact_pix_around.GetBinContent(row_ID, col_ID);
294 int leff_pixels_around=0;
295 auto legit_pixels=getLegitPixels(input, col_ID,row_ID);
296 for(auto legpix : legit_pixels)
297 {
298 if(tmp_found_leff.GetBinContent(legpix.first,legpix.second)>0) leff_pixels_around++;
299 }
300
301 if(leff_pixels_around!=0 && inactive_pixels_around==0)
302 {
303 //if low eff. pixels around > 0, recalculate avarage
304 double re_av=0.;
305 int npix=0;
306
307 auto legit_pixels=getLegitPixels(input, col_ID,row_ID, method);
308 for(auto legpix : legit_pixels)
309 {
310 if(tmp_found_leff.GetBinContent(legpix.first,legpix.second)==0)
311 {
312 re_av+=input->GetBinContent(legpix.first,legpix.second);
313 npix++;
314 }
315 }
316
317 if(npix==0)
318 {
319 re_av=0;
320 for(auto legpix : legit_pixels)
321 {
322 re_av+=input->GetBinContent(legpix.first,legpix.second);
323 }
324
325 if(re_av/input->GetBinContent(row_ID,col_ID)<1.5 && re_av/input->GetBinContent(row_ID,col_ID)>0.5) re_av=1;
326 }
327 else
328 {
329 re_av/=npix;
330 }
331
332 double ratio_noisy = m_sensitivity+input->GetBinContent(row_ID,col_ID)/re_av;
333 double ratio_leff = input->GetBinContent(row_ID,col_ID)/re_av;
334 if(re_av==1)
335 {
336 tmp_eff_leff.SetBinContent(row_ID,col_ID,1);
337 tmp_eff_noisy.SetBinContent(row_ID,col_ID,1);
338 }
339
340 if(row_ID!=input->GetNbinsX() && row_ID!=1 && col_ID!=1 && col_ID!=input->GetNbinsY())
341 {
342 if( ((ratio_noisy>m_threshNoisy && re_av!=1) || re_av==1) )
343 {
344 tmp_found_noisy.SetBinContent(row_ID,col_ID,1);
345 if(re_av!=1) tmp_eff_noisy.SetBinContent(row_ID,col_ID,100.*ratio_noisy);
346 }
347 if( ((ratio_leff<m_threshLEff && re_av!=1) || re_av==1) )
348 {
349 tmp_found_leff.SetBinContent(row_ID,col_ID,1);
350 if(re_av!=1) tmp_eff_leff.SetBinContent(row_ID,col_ID,100.*ratio_leff);
351 }
352 }
353 }
354 }
355 }
356
357 return;
358}
359
360
361std::vector<std::vector<double>> AFP_NoisyPixelTool::makeFits(std::shared_ptr<const TH2F> input) const
362{
363 const int nParams=7;
364 std::vector<std::vector<double>> params(nParams, std::vector<double>(input->GetNbinsX()+1, 0.));
365
366 for(int row_ID=0; row_ID<=input->GetNbinsX(); row_ID++)
367 {
368 std::string srow_ID = std::to_string(row_ID);
369 std::unique_ptr<TH1F> hist(new TH1F(("ROW"+srow_ID).c_str(), ("ROW"+srow_ID).c_str(), input->GetNbinsY(), 0.5, input->GetNbinsX()+0.5));
370
371 for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
372 {
373 hist->Fill(col_ID, 1.0*input->GetBinContent(row_ID,col_ID));
374 }
375
376 //Set double gaussian fit with starting parameters
377 std::shared_ptr<TF1> FIT_2(new TF1("gaus","gaus(0)+gaus(3)+[6]",0,input->GetNbinsY(),nParams));
378 if(row_ID==1)
379 {
380 FIT_2->SetParameters(25000,hist->GetMean(),hist->GetStdDev()/std::sqrt(2),
381 5000,hist->GetMean(),hist->GetStdDev()/std::sqrt(2),
382 std::max(hist->GetBinContent(2),std::max(hist->GetBinContent(3),hist->GetBinContent(4))) );
383 }
384 else
385 {
386 for(int par=0;par<nParams;++par)
387 {
388 FIT_2->SetParameter(par,params.at(par).at(row_ID-1));
389 }
390 }
391
392 //fit
393 hist->Fit(FIT_2.get(),"Q");
394
395 for(int par=0;par<nParams;++par)
396 {
397 params.at(par).at(row_ID-1)=FIT_2->GetParameter(par);
398 }
399 }
400
401 return params;
402}
403
404
std::tuple< TH2F, TH2F, TH2F, TH2F > findLEffAndNoisyPixels(std::shared_ptr< const TH2F > input, const TH2F &inact, const std::string &method) const
std::vector< std::vector< double > > makeFits(std::shared_ptr< const TH2F > input) const
TH2F countInactivePixelsAround(std::shared_ptr< const TH2F > input) const
double getNeighbours(std::shared_ptr< const TH2F > input, int row_ID, int col_ID) const
int Identify(std::shared_ptr< const TH2F > input, std::vector< TH2F > &output) const override
void filterLEffPixelsAroundNoisy(std::shared_ptr< const TH2F > input, const TH2F &inact, TH2F &fleff, TH2F &fnoisy, TH2F &eleff, TH2F &enoisy, const std::string &method) const
std::vector< std::string > m_methods
std::vector< std::pair< int, int > > getLegitPixels(std::shared_ptr< const TH2F > input, const int col_ID, const int row_ID, const std::string &method) const
int r
Definition globals.cxx:22