CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripCorrelateBadStripAndNoise.cc
Go to the documentation of this file.
2 
5 
6 
7 
9  cacheID_quality(0xFFFFFFFF),cacheID_noise(0xFFFFFFFF)
10 {
11  //now do what ever initialization is needed
12  if(!edm::Service<SiStripDetInfoFileReader>().isAvailable()){
13  edm::LogError("TkLayerMap") <<
14  "\n------------------------------------------"
15  "\nUnAvailable Service SiStripDetInfoFileReader: please insert in the configuration file an instance like"
16  "\n\tprocess.SiStripDetInfoFileReader = cms.Service(\"SiStripDetInfoFileReader\")"
17  "\n------------------------------------------";
18  }
19 
21  file = new TFile("correlTest.root","RECREATE");
22  tkmap = new TrackerMap();
23 }
24 
25 
27 {}
28 
29 //
30 
31 void
33 
35  return;
38 
39  edm::LogInfo("") << "[SiStripCorrelateBadStripAndNoise::beginRun] cacheID_quality " << cacheID_quality << " cacheID_noise " << cacheID_noise << std::endl;
40 
43 
44  DoAnalysis(es);
45 }
46 
47 void
49 
50  //Retrieve tracker topology from geometry
52  es.get<TrackerTopologyRcd>().get(tTopoHandle);
53  const TrackerTopology* const tTopo = tTopoHandle.product();
54 
55  //Loop on quality bad stirps
56  //for each strip, look at the noise
57  // evalaute the mean apv noise and the ratio among strip noise and meanApvNoise
58  // put the value in the histo in terms of ratio Vs percentage of badStrips per APV
59 
60  //Fill an histo per subdet and layer (and plus && minus for TEC/TID)
61  edm::LogInfo("") << "[Doanalysis]";
62  iterateOnDets(tTopo);
63 
64 }
65 
66 void
68 
69  SiStripQuality::RegistryIterator rbegin = qualityHandle_->getRegistryVectorBegin();
70  SiStripQuality::RegistryIterator rend = qualityHandle_->getRegistryVectorEnd();
71 
72 
73  for (SiStripBadStrip::RegistryIterator rp=rbegin; rp != rend; ++rp) {
74  const uint32_t detid=rp->detid;
75 
76  SiStripQuality::Range sqrange = SiStripQuality::Range( qualityHandle_->getDataVectorBegin()+rp->ibegin , qualityHandle_->getDataVectorBegin()+rp->iend );
77  iterateOnBadStrips(detid,tTopo,sqrange);
78  }
79 }
80 
81 void
83 
84  float percentage=0;
85  for(int it=0;it<sqrange.second-sqrange.first;it++){
86  unsigned int firstStrip=qualityHandle_->decode( *(sqrange.first+it) ).firstStrip;
87  unsigned int range=qualityHandle_->decode( *(sqrange.first+it) ).range;
88 
89  correlateWithNoise(detid,tTopo,firstStrip,range);
90 
91  edm::LogInfo("range")<< range;
92  percentage+=range;
93  }
94  if(percentage!=0)
95  percentage/=128.*fr->getNumberOfApvsAndStripLength(detid).first;
96  if(percentage>1)
97  edm::LogError("SiStripQualityStatistics") << "PROBLEM detid " << detid << " value " << percentage<< std::endl;
98 
99  //------- Global Statistics on percentage of bad components along the IOVs ------//
100  if(percentage!=0)
101  edm::LogInfo("")<< "percentage " << detid << " " << percentage;
102 }
103 
104 void
105 SiStripCorrelateBadStripAndNoise::correlateWithNoise(const uint32_t & detid, const TrackerTopology* tTopo, const uint32_t & firstStrip, const uint32_t & range){
106 
107  std::vector<TH2F *>histos;
108 
109  SiStripNoises::Range noiseRange = noiseHandle_->getRange(detid);
110  edm::LogInfo("Domenico") << "detid " << detid << " first " << firstStrip << " range " << range;
111  float meanAPVNoise=getMeanNoise(noiseRange,firstStrip/128,128);
112 
113  //float meanNoiseHotStrips=getMeanNoise(noiseRange,firstStrip,range);
114  for (size_t theStrip=firstStrip;theStrip<firstStrip+range;theStrip++){
115  float meanNoiseHotStrips=getMeanNoise(noiseRange,theStrip,1);
116 
117  //Get the histogram for this detid
118  getHistos(detid, tTopo, histos);
119  float yvalue=range<21?1.*range:21;
120 
121  for(size_t i=0;i<histos.size();++i)
122  histos[i]->Fill(meanNoiseHotStrips/meanAPVNoise-1.,yvalue);
123 
124  if(meanNoiseHotStrips/meanAPVNoise-1.<-0.3)
125  tkmap->fillc(detid,0xFF0000);
126  else
127  tkmap->fillc(detid,0x0000FF);
128  }
129 }
130 
131 float
132 SiStripCorrelateBadStripAndNoise::getMeanNoise(const SiStripNoises::Range& noiseRange,const uint32_t& firstStrip, const uint32_t& range){
133 
134  float mean=0;
135  for (size_t istrip=firstStrip;istrip<firstStrip+range;istrip++){
136  mean+=noiseHandle_->getNoise(istrip,noiseRange);
137  }
138  return mean/(1.*range);
139 }
140 
141 void
142 SiStripCorrelateBadStripAndNoise::getHistos(const uint32_t& detid, const TrackerTopology* tTopo, std::vector<TH2F*>& histos){
143 
144  histos.clear();
145 
146  int subdet=-999; int component=-999;
147  SiStripDetId a(detid);
148  if ( a.subdetId() == 3 ){
149  subdet=0;
150  component=tTopo->tibLayer(detid);
151  } else if ( a.subdetId() == 4 ) {
152  subdet=1;
153  component=tTopo->tidSide(detid)==2?tTopo->tidWheel(detid):tTopo->tidWheel(detid)+3;
154  } else if ( a.subdetId() == 5 ) {
155  subdet=2;
156  component=tTopo->tobLayer(detid);
157  } else if ( a.subdetId() == 6 ) {
158  subdet=3;
159  component=tTopo->tecSide(detid)==2?tTopo->tecWheel(detid):tTopo->tecWheel(detid)+9;
160  }
161 
162  int index=100+subdet*100+component;
163 
164 
165  histos.push_back(getHisto(subdet));
166  histos.push_back(getHisto(index));
167 
168 }
169 
170 TH2F*
172  if(vTH2.size()<index+1)
173  vTH2.resize(index+1,0);
174 
175  if(vTH2[index]==0){
176  char name[128];
177  sprintf(name,"%lu",index);
178  edm::LogInfo("")<<"[getHisto] creating index " << index << std::endl;
179  vTH2[index]=new TH2F(name,name,50,-2.,2.,21,0.5,21.5);
180  }
181 
182  return vTH2[index];
183 }
184 
185 void
187  for(size_t i=0;i<vTH2.size();i++)
188  if(vTH2[i]!=0)
189  vTH2[i]->Write();
190 
191  file->Write();
192  file->Close();
193 
194  tkmap->save(true,0,0,"testTkMap.png");
195 
196 }
197 
int i
Definition: DBlmapReader.cc:9
unsigned int tibLayer(const DetId &id) const
TH2F * getHisto(const long unsigned int &index)
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
void iterateOnBadStrips(const uint32_t &detid, const TrackerTopology *tTopo, SiStripQuality::Range &sqrange)
unsigned int tidWheel(const DetId &id) const
void getHistos(const uint32_t &detid, const TrackerTopology *tTopo, std::vector< TH2F * > &histos)
Registry::const_iterator RegistryIterator
unsigned long long getQualityCache(const edm::EventSetup &eSetup)
unsigned int tidSide(const DetId &id) const
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void iterateOnDets(const TrackerTopology *tTopo)
void save(bool print_total=true, float minval=0., float maxval=0., std::string s="svgmap.svg", int width=1500, int height=800)
Definition: TrackerMap.cc:699
unsigned long long getNoiseCache(const edm::EventSetup &eSetup)
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
edm::ESHandle< SiStripQuality > qualityHandle_
SiStripCorrelateBadStripAndNoise(const edm::ParameterSet &)
void fillc(int idmod, int RGBcode)
Definition: TrackerMap.h:104
float getMeanNoise(const SiStripNoises::Range &noiseRange, const uint32_t &first, const uint32_t &range)
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
virtual void beginRun(const edm::Run &run, const edm::EventSetup &es)
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
double a
Definition: hdecay.h:121
std::pair< ContainerIterator, ContainerIterator > Range
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:48
unsigned int tecWheel(const DetId &id) const
unsigned int tobLayer(const DetId &id) const
Definition: Run.h:43
void correlateWithNoise(const uint32_t &detid, const TrackerTopology *tTopo, const uint32_t &firstStrip, const uint32_t &range)
unsigned int tecSide(const DetId &id) const