ATLAS Offline Software
AFP_DeadPixelTool.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 
7 int AFP_DeadPixelTool::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  TH2F tmp_dead_pixels(*static_cast<TH2F*>(input->Clone("tmp_dead_pixels")));
15  tmp_dead_pixels.Reset();
16 
17  TH2F& deadpixels_output = output.at(0);
18  deadpixels_output.SetName(Form("deadpixels_%s",deadpixels_output.GetName()));
19  deadpixels_output.SetTitle(Form("dead pixels, %s",deadpixels_output.GetTitle()));
20  deadpixels_output.Reset();
21 
22  auto maxHitValue=input->GetMaximum();
23 
24  // first iteration (find singular dead pixels)
25  for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
26  {
27  for(int row_ID=1; row_ID<=input->GetNbinsX(); row_ID++)
28  {
29  if(input->GetBinContent(row_ID, col_ID)==0 )
30  {
31  double nNeighbours=getNeighbours(input, row_ID,col_ID);
32 
33  int sum_neighbours=0;
34  int inactive_pixels_around=0;
35 
36  auto legit_pixels=getLegitPixels(input, col_ID, row_ID);
37  for(auto legpix : legit_pixels)
38  {
39  sum_neighbours+=input->GetBinContent(legpix.first, legpix.second);
40  if(input->GetBinContent(legpix.first, legpix.second)<1) ++inactive_pixels_around;
41  }
42 
43  if(inactive_pixels_around<4 && round(sum_neighbours/nNeighbours)>=m_range*maxHitValue)
44  {
45  deadpixels_output.Fill(row_ID, col_ID);
46  tmp_dead_pixels.SetBinContent(row_ID,col_ID, 1);
47  }
48  }
49  }
50  }
51 
52  // second iteration (find groups of dead pixels)
53  bool newDEAD=true;
54  while(newDEAD)
55  {
56  newDEAD=false;
57  for(int col_ID=1; col_ID<=input->GetNbinsY(); col_ID++)
58  {
59  for(int row_ID=1; row_ID<=input->GetNbinsX(); row_ID++)
60  {
61  if(tmp_dead_pixels.GetBinContent(row_ID,col_ID)==0 && input->GetBinContent(row_ID, col_ID)==0)
62  {
63  double nNeighbours=getNeighbours(input, row_ID,col_ID);
64 
65  int sum_neighbours=0;
66  int inactive_pixels_around=0;
67  int dead_pixels_around=0;
68 
69  auto legit_pixels=getLegitPixels(input, col_ID, row_ID);
70  for(auto legpix : legit_pixels)
71  {
72  sum_neighbours+=input->GetBinContent(legpix.first, legpix.second);
73  if(input->GetBinContent(legpix.first, legpix.second)<1) ++inactive_pixels_around;
74 
75  if(tmp_dead_pixels.GetBinContent(legpix.first, legpix.second)>0) dead_pixels_around++;
76  }
77 
78 
79  if(dead_pixels_around>0 && inactive_pixels_around>3 && round(sum_neighbours/nNeighbours)>=m_range*maxHitValue)
80  {
81  deadpixels_output.Fill(row_ID, col_ID);
82  tmp_dead_pixels.SetBinContent(row_ID,col_ID, 1);
83  newDEAD=true;
84  }
85  }
86  }
87  }
88  }
89 
90  return 0;
91 }
92 
93 
94 std::vector<std::pair<int,int>> AFP_DeadPixelTool::getLegitPixels(std::shared_ptr<const TH2F> input, const int col_ID, const int row_ID) const
95 {
96  std::vector<std::pair<int,int>> legit_pixels={};
97 
98  for(int r=-1; r<2; r++)
99  {
100  for(int c=-1; c<2; c++)
101  {
102  if( (row_ID+r)>=1 && (col_ID+c)>=1 && (col_ID+c)<=input->GetNbinsY() && (row_ID+r)<=input->GetNbinsX() && (r!=0 || c!=0))
103  {
104  legit_pixels.push_back(std::pair<int,int>(row_ID+r,col_ID+c));
105  }
106  }
107  }
108 
109  return legit_pixels;
110 }
111 
112 
113 double AFP_DeadPixelTool::getNeighbours(std::shared_ptr<const TH2F> input, int row_ID, int col_ID) const
114 {
115  if( row_ID!=input->GetNbinsX() && row_ID!=1 && col_ID!=1 && col_ID!=input->GetNbinsY()) return 8.; // inner pixels
116  else if( (row_ID==input->GetNbinsX() || row_ID==1) && (col_ID==1 || col_ID==input->GetNbinsY())) return 3.; // corner pixels
117  else return 5; // pixels along edges
118 }
beamspotman.r
def r
Definition: beamspotman.py:676
TH2F
Definition: rootspy.cxx:420
MuonGM::round
float round(const float toRound, const unsigned int decimals)
Definition: Mdt.cxx:27
AFP_DeadPixelTool::m_range
float m_range
Definition: AFP_DeadPixelTool.h:31
PlotPulseshapeFromCool.input
input
Definition: PlotPulseshapeFromCool.py:106
AFP_DeadPixelTool::getLegitPixels
std::vector< std::pair< int, int > > getLegitPixels(std::shared_ptr< const TH2F > input, const int col_ID, const int row_ID) const
Definition: AFP_DeadPixelTool.cxx:94
TH2F::GetBinContent
double GetBinContent(int) const
Definition: rootspy.cxx:425
TH2F::SetBinContent
void SetBinContent(int, double)
Definition: rootspy.cxx:426
merge.output
output
Definition: merge.py:17
AFP_DeadPixelTool::getNeighbours
double getNeighbours(std::shared_ptr< const TH2F > input, int row_ID, int col_ID) const
Definition: AFP_DeadPixelTool.cxx:113
AFP_DeadPixelTool::Identify
int Identify(std::shared_ptr< const TH2F > input, std::vector< TH2F > &output) const override
Definition: AFP_DeadPixelTool.cxx:7
python.compressB64.c
def c
Definition: compressB64.py:93
AFP_DeadPixelTool.h
ActsTrk::nNeighbours
@ nNeighbours
Definition: StripInformationHelper.h:13