CMS 3D CMS Logo

SiStripPayloadInspectorHelper.h
Go to the documentation of this file.
1 #ifndef CONDCORE_SISTRIPPLUGINS_SISTRIPPAYLOADINSPECTORHELPER_H
2 #define CONDCORE_SISTRIPPLUGINS_SISTRIPPAYLOADINSPECTORHELPER_H
3 
4 #include <vector>
5 #include <numeric>
6 #include <string>
7 #include "TH1.h"
8 #include "TPaveText.h"
10 
11 namespace SiStripPI {
12 
13  enum TrackerRegion {
14  TIB1r = 1010, TIB1s = 1011,
15  TIB2r = 1020, TIB2s = 1021,
16  TIB3r = 1030,
17  TIB4r = 1040,
18  TOB1r = 2010, TOB1s = 2011,
19  TOB2r = 2020, TOB2s = 2021,
20  TOB3r = 2030,
21  TOB4r = 2040,
22  TOB5r = 2050,
23  TOB6r = 2060,
24  TEC1r = 3010, TEC1s = 3011,
25  TEC2r = 3020, TEC2s = 3021,
26  TEC3r = 3030, TEC3s = 3031,
27  TEC4r = 3040, TEC4s = 3041,
28  TEC5r = 3050, TEC5s = 3051,
29  TEC6r = 3060, TEC6s = 3061,
30  TEC7r = 3070, TEC7s = 3071,
31  TEC8r = 3080, TEC8s = 3081,
32  TEC9r = 3090, TEC9s = 3091,
33  TID1r = 4010, TID1s = 4011,
34  TID2r = 4020, TID2s = 4021,
35  TID3r = 4030, TID3s = 4031,
37  };
38 
39  /*--------------------------------------------------------------------*/
40  const char * regionType(int index)
41  /*--------------------------------------------------------------------*/
42  {
43 
44  auto region = static_cast<std::underlying_type_t<SiStripPI::TrackerRegion> >(index);
45 
46  switch(region){
47  case SiStripPI::TIB1r: return "TIB L1 r-#varphi";
48  case SiStripPI::TIB1s: return "TIB L1 stereo";
49  case SiStripPI::TIB2r: return "TIB L2 r-#varphi";
50  case SiStripPI::TIB2s: return "TIB L2 stereo";
51  case SiStripPI::TIB3r: return "TIB L3";
52  case SiStripPI::TIB4r: return "TIB L4";
53  case SiStripPI::TOB1r: return "TOB L1 r-#varphi";
54  case SiStripPI::TOB1s: return "TOB L1 stereo";
55  case SiStripPI::TOB2r: return "TOB L2 r-#varphi";
56  case SiStripPI::TOB2s: return "TOB L2 stereo";
57  case SiStripPI::TOB3r: return "TOB L3 r-#varphi";
58  case SiStripPI::TOB4r: return "TOB L4";
59  case SiStripPI::TOB5r: return "TOB L5";
60  case SiStripPI::TOB6r: return "TOB L6";
61  case SiStripPI::TEC1r: return "TEC D1 r-#varphi";
62  case SiStripPI::TEC1s: return "TEC D1 stereo";
63  case SiStripPI::TEC2r: return "TEC D2 r-#varphi";
64  case SiStripPI::TEC2s: return "TEC D2 stereo";
65  case SiStripPI::TEC3r: return "TEC D3 r-#varphi";
66  case SiStripPI::TEC3s: return "TEC D3 stereo";
67  case SiStripPI::TEC4r: return "TEC D4 r-#varphi";
68  case SiStripPI::TEC4s: return "TEC D4 stereo";
69  case SiStripPI::TEC5r: return "TEC D5 r-#varphi";
70  case SiStripPI::TEC5s: return "TEC D5 stereo";
71  case SiStripPI::TEC6r: return "TEC D6 r-#varphi";
72  case SiStripPI::TEC6s: return "TEC D6 stereo";
73  case SiStripPI::TEC7r: return "TEC D7 r-#varphi";
74  case SiStripPI::TEC7s: return "TEC D7 stereo";
75  case SiStripPI::TEC8r: return "TEC D8 r-#varphi";
76  case SiStripPI::TEC8s: return "TEC D8 stereo";
77  case SiStripPI::TEC9r: return "TEC D9 r-#varphi";
78  case SiStripPI::TEC9s: return "TEC D9 stereo";
79  case SiStripPI::TID1r: return "TID D1 r-#varphi";
80  case SiStripPI::TID1s: return "TID D1 stereo";
81  case SiStripPI::TID2r: return "TID D2 r-#varphi";
82  case SiStripPI::TID2s: return "TID D2 stereo";
83  case SiStripPI::TID3r: return "TID D3 r-#varphi";
84  case SiStripPI::TID3s: return "TID D3 stereo";
85  case SiStripPI::END_OF_REGIONS : return "undefined";
86  default : return "should never be here";
87  }
88  }
89 
90  /*--------------------------------------------------------------------*/
91  std::pair<float,float> getTheRange(std::map<uint32_t,float> values)
92  /*--------------------------------------------------------------------*/
93  {
94  float sum = std::accumulate(std::begin(values),
95  std::end(values),
96  0.0,
98  { return value + p.second; }
99  );
100 
101  float m = sum / values.size();
102 
103  float accum = 0.0;
104  std::for_each (std::begin(values),
105  std::end(values),
107  {accum += (p.second - m) * (p.second - m);}
108  );
109 
110  float stdev = sqrt(accum / (values.size()-1));
111 
112  return std::make_pair(m-2*stdev,m+2*stdev);
113 
114  }
115 
116  /*--------------------------------------------------------------------*/
117  void drawStatBox(std::map<std::string,std::shared_ptr<TH1F>> histos, std::map<std::string,int> colormap, std::vector<std::string> legend, double X=0.15, double Y=0.93, double W=0.15, double H=0.10)
118  /*--------------------------------------------------------------------*/
119  {
120  char buffer[255];
121 
122  int i=0;
123  for ( const auto &element : legend ){
124  TPaveText* stat = new TPaveText(X,Y-(i*H), X+W, Y-(i+1)*H, "NDC");
125  i++;
126  auto Histo = histos[element];
127  sprintf(buffer,"Entries : %i\n",(int)Histo->GetEntries());
128  stat->AddText(buffer);
129 
130  sprintf(buffer,"Mean : %6.2f\n",Histo->GetMean());
131  stat->AddText(buffer);
132 
133  sprintf(buffer,"RMS : %6.2f\n",Histo->GetRMS());
134  stat->AddText(buffer);
135 
136  stat->SetFillColor(0);
137  stat->SetLineColor(colormap[element]);
138  stat->SetTextColor(colormap[element]);
139  stat->SetTextSize(0.03);
140  stat->SetBorderSize(0);
141  stat->SetMargin(0.05);
142  stat->SetTextAlign(12);
143  stat->Draw();
144  }
145  }
146 
147  /*--------------------------------------------------------------------*/
148  std::pair<float,float> getExtrema(TH1 *h1,TH1 *h2)
149  /*--------------------------------------------------------------------*/
150  {
151  float theMax(-9999.);
152  float theMin(9999.);
153  theMax = h1->GetMaximum() > h2->GetMaximum() ? h1->GetMaximum() : h2->GetMaximum();
154  theMin = h1->GetMinimum() < h2->GetMaximum() ? h1->GetMinimum() : h2->GetMinimum();
155 
156  float add_min = theMin>0. ? -0.05 : 0.05;
157  float add_max = theMax>0. ? 0.05 : -0.05;
158 
159  auto result = std::make_pair(theMin*(1+add_min),theMax*(1+add_max));
160  return result;
161 
162  }
163 
164 
165  /*--------------------------------------------------------------------*/
167  /*--------------------------------------------------------------------*/
168  {
169  hist->SetStats(kFALSE);
170  hist->SetLineWidth(2);
171  hist->GetXaxis()->CenterTitle(true);
172  hist->GetYaxis()->CenterTitle(true);
173  hist->GetXaxis()->SetTitleFont(42);
174  hist->GetYaxis()->SetTitleFont(42);
175  hist->GetXaxis()->SetTitleSize(0.05);
176  hist->GetYaxis()->SetTitleSize(0.05);
177  hist->GetXaxis()->SetTitleOffset(0.9);
178  hist->GetYaxis()->SetTitleOffset(1.3);
179  hist->GetXaxis()->SetLabelFont(42);
180  hist->GetYaxis()->SetLabelFont(42);
181  hist->GetYaxis()->SetLabelSize(.05);
182  hist->GetXaxis()->SetLabelSize(.05);
183  }
184 
185 
186  /*--------------------------------------------------------------------*/
187  void printSummary(const std::map<unsigned int, SiStripDetSummary::Values>& map)
188  /*--------------------------------------------------------------------*/
189  {
190  for (const auto &element : map){
191  int count = element.second.count;
192  double mean = count>0 ? (element.second.mean)/count : 0. ;
193  double rms = count>0 ? (element.second.rms)/count - mean*mean : 0.;
194  if(rms <= 0)
195  rms = 0;
196  else
197  rms = sqrt(rms);
198 
200 
201  switch ((element.first)/1000)
202  {
203  case 1:
204  detector = "TIB ";
205  break;
206  case 2:
207  detector = "TOB ";
208  break;
209  case 3:
210  detector = "TEC ";
211  break;
212  case 4:
213  detector = "TID ";
214  break;
215  }
216 
217  int layer = (element.first)/10 - (element.first)/1000*100;
218  int stereo = (element.first) - (layer*10) -(element.first)/1000*1000;
219 
220  std::cout<<"key of the map:"<<element.first <<" ( region: "<<regionType(element.first) <<" ) "
221  << detector<<" layer: "<<layer<<" stereo:"<<stereo
222  <<"| count:"<<count<<" mean: "<<mean<<" rms: "<<rms<<std::endl;
223 
224  }
225  }
226 };
227 
228 #endif
const char * regionType(int index)
std::pair< float, float > getTheRange(std::map< uint32_t, float > values)
#define X(str)
Definition: MuonsGrabber.cc:48
void printSummary(const std::map< unsigned int, SiStripDetSummary::Values > &map)
void drawStatBox(std::map< std::string, std::shared_ptr< TH1F >> histos, std::map< std::string, int > colormap, std::vector< std::string > legend, double X=0.15, double Y=0.93, double W=0.15, double H=0.10)
T sqrt(T t)
Definition: SSEVec.h:18
std::pair< float, float > getExtrema(TH1 *h1, TH1 *h2)
void makeNicePlotStyle(TH1 *hist)
#define end
Definition: vmac.h:37
Definition: value.py:1
#define begin
Definition: vmac.h:30
def stdev(xlist)
Definition: plotscripts.py:67