CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PVValidationHelpers.cc
Go to the documentation of this file.
3 #include <cassert>
4 #include <numeric>
5 #include <algorithm>
6 #include <iostream>
7 #include "TMath.h"
8 #include "TH1.h"
9 #include "TF1.h"
10 
11 //*************************************************************
12 void PVValHelper::add(std::map<std::string, TH1*>& h, TH1* hist)
13 //*************************************************************
14 {
15  h[hist->GetName()] = hist;
16  hist->StatOverflows(kTRUE);
17 }
18 
19 //*************************************************************
20 void PVValHelper::fill(std::map<std::string, TH1*>& h, const std::string& s, double x)
21 //*************************************************************
22 {
23  if (h.count(s) == 0) {
24  edm::LogWarning("PVValidationHelpers") << "Trying to fill non-existing Histogram named " << s << std::endl;
25  return;
26  }
27  h[s]->Fill(x);
28 }
29 
30 //*************************************************************
31 void PVValHelper::fill(std::map<std::string, TH1*>& h, const std::string& s, double x, double y)
32 //*************************************************************
33 {
34  if (h.count(s) == 0) {
35  edm::LogWarning("PVValidationHelpers") << "Trying to fill non-existing Histogram named " << s << std::endl;
36  return;
37  }
38  h[s]->Fill(x, y);
39 }
40 
41 //*************************************************************
42 void PVValHelper::fillByIndex(std::vector<TH1F*>& h, unsigned int index, double x, std::string tag)
43 //*************************************************************
44 {
45  assert(!h.empty());
46  if (index < h.size()) {
47  h[index]->Fill(x);
48  } else {
49  edm::LogWarning("PVValidationHelpers") << "Trying to fill non-existing Histogram with index " << index
50  << " for array with size: " << h.size() << " tag: " << tag << std::endl;
51  return;
52  }
53 }
54 
55 //*************************************************************
56 void PVValHelper::shrinkHistVectorToFit(std::vector<TH1F*>& h, unsigned int desired_size)
57 //*************************************************************
58 {
59  h.erase(h.begin() + desired_size, h.end());
60 }
61 
62 //*************************************************************
64 //*************************************************************
65 {
67  switch (type) {
68  // absoulte
69 
70  case PVValHelper::dxy:
71  returnType = std::make_tuple("dxy", "d_{xy}", "[#mum]");
72  break;
73  case PVValHelper::dx:
74  returnType = std::make_tuple("dx", "d_{x}", "[#mum]");
75  break;
76  case PVValHelper::dy:
77  returnType = std::make_tuple("dy", "d_{y}", "[#mum]");
78  break;
79  case PVValHelper::dz:
80  returnType = std::make_tuple("dz", "d_{z}", "[#mum]");
81  break;
82  case PVValHelper::IP2D:
83  returnType = std::make_tuple("IP2D", "IP_{2D}", "[#mum]");
84  break;
85  case PVValHelper::resz:
86  returnType = std::make_tuple("resz", "z_{trk}-z_{vtx}", "[#mum]");
87  break;
88  case PVValHelper::IP3D:
89  returnType = std::make_tuple("IP3D", "IP_{3D}", "[#mum]");
90  break;
91  case PVValHelper::d3D:
92  returnType = std::make_tuple("d3D", "d_{3D}", "[#mum]");
93  break;
94 
95  // normalized
96 
98  returnType = std::make_tuple("norm_dxy", "d_{xy}/#sigma_{d_{xy}}", "");
99  break;
101  returnType = std::make_tuple("norm_dx", "d_{x}/#sigma_{d_{x}}", "");
102  break;
104  returnType = std::make_tuple("norm_dy", "d_{y}/#sigma_{d_{y}}", "");
105  break;
107  returnType = std::make_tuple("norm_dz", "d_{z}/#sigma_{d_{z}}", "");
108  break;
110  returnType = std::make_tuple("norm_IP2D", "IP_{2D}/#sigma_{IP_{2D}}", "");
111  break;
113  returnType = std::make_tuple("norm_resz", "z_{trk}-z_{vtx}/#sigma_{res_{z}}", "");
114  break;
116  returnType = std::make_tuple("norm_IP3D", "IP_{3D}/#sigma_{IP_{3D}}", "");
117  break;
119  returnType = std::make_tuple("norm_d3D", "d_{3D}/#sigma_{d_{3D}}", "");
120  break;
121 
122  default:
123  edm::LogWarning("PVValidationHelpers") << " getTypeString() unknown residual type: " << type << std::endl;
124  }
125 
126  return returnType;
127 }
128 
129 //*************************************************************
131 //*************************************************************
132 {
133  PVValHelper::plotLabels returnVar;
134  switch (var) {
135  case PVValHelper::phi:
136  returnVar = std::make_tuple("phi", "#phi", "[rad]");
137  break;
138  case PVValHelper::eta:
139  returnVar = std::make_tuple("eta", "#eta", "");
140  break;
141  case PVValHelper::pT:
142  returnVar = std::make_tuple("pT", "p_{T}", "[GeV]");
143  break;
145  returnVar = std::make_tuple("pTCentral", "p_{T} |#eta|<1.", "[GeV]");
146  break;
147  case PVValHelper::ladder:
148  returnVar = std::make_tuple("ladder", "ladder number", "");
149  break;
150  case PVValHelper::modZ:
151  returnVar = std::make_tuple("modZ", "module number", "");
152  break;
153  default:
154  edm::LogWarning("PVValidationHelpers") << " getVarString() unknown plot variable: " << var << std::endl;
155  }
156 
157  return returnVar;
158 }
159 
160 //*************************************************************
161 std::vector<float> PVValHelper::generateBins(int n, float start, float range)
162 //*************************************************************
163 {
164  std::vector<float> v(n);
165  float interval = range / (n - 1);
166  std::iota(v.begin(), v.end(), 1.);
167 
168  /*
169  std::cout<<" interval:"<<interval<<std::endl;
170  for(float &a : v) { std::cout<< a << " "; }
171  std::cout<< "\n";
172  */
173 
174  std::for_each(begin(v), end(v), [&](float& a) { a = start + ((a - 1) * interval); });
175 
176  return v;
177 }
178 
179 //*************************************************************
181 //*************************************************************
182 {
183  double median = 0.;
184  double q = 0.5; // 0.5 quantile for "median"
185  // protect against empty histograms
186  if (histo->Integral() != 0) {
187  histo->GetQuantiles(1, &median, &q);
188  }
189 
190  Measurement1D result(median, median / TMath::Sqrt(histo->GetEntries()));
191 
192  return result;
193 }
194 
195 //*************************************************************
197 //*************************************************************
198 {
199  int nbins = histo->GetNbinsX();
200  double median = getMedian(histo).value();
201  double x_lastBin = histo->GetBinLowEdge(nbins + 1);
202  const char* HistoName = histo->GetName();
203  std::string Finalname = "resMed_";
204  Finalname.append(HistoName);
205  TH1F* newHisto = new TH1F(Finalname.c_str(), Finalname.c_str(), nbins, 0., x_lastBin);
206  double* residuals = new double[nbins];
207  const float* weights = histo->GetArray();
208 
209  for (int j = 0; j < nbins; j++) {
210  residuals[j] = std::abs(median - histo->GetBinCenter(j + 1));
211  newHisto->Fill(residuals[j], weights[j]);
212  }
213 
214  double theMAD = (PVValHelper::getMedian(newHisto).value()) * 1.4826;
215 
216  delete[] residuals;
217  residuals = nullptr;
218  newHisto->Delete("");
219 
220  Measurement1D result(theMAD, theMAD / histo->GetEntries());
221  return result;
222 }
223 
224 //*************************************************************
225 std::pair<Measurement1D, Measurement1D> PVValHelper::fitResiduals(TH1* hist)
226 //*************************************************************
227 {
228  //float fitResult(9999);
229  if (hist->GetEntries() < 1) {
230  return std::make_pair(Measurement1D(0., 0.), Measurement1D(0., 0.));
231  };
232 
233  float mean = hist->GetMean();
234  float sigma = hist->GetRMS();
235 
236  TF1 func("tmp", "gaus", mean - 1.5 * sigma, mean + 1.5 * sigma);
237  if (0 == hist->Fit(&func, "QNR")) { // N: do not blow up file by storing fit!
238  mean = func.GetParameter(1);
239  sigma = func.GetParameter(2);
240  // second fit: three sigma of first fit around mean of first fit
241  func.SetRange(mean - 2 * sigma, mean + 2 * sigma);
242  // I: integral gives more correct results if binning is too wide
243  // L: Likelihood can treat empty bins correctly (if hist not weighted...)
244  if (0 == hist->Fit(&func, "Q0LR")) {
245  if (hist->GetFunction(func.GetName())) { // Take care that it is later on drawn:
246  hist->GetFunction(func.GetName())->ResetBit(TF1::kNotDraw);
247  }
248  }
249  }
250 
251  float res_mean = func.GetParameter(1);
252  float res_width = func.GetParameter(2);
253 
254  float res_mean_err = func.GetParError(1);
255  float res_width_err = func.GetParError(2);
256 
257  Measurement1D resultM(res_mean, res_mean_err);
258  Measurement1D resultW(res_width, res_width_err);
259 
260  std::pair<Measurement1D, Measurement1D> result;
261 
262  result = std::make_pair(resultM, resultW);
263  return result;
264 }
Measurement1D getMedian(TH1F *histo)
void fillByIndex(std::vector< TH1F * > &h, unsigned int index, double x, std::string tag="")
std::pair< Measurement1D, Measurement1D > fitResiduals(TH1 *hist)
tuple interval
Definition: MergeJob_cfg.py:20
std::vector< float > generateBins(int n, float start, float range)
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t Func __host__ __device__ V int Func func
assert(be >=bs)
Measurement1D getMAD(TH1F *histo)
plotLabels getVarString(plotVariable var)
const uint16_t range(const Frame &aFrame)
tuple result
Definition: mps_fire.py:311
edm::TypeWithDict returnType(const edm::FunctionWithDict &)
Definition: returnType.cc:16
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
__shared__ Hist hist
void add(std::map< std::string, TH1 * > &h, TH1 *hist)
plotLabels getTypeString(residualType type)
void fill(std::map< std::string, TH1 * > &h, const std::string &s, double x)
void shrinkHistVectorToFit(std::vector< TH1F * > &h, unsigned int desired_size)
double value() const
Definition: Measurement1D.h:25
std::string HistoName
T median(std::vector< T > values)
Definition: median.h:16
double a
Definition: hdecay.h:119
string end
Definition: dataset.py:937
float x
std::tuple< std::string, std::string, std::string > plotLabels
Log< level::Warning, false > LogWarning
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4