CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripHotStripAlgorithmFromClusterOccupancy.cc
Go to the documentation of this file.
4 
5 
6 
7 
9  prob_(1.E-7),
10  ratio_(1.5),
11  MinNumEntries_(0),
12  MinNumEntriesPerStrip_(0),
13  Nevents_(0),
14  occupancy_(0),
15  OutFileName_("Occupancy.root"),
16  tTopo(theTopo),
17  UseInputDB_(iConfig.getUntrackedParameter<bool>("UseInputDB",false))
18  {
20  }
21 
22 
24  LogTrace("SiStripHotStripAlgorithmFromClusterOccupancy")<<"[SiStripHotStripAlgorithmFromClusterOccupancy::~SiStripHotStripAlgorithmFromClusterOccupancy] "<<std::endl;
25 }
26 
28 
29  LogTrace("SiStripHotStripAlgorithmFromClusterOccupancy")<<"[SiStripHotStripAlgorithmFromClusterOccupancy::extractBadStrips] "<<std::endl;
30 
31 
32 
33  if (WriteOutputFile_==true){
34  f = new TFile(OutFileName_.c_str(),"RECREATE");
35  f->cd();
36 
37  striptree = new TTree("stripOccupancy","tree");
38 
39  striptree->Branch("DetRawId", &detrawid, "DetRawId/I");
40  striptree->Branch("SubDetId", &subdetid, "SubDetId/I");
41  striptree->Branch("Layer_Ring", &layer_ring, "Layer_Ring/I");
42  striptree->Branch("Disc", &disc, "Disc/I");
43  striptree->Branch("IsBack", &isback, "IsBack/I");
44  striptree->Branch("IsExternalString", &isexternalstring, "IsExternalString/I");
45  striptree->Branch("IsZMinusSide", &iszminusside, "IsZMinusSide/I");
46  striptree->Branch("RodStringPetal", &rodstringpetal, "RodStringPetal/I");
47  striptree->Branch("IsStereo", &isstereo, "IsStereo/I");
48  striptree->Branch("ModulePosition", &module_position, "ModulePosition/I");
49  striptree->Branch("NumberOfStrips", &number_strips, "NumberOfStrips/I");
50  striptree->Branch("StripNumber", &strip_number, "StripNumber/I");
51  striptree->Branch("APVChannel", &apv_channel, "APVChannel/I");
52  striptree->Branch("StripGlobalPositionX", &global_position_x, "StripGlobalPositionX/F");
53  striptree->Branch("StripGlobalPositionY", &global_position_y, "StripGlobalPositionY/F");
54  striptree->Branch("StripGlobalPositionZ", &global_position_z, "StripGlobalPositionZ/F");
55  striptree->Branch("IsHot", &isHot, "IsHot/I");
56  striptree->Branch("HotStripsPerAPV", &hotStripsPerAPV, "HotStripsPerAPV/I");
57  striptree->Branch("HotStripsPerModule", &hotStripsPerModule,"HotStripsPerModule/I");
58  striptree->Branch("StripOccupancy", &stripOccupancy, "StripOccupancy/D");
59  striptree->Branch("StripHits", &stripHits, "StripHits/I");
60  striptree->Branch("PoissonProb", &poissonProb, "PoissonProb/D");
61  striptree->Branch("MedianAPVHits", &medianAPVHits, "MedianAPVHits/D");
62  striptree->Branch("AvgAPVHits", &avgAPVHits, "AvgAPVHits/D");
63 }
64 
65 
66  HistoMap::iterator it=DM.begin();
67  HistoMap::iterator itEnd=DM.end();
68  std::vector<unsigned int> badStripList;
69  uint32_t detid;
70  for (;it!=itEnd;++it){
71  pHisto phisto;
72  detid=it->first;
73 
74  DetId detectorId=DetId(detid);
75  phisto._SubdetId=detectorId.subdetId();
76 
77  if (edm::isDebugEnabled())
78  LogTrace("SiStripHotStrip") << "Analyzing detid " << detid<< std::endl;
79 
80  int numberAPVs = (int)(it->second.get())->GetNbinsX()/128;
81 
82  // Set the values for the tree:
83 
84  detrawid = detid;
85  subdetid = detectorId.subdetId();
86  number_strips = (int)(it->second.get())->GetNbinsX();
87  switch (detectorId.subdetId())
88  {
91  disc = -1;
93  isback = -1;
95  else isexternalstring = 0;
97  else iszminusside = 0;
100  break;
101 
102  case StripSubdetector::TID :
106  if (tTopo->tidIsBackRing(detrawid)) isback = 1;
107  else isback = 0;
109  else iszminusside = 0;
110  isexternalstring = -1;
111  rodstringpetal = -1;
113  break;
114 
115  case StripSubdetector::TOB :
117  disc = -1;
119  isback = -1;
121  else iszminusside = 0;
122  isexternalstring = -1;
125  break;
126 
127  case StripSubdetector::TEC :
131  if (tTopo->tecIsBackPetal(detrawid)) isback = 1;
132  else isback = 0;
134  else iszminusside = 0;
135  isexternalstring = -1;
138  break;
139 
140  default :
141  std::cout << "### Detector does not belong to TIB, TID, TOB or TEC !? ###" << std::endl;
142  std::cout << "### DetRawId: " << detrawid << " ###" << std::endl;
143  }
144 
145  // End: Set the values for the tree.
146 
147 
148  pQuality=OutSiStripQuality;
149  badStripList.clear();
150 
151  for (int i=0; i<768; i++){
152  ishot[i] = 0;
153  stripoccupancy[i] = 0;
154  striphits[i] = 0;
155  poissonprob[i] = 0;
156  hotstripsperapv[i/128] = 0;
157  }
158 
159  hotstripspermodule = 0;
160 
161  for (int apv=0; apv<numberAPVs; apv++){
162  if(UseInputDB_){
163  if(InSiStripQuality->IsApvBad(detid,apv) ){
164  if(edm::isDebugEnabled())
165  LogTrace("SiStripHotStrip")<<"(Module and Apv number) "<<detid<<" , "<<apv<<" excluded by input ESetup."<<std::endl;
166  continue;//if the apv is already flagged as bad, continue.
167  }
168  else{
169  if(edm::isDebugEnabled())
170  LogTrace("SiStripHotStrip")<<"(Module and Apv number) "<<detid<<" , "<<apv<<" good by input ESetup."<<std::endl;
171  }
172  }
173 
174  phisto._th1f = new TH1F("tmp","tmp",128,0.5,128.5);
175  int NumberEntriesPerAPV=0;
176 
177  for (int strip=0; strip<128; strip++){
178  phisto._th1f->SetBinContent(strip+1,(it->second.get())->GetBinContent((apv*128)+strip+1));
179  NumberEntriesPerAPV += (int)(it->second.get())->GetBinContent((apv*128)+strip+1);
180  }
181 
182  phisto._th1f->SetEntries(NumberEntriesPerAPV);
183  phisto._NEntries=(int)phisto._th1f->GetEntries();
184  phisto._NEmptyBins=0;
185 
186  LogTrace("SiStripHotStrip") << "Number of clusters in APV " << apv << ": " << NumberEntriesPerAPV << std::endl;
187 
188  iterativeSearch(phisto,badStripList,apv);
189 
190  delete phisto._th1f;
191  }
192 
193  const StripGeomDetUnit* theStripDet = dynamic_cast<const StripGeomDetUnit*>( (TkGeom->idToDet(detectorId)) );
194  const StripTopology* theStripTopol = dynamic_cast<const StripTopology*>( &(theStripDet->specificTopology()) );
195 
196  for (int strip=0; strip<number_strips; strip++)
197  {
198  strip_number = strip+1;
199  apv_channel = (strip%128)+1;
200  isHot = ishot[strip];
202  stripHits = striphits[strip];
203  poissonProb = poissonprob[strip];
204  medianAPVHits = medianapvhits[strip/128];
205  avgAPVHits = avgapvhits[strip/128];
206 
208  hotStripsPerAPV = hotstripsperapv[strip/128];
209 
210  LocalPoint pos_strip_local = theStripTopol->localPosition(strip);
211  GlobalPoint pos_strip_global = (TkGeom->idToDet(detectorId))->surface().toGlobal(pos_strip_local);
212 
213  global_position_x = pos_strip_global.x();
214  global_position_y = pos_strip_global.y();
215  global_position_z = pos_strip_global.z();
216 
217  if (WriteOutputFile_==true) striptree->Fill();
218  }
219 
220  if (badStripList.begin()==badStripList.end())
221  continue;
222 
223  OutSiStripQuality->compact(detid,badStripList);
224 
225  SiStripQuality::Range range(badStripList.begin(),badStripList.end());
226  if ( ! OutSiStripQuality->put(detid,range) )
227  edm::LogError("SiStripHotStrip")<<"[SiStripHotStripAlgorithmFromClusterOccupancy::extractBadStrips] detid already exists"<<std::endl;
228  }
229  OutSiStripQuality->fillBadComponents();
230 
231  if (WriteOutputFile_==true){
232  f->cd();
233  striptree->Write();
234  f->Close();
235  }
236 
237  LogTrace("SiStripHotStrip") << ss.str() << std::endl;
238 }
239 
240 
241 void SiStripHotStripAlgorithmFromClusterOccupancy::iterativeSearch(pHisto& histo,std::vector<unsigned int>& vect, int apv){
242  if (!histo._NEntries || histo._NEntries <=MinNumEntries_ || histo._NEntries <= minNevents_)
243  return;
244 
245  size_t startingSize=vect.size();
246  long double diff=1.-prob_;
247 
248  size_t Nbins = histo._th1f->GetNbinsX();
249  size_t ibinStart = 1;
250  size_t ibinStop = Nbins+1;
251  int MaxEntry = (int)histo._th1f->GetMaximum();
252 
253  std::vector<long double> vPoissonProbs(MaxEntry+1,0);
254  long double meanVal=1.*histo._NEntries/(1.*Nbins-histo._NEmptyBins);
255  evaluatePoissonian(vPoissonProbs,meanVal);
256 
257  // Find median occupancy, taking into account only good strips
258  unsigned int goodstripentries[128];
259  int nGoodStrips = 0;
260  for (size_t i=ibinStart; i<ibinStop; ++i){
261  if (ishot[(apv*128)+i-1]==0){
262  goodstripentries[nGoodStrips] = (unsigned int)histo._th1f->GetBinContent(i);
263  nGoodStrips++;
264  }
265  }
266  double median = TMath::Median(nGoodStrips,goodstripentries);
267 
268  for (size_t i=ibinStart; i<ibinStop; ++i){
269  unsigned int entries= (unsigned int)histo._th1f->GetBinContent(i);
270 
271  if (ishot[(apv*128)+i-1]==0){
272  stripoccupancy[(apv*128)+i-1] = entries/(double) Nevents_;
273  striphits[(apv*128)+i-1] = entries;
274  poissonprob[(apv*128)+i-1] = 1-vPoissonProbs[entries];
275  medianapvhits[apv] = median;
276  avgapvhits[apv] = meanVal;
277  }
278  if (entries<=MinNumEntriesPerStrip_ || entries <= minNevents_ || entries / median < ratio_) continue;
279 
280  if(diff<vPoissonProbs[entries]){
281  ishot[(apv*128)+i-1] = 1;
283  hotstripsperapv[apv]++;
284  histo._th1f->SetBinContent(i,0.);
285  histo._NEntries-=entries;
286  histo._NEmptyBins++;
287  if (edm::isDebugEnabled())
288  LogTrace("SiStripHotStrip")<< " rejecting strip " << (apv*128)+i-1 << " value " << entries << " diff " << diff << " prob " << vPoissonProbs[entries]<< std::endl;
289  vect.push_back(pQuality->encode((apv*128)+i-1,1,0));
290  }
291 
292  }
293  if (edm::isDebugEnabled())
294  LogTrace("SiStripHotStrip") << " [SiStripHotStripAlgorithmFromClusterOccupancy::iterativeSearch] Nbins="<< Nbins << " MaxEntry="<<MaxEntry << " meanVal=" << meanVal << " NEmptyBins="<<histo._NEmptyBins<< " NEntries=" << histo._NEntries << " thEntries " << histo._th1f->GetEntries()<< " startingSize " << startingSize << " vector.size " << vect.size() << std::endl;
295 
296  if (vect.size()!=startingSize)
297  iterativeSearch(histo,vect,apv);
298 }
299 
300 void SiStripHotStripAlgorithmFromClusterOccupancy::evaluatePoissonian(std::vector<long double>& vPoissonProbs, long double& meanVal){
301  for(size_t i=0;i<vPoissonProbs.size();++i){
302  vPoissonProbs[i]= (i==0)?TMath::Poisson(i,meanVal):vPoissonProbs[i-1]+TMath::Poisson(i,meanVal);
303  }
304 }
305 
307  Nevents_=Nevents;
309  if (edm::isDebugEnabled())
310  LogTrace("SiStripHotStrip")<<" [SiStripHotStripAlgorithmFromClusterOccupancy::setNumberOfEvents] minNumber of Events per strip used to consider a strip bad" << minNevents_ << " for occupancy " << occupancy_ << std::endl;
311 }
bool isDebugEnabled()
int i
Definition: DBlmapReader.cc:9
SiStripHotStripAlgorithmFromClusterOccupancy(const edm::ParameterSet &, const TrackerTopology *)
void extractBadStrips(SiStripQuality *, HistoMap &, edm::ESHandle< SiStripQuality > &)
unsigned int tibLayer(const DetId &id) const
unsigned int tibString(const DetId &id) const
unsigned int tidRing(const DetId &id) const
bool tobIsStereo(const DetId &id) const
unsigned int tecRing(const DetId &id) const
ring id
T y() const
Definition: PV3DBase.h:63
void evaluatePoissonian(std::vector< long double > &, long double &meanVal)
unsigned int tidWheel(const DetId &id) const
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
bool tecIsZMinusSide(const DetId &id) const
bool tidIsStereo(const DetId &id) const
bool tidIsZMinusSide(const DetId &id) const
static const char occupancy_[]
bool tecIsStereo(const DetId &id) const
bool tibIsExternalString(const DetId &id) const
bool tibIsZMinusSide(const DetId &id) const
unsigned int tidModule(const DetId &id) const
bool tidIsBackRing(const DetId &id) const
T z() const
Definition: PV3DBase.h:64
void compact(unsigned int &, std::vector< unsigned int > &)
bool tobIsZMinusSide(const DetId &id) const
void fillBadComponents()
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
#define LogTrace(id)
unsigned int tibModule(const DetId &id) const
unsigned int tecModule(const DetId &id) const
bool tecIsBackPetal(const DetId &id) const
Definition: DetId.h:18
void iterativeSearch(pHisto &, std::vector< unsigned int > &, int)
bool tibIsStereo(const DetId &id) const
unsigned int tobModule(const DetId &id) const
std::pair< ContainerIterator, ContainerIterator > Range
virtual LocalPoint localPosition(float strip) const =0
unsigned int tecPetalNumber(const DetId &id) const
tuple cout
Definition: gather_cfg.py:121
volatile std::atomic< bool > shutdown_flag false
unsigned int tobRod(const DetId &id) const
bool put(const uint32_t &detID, const InputVector &vect)
T x() const
Definition: PV3DBase.h:62
unsigned int tecWheel(const DetId &id) const
unsigned int encode(const unsigned short &first, const unsigned short &NconsecutiveBadStrips, const unsigned short &flag=0)
unsigned int tobLayer(const DetId &id) const
virtual const TrackerGeomDet * idToDet(DetId) const