CMS 3D CMS Logo

SiStripPedestals_PayloadInspector.cc
Go to the documentation of this file.
1 
12 
13 // the data format of the condition to be inspected
19 
20 // needed for the tracker map
22 
23 // auxilliary functions
26 
27 #include <memory>
28 #include <sstream>
29 #include <iostream>
30 
31 // include ROOT
32 #include "TH2F.h"
33 #include "TLegend.h"
34 #include "TCanvas.h"
35 #include "TLine.h"
36 #include "TStyle.h"
37 #include "TLatex.h"
38 #include "TPave.h"
39 #include "TPaveStats.h"
40 
41 namespace {
42 
43  /************************************************
44  test class
45  *************************************************/
46 
47  class SiStripPedestalsTest : public cond::payloadInspector::Histogram1D<SiStripPedestals> {
48 
49  public:
50  SiStripPedestalsTest() : cond::payloadInspector::Histogram1D<SiStripPedestals>("SiStrip Pedestals test",
51  "SiStrip Pedestals test", 10,0.0,10.0),
52  m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())}
53  {
54  Base::setSingleIov( true );
55  }
56 
57  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
58  for ( auto const & iov: iovs) {
59  std::shared_ptr<SiStripPedestals> payload = Base::fetchPayload( std::get<1>(iov) );
60  if( payload.get() ){
61 
62  fillWithValue(1.);
63 
64  std::stringstream ss;
65  ss << "Summary of strips pedestals:" << std::endl;
66 
67  //payload->printDebug(ss);
68  payload->printSummary(ss,&m_trackerTopo);
69 
70  std::vector<uint32_t> detid;
71  payload->getDetIds(detid);
72 
73  // for (const auto & d : detid) {
74  // int nstrip=0;
75  // SiStripPedestals::Range range=payload->getRange(d);
76  // for( int it=0; it < (range.second-range.first)*8/10; ++it ){
77  // auto ped = payload->getPed(it,range);
78  // nstrip++;
79  // ss << "DetId="<< d << " Strip=" << nstrip <<": "<< ped << std::endl;
80  // } // end of loop on strips
81  // } // end of loop on detIds
82 
83  std::cout<<ss.str()<<std::endl;
84 
85  }// payload
86  }// iovs
87  return true;
88  }// fill
89  private:
90  TrackerTopology m_trackerTopo;
91  };
92 
93  /************************************************
94  1d histogram of SiStripNoises of 1 IOV
95  *************************************************/
96 
97  // inherit from one of the predefined plot class: Histogram1D
98  class SiStripPedestalsValue : public cond::payloadInspector::Histogram1D<SiStripPedestals> {
99 
100  public:
101  SiStripPedestalsValue() : cond::payloadInspector::Histogram1D<SiStripPedestals>("SiStrip Pedestals values",
102  "SiStrip Pedestals values", 300,0.0,300.0){
103  Base::setSingleIov( true );
104  }
105 
106  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
107  for ( auto const & iov: iovs) {
108  std::shared_ptr<SiStripPedestals> payload = Base::fetchPayload( std::get<1>(iov) );
109  if( payload.get() ){
110 
111  std::vector<uint32_t> detid;
112  payload->getDetIds(detid);
113 
114  for (const auto & d : detid) {
115  SiStripPedestals::Range range=payload->getRange(d);
116  for( int it=0; it < (range.second-range.first)*8/10; ++it ){
117  auto ped = payload->getPed(it,range);
118  //to be used to fill the histogram
119  fillWithValue(ped);
120  }// loop over APVs
121  } // loop over detIds
122  }// payload
123  }// iovs
124  return true;
125  }// fill
126  };
127 
128  /************************************************
129  Tracker Map of SiStrip Pedestals
130  *************************************************/
131 
132  template<SiStripPI::estimator est> class SiStripPedestalsTrackerMap : public cond::payloadInspector::PlotImage<SiStripPedestals> {
133 
134  public:
135  SiStripPedestalsTrackerMap() : cond::payloadInspector::PlotImage<SiStripPedestals> ( "Tracker Map of SiStripPedestals "+estimatorType(est)+" per module" )
136  {
137  setSingleIov( true );
138  }
139 
140  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override {
141  auto iov = iovs.front();
142  std::shared_ptr<SiStripPedestals> payload = fetchPayload( std::get<1>(iov) );
143 
144  std::string titleMap = "Tracker Map of SiStrip Pedestals "+estimatorType(est)+" per module (payload : "+std::get<1>(iov)+")";
145 
146  std::unique_ptr<TrackerMap> tmap = std::unique_ptr<TrackerMap>(new TrackerMap("SiStripPedestals"));
147  tmap->setTitle(titleMap);
148  tmap->setPalette(1);
149 
150  std::vector<uint32_t> detid;
151  payload->getDetIds(detid);
152 
153  std::map<unsigned int,float> info_per_detid;
154 
155  for (const auto & d : detid) {
156  int nstrips=0;
157  double mean(0.),rms(0.),min(10000.), max(0.);
158  SiStripPedestals::Range range = payload->getRange(d);
159  for( int it=0; it < (range.second-range.first)*8/10; ++it ){
160  nstrips++;
161  auto ped = payload->getPed(it,range);
162  mean+=ped;
163  rms+=(ped*ped);
164  if(ped<min) min=ped;
165  if(ped>max) max=ped;
166  } // end of loop on strips
167 
168  mean/=nstrips;
169  if((rms/nstrips-mean*mean)>0.){
170  rms = sqrt(rms/nstrips-mean*mean);
171  } else {
172  rms=0.;
173  }
174 
175  switch(est){
176  case SiStripPI::min:
177  info_per_detid[d]=min;
178  break;
179  case SiStripPI::max:
180  info_per_detid[d]=max;
181  break;
182  case SiStripPI::mean:
183  info_per_detid[d]=mean;
184  break;
185  case SiStripPI::rms:
186  info_per_detid[d]=rms;
187  break;
188  default:
189  edm::LogWarning("LogicError") << "Unknown estimator: " << est;
190  break;
191  }
192  } // end of loop on detids
193 
194  for(const auto & d : detid){
195  tmap->fill(d,info_per_detid[d]);
196  }
197 
198  std::string fileName(m_imageFileName);
199  tmap->save(true,0.,0.,fileName);
200 
201  return true;
202  }
203  };
204 
205  typedef SiStripPedestalsTrackerMap<SiStripPI::min> SiStripPedestalsMin_TrackerMap;
206  typedef SiStripPedestalsTrackerMap<SiStripPI::max> SiStripPedestalsMax_TrackerMap;
207  typedef SiStripPedestalsTrackerMap<SiStripPI::mean> SiStripPedestalsMean_TrackerMap;
208  typedef SiStripPedestalsTrackerMap<SiStripPI::rms> SiStripPedestalsRMS_TrackerMap;
209 
210  /************************************************
211  Tracker Map of SiStrip Pedestals Summaries
212  *************************************************/
213 
214  template<SiStripPI::estimator est> class SiStripPedestalsByRegion : public cond::payloadInspector::PlotImage<SiStripPedestals> {
215  public:
216  SiStripPedestalsByRegion() : cond::payloadInspector::PlotImage<SiStripPedestals>( "SiStrip Pedestals "+estimatorType(est)+" by Region" ),
217  m_trackerTopo{StandaloneTrackerTopology::fromTrackerParametersXMLFile(edm::FileInPath("Geometry/TrackerCommonData/data/trackerParameters.xml").fullPath())}
218  {
219  setSingleIov( true );
220  }
221 
222  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override {
223  auto iov = iovs.front();
224  std::shared_ptr<SiStripPedestals> payload = fetchPayload( std::get<1>(iov) );
225 
226  SiStripDetSummary summaryPedestals{&m_trackerTopo};
227  std::vector<uint32_t> detid;
228  payload->getDetIds(detid);
229 
230  for (const auto & d : detid) {
231  int nstrips=0;
232  double mean(0.),rms(0.),min(10000.), max(0.);
233  SiStripPedestals::Range range = payload->getRange(d);
234  for( int it=0; it < (range.second-range.first)*8/10; ++it ){
235  nstrips++;
236  auto ped = payload->getPed(it,range);
237  mean+=ped;
238  rms+=(ped*ped);
239  if(ped<min) min=ped;
240  if(ped>max) max=ped;
241  } // end of loop on strips
242 
243  mean/=nstrips;
244  if((rms/nstrips-mean*mean)>0.){
245  rms = sqrt(rms/nstrips-mean*mean);
246  } else {
247  rms=0.;
248  }
249 
250  switch(est){
251  case SiStripPI::min:
252  summaryPedestals.add(d,min);
253  break;
254  case SiStripPI::max:
255  summaryPedestals.add(d,max);
256  break;
257  case SiStripPI::mean:
258  summaryPedestals.add(d,mean);
259  break;
260  case SiStripPI::rms:
261  summaryPedestals.add(d,rms);
262  break;
263  default:
264  edm::LogWarning("LogicError") << "Unknown estimator: " << est;
265  break;
266  }
267  } // loop on the detIds
268 
269  std::map<unsigned int, SiStripDetSummary::Values> map = summaryPedestals.getCounts();
270  //=========================
271 
272  TCanvas canvas("Partion summary","partition summary",1200,1000);
273  canvas.cd();
274  auto h1 = std::unique_ptr<TH1F>(new TH1F("byRegion",Form("Average by partition of %s SiStrip Pedestals per module;;average SiStrip Pedestals %s [ADC counts]",estimatorType(est).c_str(),estimatorType(est).c_str()),map.size(),0.,map.size()));
275  h1->SetStats(false);
276  canvas.SetBottomMargin(0.18);
277  canvas.SetLeftMargin(0.17);
278  canvas.SetRightMargin(0.05);
279  canvas.Modified();
280 
281  std::vector<int> boundaries;
282  unsigned int iBin=0;
283 
285  std::string currentDetector;
286 
287  for (const auto &element : map){
288  iBin++;
289  int count = element.second.count;
290  double mean = (element.second.mean)/count;
291  double rms = (element.second.rms)/count - mean*mean;
292 
293  if(rms <= 0)
294  rms = 0;
295  else
296  rms = sqrt(rms);
297 
298  if(currentDetector.empty()) currentDetector="TIB";
299 
300  switch ((element.first)/1000)
301  {
302  case 1:
303  detector = "TIB";
304  break;
305  case 2:
306  detector = "TOB";
307  break;
308  case 3:
309  detector = "TEC";
310  break;
311  case 4:
312  detector = "TID";
313  break;
314  }
315 
316  h1->SetBinContent(iBin,mean);
317  h1->GetXaxis()->SetBinLabel(iBin,SiStripPI::regionType(element.first).second);
318  h1->GetXaxis()->LabelsOption("v");
319 
320  if(detector!=currentDetector) {
321  boundaries.push_back(iBin);
322  currentDetector=detector;
323  }
324  }
325 
326  h1->SetMarkerStyle(20);
327  h1->SetMarkerSize(1);
328  h1->SetMaximum(h1->GetMaximum()*1.1);
329  h1->Draw("HIST");
330  h1->Draw("Psame");
331 
332  canvas.Update();
333 
334  TLine l[boundaries.size()];
335  unsigned int i=0;
336  for (const auto & line : boundaries){
337  l[i] = TLine(h1->GetBinLowEdge(line),canvas.GetUymin(),h1->GetBinLowEdge(line),canvas.GetUymax());
338  l[i].SetLineWidth(1);
339  l[i].SetLineStyle(9);
340  l[i].SetLineColor(2);
341  l[i].Draw("same");
342  i++;
343  }
344 
345  TLegend legend = TLegend(0.52,0.82,0.95,0.9);
346  legend.SetHeader((std::get<1>(iov)).c_str(),"C"); // option "C" allows to center the header
347  legend.AddEntry(h1.get(),("IOV: "+std::to_string(std::get<0>(iov))).c_str(),"PL");
348  legend.SetTextSize(0.025);
349  legend.Draw("same");
350 
351  std::string fileName(m_imageFileName);
352  canvas.SaveAs(fileName.c_str());
353 
354  return true;
355  }
356  private:
357  TrackerTopology m_trackerTopo;
358  };
359 
360 
361  typedef SiStripPedestalsByRegion<SiStripPI::mean> SiStripPedestalsMeanByRegion;
362  typedef SiStripPedestalsByRegion<SiStripPI::min> SiStripPedestalsMinByRegion;
363  typedef SiStripPedestalsByRegion<SiStripPI::max> SiStripPedestalsMaxByRegion;
364  typedef SiStripPedestalsByRegion<SiStripPI::rms> SiStripPedestalsRMSByRegion;
365 
366 } // close namespace
367 
369  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsTest);
370  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsValue);
371  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMin_TrackerMap);
372  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMax_TrackerMap);
373  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMean_TrackerMap);
374  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsRMS_TrackerMap);
375  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMeanByRegion);
376  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMinByRegion);
377  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsMaxByRegion);
378  PAYLOAD_INSPECTOR_CLASS(SiStripPedestalsRMSByRegion);
379 }
std::pair< ContainerIterator, ContainerIterator > Range
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
T sqrt(T t)
Definition: SSEVec.h:18
void fillWithValue(float value, float weight=1)
T min(T a, T b)
Definition: MathUtil.h:58
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
std::string estimatorType(SiStripPI::estimator e)
std::pair< int, const char * > regionType(int index)
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
Definition: plugin.cc:24
bool fill(const std::vector< std::tuple< cond::Time_t, cond::Hash > > &iovs) override
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
def canvas(sub, attr)
Definition: svgfig.py:481