CMS 3D CMS Logo

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  //std::cout<<" interval:"<<interval<<std::endl;
169  //for(float &a : v) { std::cout<< a << " "; }
170  //std::cout<< "\n";
171 
172  std::for_each(begin(v), end(v), [&](float& a) { a = start + ((a - 1) * interval); });
173 
174  return v;
175 }
176 
177 //*************************************************************
179 //*************************************************************
180 {
181  double median = 0.;
182  double q = 0.5; // 0.5 quantile for "median"
183  // protect against empty histograms
184  if (histo->Integral() != 0) {
185  histo->GetQuantiles(1, &median, &q);
186  }
187 
188  Measurement1D result(median, median / TMath::Sqrt(histo->GetEntries()));
189 
190  return result;
191 }
192 
193 //*************************************************************
195 //*************************************************************
196 {
197  int nbins = histo->GetNbinsX();
198  double median = getMedian(histo).value();
199  double x_lastBin = histo->GetBinLowEdge(nbins + 1);
200  const char* HistoName = histo->GetName();
201  std::string Finalname = "resMed_";
202  Finalname.append(HistoName);
203  TH1F* newHisto = new TH1F(Finalname.c_str(), Finalname.c_str(), nbins, 0., x_lastBin);
204  double* residuals = new double[nbins];
205  const float* weights = histo->GetArray();
206 
207  for (int j = 0; j < nbins; j++) {
208  residuals[j] = std::abs(median - histo->GetBinCenter(j + 1));
209  newHisto->Fill(residuals[j], weights[j]);
210  }
211 
212  double theMAD = (PVValHelper::getMedian(newHisto).value()) * 1.4826;
213 
214  delete[] residuals;
215  residuals = nullptr;
216  newHisto->Delete("");
217 
218  Measurement1D result(theMAD, theMAD / histo->GetEntries());
219  return result;
220 }
221 
222 //*************************************************************
223 std::pair<Measurement1D, Measurement1D> PVValHelper::fitResiduals(TH1* hist)
224 //*************************************************************
225 {
226  //float fitResult(9999);
227  if (hist->GetEntries() < 1) {
228  return std::make_pair(Measurement1D(0., 0.), Measurement1D(0., 0.));
229  };
230 
231  float mean = hist->GetMean();
232  float sigma = hist->GetRMS();
233 
234  TF1 func("tmp", "gaus", mean - 1.5 * sigma, mean + 1.5 * sigma);
235  if (0 == hist->Fit(&func, "QNR")) { // N: do not blow up file by storing fit!
236  mean = func.GetParameter(1);
237  sigma = func.GetParameter(2);
238  // second fit: three sigma of first fit around mean of first fit
239  func.SetRange(mean - 2 * sigma, mean + 2 * sigma);
240  // I: integral gives more correct results if binning is too wide
241  // L: Likelihood can treat empty bins correctly (if hist not weighted...)
242  if (0 == hist->Fit(&func, "Q0LR")) {
243  if (hist->GetFunction(func.GetName())) { // Take care that it is later on drawn:
244  hist->GetFunction(func.GetName())->ResetBit(TF1::kNotDraw);
245  }
246  }
247  }
248 
249  float res_mean = func.GetParameter(1);
250  float res_width = func.GetParameter(2);
251 
252  float res_mean_err = func.GetParError(1);
253  float res_width_err = func.GetParError(2);
254 
255  Measurement1D resultM(res_mean, res_mean_err);
256  Measurement1D resultW(res_width, res_width_err);
257 
258  std::pair<Measurement1D, Measurement1D> result;
259 
260  result = std::make_pair(resultM, resultW);
261  return result;
262 }
PVValHelper::getVarString
plotLabels getVarString(plotVariable var)
Definition: PVValidationHelpers.cc:130
PVValHelper::resz
Definition: PVValidationHelpers.h:53
FastTimerService_cff.range
range
Definition: FastTimerService_cff.py:34
Measurement1D
Definition: Measurement1D.h:11
start
Definition: start.py:1
MessageLogger.h
SiStripPI::mean
Definition: SiStripPayloadInspectorHelper.h:169
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
PVValHelper::norm_IP2D
Definition: PVValidationHelpers.h:61
Measurement1D::value
double value() const
Definition: Measurement1D.h:25
pfDeepBoostedJetPreprocessParams_cfi.median
median
Definition: pfDeepBoostedJetPreprocessParams_cfi.py:12
PVValHelper::fill
void fill(std::map< std::string, TH1 * > &h, const std::string &s, double x)
Definition: PVValidationHelpers.cc:20
PVValHelper::getTypeString
plotLabels getTypeString(residualType type)
Definition: PVValidationHelpers.cc:63
PVValHelper::fillByIndex
void fillByIndex(std::vector< TH1F * > &h, unsigned int index, double x, std::string tag="")
Definition: PVValidationHelpers.cc:42
timingPdfMaker.histo
histo
Definition: timingPdfMaker.py:278
cms::cuda::assert
assert(be >=bs)
PVValHelper::plotVariable
plotVariable
Definition: PVValidationHelpers.h:68
findQualityFiles.v
v
Definition: findQualityFiles.py:179
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
trigObjTnPSource_cfi.var
var
Definition: trigObjTnPSource_cfi.py:21
PVValHelper::norm_resz
Definition: PVValidationHelpers.h:62
alignCSCRings.s
s
Definition: alignCSCRings.py:92
PVValHelper::norm_d3D
Definition: PVValidationHelpers.h:64
PVValHelper::eta
Definition: PVValidationHelpers.h:70
PVValHelper::pT
Definition: PVValidationHelpers.h:71
HLT_FULL_cff.weights
weights
Definition: HLT_FULL_cff.py:99207
PVValHelper::pTCentral
Definition: PVValidationHelpers.h:72
mps_fire.end
end
Definition: mps_fire.py:242
LaserClient_cfi.nbins
nbins
Definition: LaserClient_cfi.py:51
PVValHelper::modZ
Definition: PVValidationHelpers.h:74
PVValHelper::norm_dxy
Definition: PVValidationHelpers.h:57
h
PVValHelper::norm_dx
Definition: PVValidationHelpers.h:58
reco::returnType
edm::TypeWithDict returnType(const edm::FunctionWithDict &func)
Definition: returnType.cc:16
PVValHelper::getMedian
Measurement1D getMedian(TH1F *histo)
Definition: PVValidationHelpers.cc:178
a
double a
Definition: hdecay.h:119
PVValHelper::IP3D
Definition: PVValidationHelpers.h:54
PVValHelper::phi
Definition: PVValidationHelpers.h:69
makeGlobalPositionRcd_cfg.tag
tag
Definition: makeGlobalPositionRcd_cfg.py:6
PVValidationHelpers.h
type
type
Definition: SiPixelVCal_PayloadInspector.cc:39
gpuVertexFinder::hist
__shared__ Hist hist
Definition: gpuClusterTracksDBSCAN.h:48
PVValHelper::add
void add(std::map< std::string, TH1 * > &h, TH1 *hist)
Definition: PVValidationHelpers.cc:12
PVValHelper::norm_IP3D
Definition: PVValidationHelpers.h:63
PVValHelper::residualType
residualType
Definition: PVValidationHelpers.h:47
cscdqm::HistoName
std::string HistoName
Definition: CSCDQM_HistoDef.h:32
readEcalDQMStatus.interval
interval
Definition: readEcalDQMStatus.py:18
PVValHelper::dy
Definition: PVValidationHelpers.h:50
submitPVResolutionJobs.q
q
Definition: submitPVResolutionJobs.py:84
PVValHelper::norm_dy
Definition: PVValidationHelpers.h:59
PVValHelper::generateBins
std::vector< float > generateBins(int n, float start, float range)
Definition: PVValidationHelpers.cc:161
TrackCollections2monitor_cff.func
func
Definition: TrackCollections2monitor_cff.py:359
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PVValHelper::d3D
Definition: PVValidationHelpers.h:55
PVValHelper::dxy
Definition: PVValidationHelpers.h:48
PVValHelper::dz
Definition: PVValidationHelpers.h:51
genVertex_cff.x
x
Definition: genVertex_cff.py:13
PVValHelper::getMAD
Measurement1D getMAD(TH1F *histo)
Definition: PVValidationHelpers.cc:194
PVValHelper::ladder
Definition: PVValidationHelpers.h:73
detailsBasic3DVector::y
float float y
Definition: extBasic3DVector.h:14
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
mps_fire.result
result
Definition: mps_fire.py:311
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
PVValHelper::dx
Definition: PVValidationHelpers.h:49
PVValHelper::plotLabels
std::tuple< std::string, std::string, std::string > plotLabels
Definition: PVValidationHelpers.h:123
PVValHelper::IP2D
Definition: PVValidationHelpers.h:52
PVValHelper::shrinkHistVectorToFit
void shrinkHistVectorToFit(std::vector< TH1F * > &h, unsigned int desired_size)
Definition: PVValidationHelpers.cc:56
PVValHelper::norm_dz
Definition: PVValidationHelpers.h:60
PVValHelper::fitResiduals
std::pair< Measurement1D, Measurement1D > fitResiduals(TH1 *hist)
Definition: PVValidationHelpers.cc:223