CMS 3D CMS Logo

ResolutionAnalyzer.cc
Go to the documentation of this file.
1 #ifndef RESOLUTIONANALYZER_CC
2 #define RESOLUTIONANALYZER_CC
3 
4 #include "ResolutionAnalyzer.h"
7 
8 //
9 // constants, enums and typedefs
10 //
11 
12 //
13 // static data member definitions
14 //
15 
16 //
17 // constructors and destructor
18 //
20  : theMuonLabel_(iConfig.getParameter<edm::InputTag>("MuonLabel")),
21  theMuonType_(iConfig.getParameter<int>("MuonType")),
22  theRootFileName_(iConfig.getUntrackedParameter<std::string>("OutputFileName")),
23  theCovariancesRootFileName_(iConfig.getUntrackedParameter<std::string>("InputFileName")),
24  debug_(iConfig.getUntrackedParameter<bool>("Debug")),
25  resonance_(iConfig.getParameter<bool>("Resonance")),
26  readCovariances_(iConfig.getParameter<bool>("ReadCovariances")),
27  treeFileName_(iConfig.getParameter<std::string>("InputTreeName")),
28  maxEvents_(iConfig.getParameter<int>("MaxEvents")),
29  ptMax_(iConfig.getParameter<double>("PtMax")) {
30  //now do what ever initialization is needed
31 
32  // Initial parameters values
33  // -------------------------
34  int resolFitType = iConfig.getParameter<int>("ResolFitType");
35  MuScleFitUtils::ResolFitType = resolFitType;
36  // MuScleFitUtils::resolutionFunction = resolutionFunctionArray[resolFitType];
38  // MuScleFitUtils::resolutionFunctionForVec = resolutionFunctionArrayForVec[resolFitType];
40 
41  MuScleFitUtils::parResol = iConfig.getParameter<std::vector<double> >("parResol");
42 
43  MuScleFitUtils::resfind = iConfig.getParameter<std::vector<int> >("ResFind");
44 
45  outputFile_ = new TFile(theRootFileName_.c_str(), "RECREATE");
46  outputFile_->cd();
47  fillHistoMap();
48 
49  eventCounter_ = 0;
50 }
51 
53  outputFile_->cd();
54  writeHistoMap();
55  outputFile_->Close();
56  std::cout << "Total analyzed events = " << eventCounter_ << std::endl;
57 }
58 
59 //
60 // member functions
61 //
62 
63 // ------------ method called to for each event ------------
65  using namespace edm;
66 
67  std::cout << "starting" << std::endl;
68 
69  lorentzVector nullLorentzVector(0, 0, 0, 0);
70 
71  RootTreeHandler rootTreeHandler;
72  typedef std::vector<std::pair<lorentzVector, lorentzVector> > MuonPairVector;
73  MuonPairVector savedPairVector;
74  MuonPairVector genPairVector;
75 
76  std::vector<std::pair<unsigned int, unsigned long long> > evtRun;
77  rootTreeHandler.readTree(maxEvents_, treeFileName_, &savedPairVector, 0, &evtRun, &genPairVector);
78  MuonPairVector::iterator savedPair = savedPairVector.begin();
79  MuonPairVector::iterator genPair = genPairVector.begin();
80  std::cout << "Starting loop on " << savedPairVector.size() << " muons" << std::endl;
81  for (; savedPair != savedPairVector.end(); ++savedPair, ++genPair) {
82  ++eventCounter_;
83 
84  if ((eventCounter_ % 10000) == 0) {
85  std::cout << "event = " << eventCounter_ << std::endl;
86  }
87 
88  lorentzVector recMu1(savedPair->first);
89  lorentzVector recMu2(savedPair->second);
90 
91  if (resonance_) {
92  // Histograms with genParticles characteristics
93  // --------------------------------------------
94 
95  reco::Particle::LorentzVector genMother(genPair->first + genPair->second);
96 
97  mapHisto_["GenMother"]->Fill(genMother);
98  mapHisto_["DeltaGenMotherMuons"]->Fill(genPair->first, genPair->second);
99  mapHisto_["GenMotherMuons"]->Fill(genPair->first);
100  mapHisto_["GenMotherMuons"]->Fill(genPair->second);
101 
102  // Match the reco muons with the gen and sim tracks
103  // ------------------------------------------------
104  if (checkDeltaR(genPair->first, recMu1)) {
105  mapHisto_["PtResolutionGenVSMu"]->Fill(recMu1, (-genPair->first.Pt() + recMu1.Pt()) / genPair->first.Pt(), -1);
106  mapHisto_["ThetaResolutionGenVSMu"]->Fill(recMu1, (-genPair->first.Theta() + recMu1.Theta()), -1);
107  mapHisto_["CotgThetaResolutionGenVSMu"]->Fill(
108  recMu1,
109  (-cos(genPair->first.Theta()) / sin(genPair->first.Theta()) + cos(recMu1.Theta()) / sin(recMu1.Theta())),
110  -1);
111  mapHisto_["EtaResolutionGenVSMu"]->Fill(recMu1, (-genPair->first.Eta() + recMu1.Eta()), -1);
112  // mapHisto_["PhiResolutionGenVSMu"]->Fill(recMu1,(-genPair->first.Phi()+recMu1.Phi()),-1);
113  mapHisto_["PhiResolutionGenVSMu"]->Fill(
114  recMu1, MuScleFitUtils::deltaPhiNoFabs(recMu1.Phi(), genPair->first.Phi()), -1);
115  recoPtVsgenPt_->Fill(genPair->first.Pt(), recMu1.Pt());
116  deltaPtOverPt_->Fill((recMu1.Pt() - genPair->first.Pt()) / genPair->first.Pt());
117  if (fabs(recMu1.Eta()) > 1 && fabs(recMu1.Eta()) < 1.2) {
118  recoPtVsgenPtEta12_->Fill(genPair->first.Pt(), recMu1.Pt());
119  deltaPtOverPtForEta12_->Fill((recMu1.Pt() - genPair->first.Pt()) / genPair->first.Pt());
120  }
121  }
122  if (checkDeltaR(genPair->second, recMu2)) {
123  mapHisto_["PtResolutionGenVSMu"]->Fill(
124  recMu2, (-genPair->second.Pt() + recMu2.Pt()) / genPair->second.Pt(), +1);
125  mapHisto_["ThetaResolutionGenVSMu"]->Fill(recMu2, (-genPair->second.Theta() + recMu2.Theta()), +1);
126  mapHisto_["CotgThetaResolutionGenVSMu"]->Fill(
127  recMu2,
128  (-cos(genPair->second.Theta()) / sin(genPair->second.Theta()) + cos(recMu2.Theta()) / sin(recMu2.Theta())),
129  +1);
130  mapHisto_["EtaResolutionGenVSMu"]->Fill(recMu2, (-genPair->second.Eta() + recMu2.Eta()), +1);
131  // mapHisto_["PhiResolutionGenVSMu"]->Fill(recMu2,(-genPair->second.Phi()+recMu2.Phi()),+1);
132  mapHisto_["PhiResolutionGenVSMu"]->Fill(
133  recMu2, MuScleFitUtils::deltaPhiNoFabs(recMu2.Phi(), genPair->second.Phi()), +1);
134  recoPtVsgenPt_->Fill(genPair->second.Pt(), recMu2.Pt());
135  deltaPtOverPt_->Fill((recMu2.Pt() - genPair->second.Pt()) / genPair->second.Pt());
136  if (fabs(recMu2.Eta()) > 1 && fabs(recMu2.Eta()) < 1.2) {
137  recoPtVsgenPtEta12_->Fill(genPair->second.Pt(), recMu2.Pt());
138  deltaPtOverPtForEta12_->Fill((recMu2.Pt() - genPair->second.Pt()) / genPair->second.Pt());
139  }
140  }
141 
142  // Fill the mass resolution histograms
143  // -----------------------------------
144  // check if the recoMuons match the genMuons
145  // if( MuScleFitUtils::ResFound && checkDeltaR(simMu.first,recMu1) && checkDeltaR(simMu.second,recMu2) ) {
146  if (genPair->first != nullLorentzVector && genPair->second != nullLorentzVector &&
147  checkDeltaR(genPair->first, recMu1) && checkDeltaR(genPair->second, recMu2)) {
148  double recoMass = (recMu1 + recMu2).mass();
149  double genMass = (genPair->first + genPair->second).mass();
150  // first is always mu-, second is always mu+
151  mapHisto_["MassResolution"]->Fill(recMu1, -1, genPair->first, recMu2, +1, genPair->second, recoMass, genMass);
152 
153  // Fill the reconstructed resonance
154  reco::Particle::LorentzVector recoResonance(recMu1 + recMu2);
155  mapHisto_["RecoResonance"]->Fill(recoResonance);
156  mapHisto_["DeltaRecoResonanceMuons"]->Fill(recMu1, recMu2);
157  mapHisto_["RecoResonanceMuons"]->Fill(recMu1);
158  mapHisto_["RecoResonanceMuons"]->Fill(recMu2);
159 
160  // Fill the mass resolution (computed from MC), we use the covariance class to compute the variance
161  if (genMass != 0) {
162  // double diffMass = (recoMass - genMass)/genMass;
163  double diffMass = recoMass - genMass;
164  // Fill if for both muons
165  double pt1 = recMu1.pt();
166  double eta1 = recMu1.eta();
167  double pt2 = recMu2.pt();
168  double eta2 = recMu2.eta();
169  // This is to avoid nan
170  if (diffMass == diffMass) {
171  massResolutionVsPtEta_->Fill(pt1, eta1, diffMass, diffMass);
172  massResolutionVsPtEta_->Fill(pt2, eta2, diffMass, diffMass);
173  } else {
174  std::cout << "Error, there is a nan: recoMass = " << recoMass << ", genMass = " << genMass << std::endl;
175  }
176  // Fill with mass resolution from resolution function
177  double massRes = MuScleFitUtils::massResolution(recMu1, recMu2, MuScleFitUtils::parResol);
178  // The value given by massRes is already divided by the mass, since the derivative functions have mass at the denominator.
179  // We alos take the squared value, since var = sigma^2.
180  mapHisto_["hFunctionResolMass"]->Fill(recMu1, std::pow(massRes, 2), -1);
181  mapHisto_["hFunctionResolMass"]->Fill(recMu2, std::pow(massRes, 2), +1);
182  }
183 
184  // Fill resolution functions for the muons (fill the squared value to make it comparable with the variance)
185  mapHisto_["hFunctionResolPt"]->Fill(
186  recMu1,
187  MuScleFitUtils::resolutionFunctionForVec->sigmaPt(recMu1.Pt(), recMu1.Eta(), MuScleFitUtils::parResol),
188  -1);
189  mapHisto_["hFunctionResolCotgTheta"]->Fill(
190  recMu1,
191  MuScleFitUtils::resolutionFunctionForVec->sigmaCotgTh(recMu1.Pt(), recMu1.Eta(), MuScleFitUtils::parResol),
192  -1);
193  mapHisto_["hFunctionResolPhi"]->Fill(
194  recMu1,
195  MuScleFitUtils::resolutionFunctionForVec->sigmaPhi(recMu1.Pt(), recMu1.Eta(), MuScleFitUtils::parResol),
196  -1);
197  mapHisto_["hFunctionResolPt"]->Fill(
198  recMu2,
199  MuScleFitUtils::resolutionFunctionForVec->sigmaPt(recMu2.Pt(), recMu2.Eta(), MuScleFitUtils::parResol),
200  +1);
201  mapHisto_["hFunctionResolCotgTheta"]->Fill(
202  recMu2,
203  MuScleFitUtils::resolutionFunctionForVec->sigmaCotgTh(recMu2.Pt(), recMu2.Eta(), MuScleFitUtils::parResol),
204  +1);
205  mapHisto_["hFunctionResolPhi"]->Fill(
206  recMu2,
207  MuScleFitUtils::resolutionFunctionForVec->sigmaPhi(recMu2.Pt(), recMu2.Eta(), MuScleFitUtils::parResol),
208  +1);
209 
210  if (readCovariances_) {
211  // Compute mass error terms
212  // ------------------------
213  double mass = (recMu1 + recMu2).mass();
214  double pt1 = recMu1.Pt();
215  double phi1 = recMu1.Phi();
216  double eta1 = recMu1.Eta();
217  double theta1 = 2 * atan(exp(-eta1));
218  double pt2 = recMu2.Pt();
219  double phi2 = recMu2.Phi();
220  double eta2 = recMu2.Eta();
221  double theta2 = 2 * atan(exp(-eta2));
222  // Derivatives
223  double mMu2 = MuScleFitUtils::mMu2;
224  double dmdpt1 = (pt1 / std::pow(sin(theta1), 2) *
225  sqrt((std::pow(pt2 / sin(theta2), 2) + mMu2) / (std::pow(pt1 / sin(theta1), 2) + mMu2)) -
226  pt2 * (cos(phi1 - phi2) + cos(theta1) * cos(theta2) / (sin(theta1) * sin(theta2)))) /
227  mass;
228  double dmdpt2 = (pt2 / std::pow(sin(theta2), 2) *
229  sqrt((std::pow(pt1 / sin(theta1), 2) + mMu2) / (std::pow(pt2 / sin(theta2), 2) + mMu2)) -
230  pt1 * (cos(phi2 - phi1) + cos(theta2) * cos(theta1) / (sin(theta2) * sin(theta1)))) /
231  mass;
232  double dmdphi1 = pt1 * pt2 / mass * sin(phi1 - phi2);
233  double dmdphi2 = pt2 * pt1 / mass * sin(phi2 - phi1);
234  double dmdcotgth1 =
235  (pt1 * pt1 * cos(theta1) / sin(theta1) *
236  sqrt((std::pow(pt2 / sin(theta2), 2) + mMu2) / (std::pow(pt1 / sin(theta1), 2) + mMu2)) -
237  pt1 * pt2 * cos(theta2) / sin(theta2)) /
238  mass;
239  double dmdcotgth2 =
240  (pt2 * pt2 * cos(theta2) / sin(theta2) *
241  sqrt((std::pow(pt1 / sin(theta1), 2) + mMu2) / (std::pow(pt2 / sin(theta2), 2) + mMu2)) -
242  pt2 * pt1 * cos(theta1) / sin(theta1)) /
243  mass;
244 
245  // Multiplied by the pt here
246  // -------------------------
247  double dmdpt[2] = {dmdpt1 * recMu1.Pt(), dmdpt2 * recMu2.Pt()};
248  double dmdphi[2] = {dmdphi1, dmdphi2};
249  double dmdcotgth[2] = {dmdcotgth1, dmdcotgth2};
250 
251  // Evaluate the single terms in the mass error expression
252 
253  reco::Particle::LorentzVector* recMu[2] = {&recMu1, &recMu2};
254  int charge[2] = {-1, +1};
255 
256  double fullMassRes = 0.;
257  double massRes1 = 0.;
258  double massRes2 = 0.;
259  double massRes3 = 0.;
260  double massRes4 = 0.;
261  double massRes5 = 0.;
262  double massRes6 = 0.;
263  double massResPtAndPt12 = 0.;
264 
265  for (int i = 0; i < 2; ++i) {
266  double ptVariance = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "Pt");
267  double cotgThetaVariance = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "CotgTheta");
268  double phiVariance = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "Phi");
269  double pt_cotgTheta = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "Pt-CotgTheta");
270  double pt_phi = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "Pt-Phi");
271  double cotgTheta_phi = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "CotgTheta-Phi");
272 
273  double pt1_pt2 = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "Pt1-Pt2");
274  double cotgTheta1_cotgTheta2 = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "CotgTheta1-CotgTheta2");
275  double phi1_phi2 = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "Phi1-Phi2");
276  double pt12_cotgTheta21 = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "Pt12-CotgTheta21");
277  double pt12_phi21 = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "Pt12-Phi21");
278  double cotgTheta12_phi21 = mapHisto_["ReadCovariances"]->Get(*(recMu[i]), "CotgTheta12-Phi21");
279 
280  // ATTENTION: Pt covariance terms are multiplied by Pt, since DeltaPt/Pt was used to compute them
281  mapHisto_["MassResolutionPt"]->Fill(*(recMu[i]), ptVariance * std::pow(dmdpt[i], 2), charge[i]);
282  mapHisto_["MassResolutionCotgTheta"]->Fill(
283  *(recMu[i]), cotgThetaVariance * std::pow(dmdcotgth[i], 2), charge[i]);
284  mapHisto_["MassResolutionPhi"]->Fill(*(recMu[i]), phiVariance * std::pow(dmdphi[i], 2), charge[i]);
285  mapHisto_["MassResolutionPt-CotgTheta"]->Fill(
286  *(recMu[i]), 2 * pt_cotgTheta * dmdpt[i] * dmdcotgth[i], charge[i]);
287  mapHisto_["MassResolutionPt-Phi"]->Fill(*(recMu[i]), 2 * pt_phi * dmdpt[i] * dmdphi[i], charge[i]);
288  mapHisto_["MassResolutionCotgTheta-Phi"]->Fill(
289  *(recMu[i]), 2 * cotgTheta_phi * dmdcotgth[i] * dmdphi[i], charge[i]);
290 
291  mapHisto_["MassResolutionPt1-Pt2"]->Fill(*(recMu[i]), pt1_pt2 * dmdpt[0] * dmdpt[1], charge[i]);
292  mapHisto_["MassResolutionCotgTheta1-CotgTheta2"]->Fill(
293  *(recMu[i]), cotgTheta1_cotgTheta2 * dmdcotgth[0] * dmdcotgth[1], charge[i]);
294  mapHisto_["MassResolutionPhi1-Phi2"]->Fill(*(recMu[i]), phi1_phi2 * dmdphi[0] * dmdphi[1], charge[i]);
295  // This must be filled for both configurations: 12 and 21
296  mapHisto_["MassResolutionPt12-CotgTheta21"]->Fill(
297  *(recMu[i]), pt12_cotgTheta21 * dmdpt[0] * dmdcotgth[1], charge[i]);
298  mapHisto_["MassResolutionPt12-CotgTheta21"]->Fill(
299  *(recMu[i]), pt12_cotgTheta21 * dmdpt[1] * dmdcotgth[0], charge[i]);
300  mapHisto_["MassResolutionPt12-Phi21"]->Fill(*(recMu[i]), pt12_phi21 * dmdpt[0] * dmdphi[1], charge[i]);
301  mapHisto_["MassResolutionPt12-Phi21"]->Fill(*(recMu[i]), pt12_phi21 * dmdpt[1] * dmdphi[0], charge[i]);
302  mapHisto_["MassResolutionCotgTheta12-Phi21"]->Fill(
303  *(recMu[i]), cotgTheta12_phi21 * dmdcotgth[0] * dmdphi[1], charge[i]);
304  mapHisto_["MassResolutionCotgTheta12-Phi21"]->Fill(
305  *(recMu[i]), cotgTheta12_phi21 * dmdcotgth[1] * dmdphi[0], charge[i]);
306 
307  // Sigmas for comparison
308  mapHisto_["sigmaPtFromVariance"]->Fill(*(recMu[i]), sqrt(ptVariance), charge[i]);
309  mapHisto_["sigmaCotgThetaFromVariance"]->Fill(*(recMu[i]), sqrt(cotgThetaVariance), charge[i]);
310  mapHisto_["sigmaPhiFromVariance"]->Fill(*(recMu[i]), sqrt(phiVariance), charge[i]);
311 
312  // Pt term from function
313  mapHisto_["MassResolutionPtFromFunction"]->Fill(
314  *(recMu[i]),
316  (recMu[i])->Pt(), (recMu[i])->Eta(), MuScleFitUtils::parResol)) *
317  std::pow(dmdpt[i], 2),
318  charge[i]);
319 
320  fullMassRes += ptVariance * std::pow(dmdpt[i], 2) + cotgThetaVariance * std::pow(dmdcotgth[i], 2) +
321  phiVariance * std::pow(dmdphi[i], 2) +
322 
323  // These are worth twice the others since there are: pt1-phi1, phi1-pt1, pt2-phi2, phi2-pt2
324  2 * pt_cotgTheta * dmdpt[i] * dmdcotgth[i] + 2 * pt_phi * dmdpt[i] * dmdphi[i] +
325  2 * cotgTheta_phi * dmdcotgth[i] * dmdphi[i] +
326 
327  pt1_pt2 * dmdpt[0] * dmdpt[1] + cotgTheta1_cotgTheta2 * dmdcotgth[0] * dmdcotgth[1] +
328  phi1_phi2 * dmdphi[0] * dmdphi[1] +
329 
330  // These are filled twice, because of the two combinations
331  pt12_cotgTheta21 * dmdpt[0] * dmdcotgth[1] + pt12_cotgTheta21 * dmdpt[1] * dmdcotgth[0] +
332  pt12_phi21 * dmdpt[0] * dmdphi[1] + pt12_phi21 * dmdpt[1] * dmdphi[0] +
333  cotgTheta12_phi21 * dmdcotgth[0] * dmdphi[1] + cotgTheta12_phi21 * dmdcotgth[1] * dmdphi[0];
334 
335  massRes1 += ptVariance * std::pow(dmdpt[i], 2);
336  massRes2 += ptVariance * std::pow(dmdpt[i], 2) + cotgThetaVariance * std::pow(dmdcotgth[i], 2);
337  massRes3 += ptVariance * std::pow(dmdpt[i], 2) + cotgThetaVariance * std::pow(dmdcotgth[i], 2) +
338  phiVariance * std::pow(dmdphi[i], 2);
339  massRes4 += ptVariance * std::pow(dmdpt[i], 2) + cotgThetaVariance * std::pow(dmdcotgth[i], 2) +
340  phiVariance * std::pow(dmdphi[i], 2) + pt1_pt2 * dmdpt[0] * dmdpt[1] +
341  2 * pt_cotgTheta * dmdpt[i] * dmdcotgth[i] + 2 * pt_phi * dmdpt[i] * dmdphi[i] +
342  2 * cotgTheta_phi * dmdcotgth[i] * dmdphi[i];
343  massRes5 += ptVariance * std::pow(dmdpt[i], 2) + cotgThetaVariance * std::pow(dmdcotgth[i], 2) +
344  phiVariance * std::pow(dmdphi[i], 2) + pt1_pt2 * dmdpt[0] * dmdpt[1] +
345  2 * pt_cotgTheta * dmdpt[i] * dmdcotgth[i] + 2 * pt_phi * dmdpt[i] * dmdphi[i] +
346  2 * cotgTheta_phi * dmdcotgth[i] * dmdphi[i] +
347  cotgTheta1_cotgTheta2 * dmdcotgth[0] * dmdcotgth[1] + phi1_phi2 * dmdphi[0] * dmdphi[1];
348  massRes6 += ptVariance * std::pow(dmdpt[i], 2) + cotgThetaVariance * std::pow(dmdcotgth[i], 2) +
349  phiVariance * std::pow(dmdphi[i], 2) + pt1_pt2 * dmdpt[0] * dmdpt[1] +
350  2 * pt_cotgTheta * dmdpt[i] * dmdcotgth[i] + 2 * pt_phi * dmdpt[i] * dmdphi[i] +
351  2 * cotgTheta_phi * dmdcotgth[i] * dmdphi[i] +
352  cotgTheta1_cotgTheta2 * dmdcotgth[0] * dmdcotgth[1] + phi1_phi2 * dmdphi[0] * dmdphi[1] +
353  pt12_cotgTheta21 * dmdpt[0] * dmdcotgth[1] + pt12_cotgTheta21 * dmdpt[1] * dmdcotgth[0] +
354  pt12_phi21 * dmdpt[0] * dmdphi[1] + pt12_phi21 * dmdpt[1] * dmdphi[0] +
355  cotgTheta12_phi21 * dmdcotgth[0] * dmdphi[1] + cotgTheta12_phi21 * dmdcotgth[1] * dmdphi[0];
356 
357  massResPtAndPt12 += ptVariance * std::pow(dmdpt[i], 2) + pt1_pt2 * dmdpt[0] * dmdpt[1];
358 
359  // Derivatives
360  mapHisto_["DerivativePt"]->Fill(*(recMu[i]), dmdpt[i], charge[i]);
361  mapHisto_["DerivativeCotgTheta"]->Fill(*(recMu[i]), dmdcotgth[i], charge[i]);
362  mapHisto_["DerivativePhi"]->Fill(*(recMu[i]), dmdphi[i], charge[i]);
363  }
364  // Fill the complete resolution function with covariance terms
365  mapHisto_["FullMassResolution"]->Fill(*(recMu[0]), fullMassRes, charge[0]);
366  mapHisto_["FullMassResolution"]->Fill(*(recMu[1]), fullMassRes, charge[1]);
367 
368  mapHisto_["MassRes1"]->Fill(*(recMu[0]), massRes1, charge[0]);
369  mapHisto_["MassRes1"]->Fill(*(recMu[1]), massRes1, charge[1]);
370  mapHisto_["MassRes2"]->Fill(*(recMu[0]), massRes2, charge[0]);
371  mapHisto_["MassRes2"]->Fill(*(recMu[1]), massRes2, charge[1]);
372  mapHisto_["MassRes3"]->Fill(*(recMu[0]), massRes3, charge[0]);
373  mapHisto_["MassRes3"]->Fill(*(recMu[1]), massRes3, charge[1]);
374  mapHisto_["MassRes4"]->Fill(*(recMu[0]), massRes4, charge[0]);
375  mapHisto_["MassRes4"]->Fill(*(recMu[1]), massRes4, charge[1]);
376  mapHisto_["MassRes5"]->Fill(*(recMu[0]), massRes5, charge[0]);
377  mapHisto_["MassRes5"]->Fill(*(recMu[1]), massRes5, charge[1]);
378  mapHisto_["MassRes6"]->Fill(*(recMu[0]), massRes6, charge[0]);
379  mapHisto_["MassRes6"]->Fill(*(recMu[1]), massRes6, charge[1]);
380  mapHisto_["MassResPtAndPt12"]->Fill(*(recMu[0]), massResPtAndPt12, charge[0]);
381  mapHisto_["MassResPtAndPt12"]->Fill(*(recMu[1]), massResPtAndPt12, charge[1]);
382  } else {
383  // Fill the covariances histograms
384  mapHisto_["Covariances"]->Fill(recMu1, genPair->first, recMu2, genPair->second);
385  }
386  }
387  } // end if resonance
388  }
389  // else {
390  //
391  // // Loop on the recMuons
392  // std::vector<reco::LeafCandidate>::const_iterator recMuon = muons.begin();
393  // for ( ; recMuon!=muons.end(); ++recMuon ) {
394  // int charge = recMuon->charge();
395  //
396  // lorentzVector recMu(recMuon->p4());
397  //
398  // // Find the matching MC muon
399  // const HepMC::GenEvent* Evt = evtMC->GetEvent();
400  // //Loop on generated particles
401  // std::map<double, lorentzVector> genAssocMap;
402  // HepMC::GenEvent::particle_const_iterator part = Evt->particles_begin();
403  // for( ; part!=Evt->particles_end(); ++part ) {
404  // if (fabs((*part)->pdg_id())==13 && (*part)->status()==1) {
405  // lorentzVector genMu = (lorentzVector((*part)->momentum().px(),(*part)->momentum().py(),
406  // (*part)->momentum().pz(),(*part)->momentum().e()));
407  //
408  // double deltaR = sqrt(MuScleFitUtils::deltaPhi(recMu.Phi(),genPair->Phi()) * MuScleFitUtils::deltaPhi(recMu.Phi(),genPair->Phi()) +
409  // ((recMu.Eta()-genPair->Eta()) * (recMu.Eta()-genPair->Eta())));
410  //
411  // // 13 for the muon (-1) and -13 for the antimuon (+1), thus pdg*charge = -13.
412  // // Only in this case we consider it matching.
413  // if( ((*part)->pdg_id())*charge == -13 ) genAssocMap.insert(std::make_pair(deltaR, genMu));
414  // }
415  // }
416  // // Take the closest in deltaR
417  // lorentzVector genMu(genAssocMap.begin()->second);
418  //
419  // // Histograms with genParticles characteristics
420  // // --------------------------------------------
421  //
422  // if(checkDeltaR(genMu,recMu)){
423  // mapHisto_["PtResolutionGenVSMu"]->Fill(genMu,(-genPair->Pt()+recMu.Pt())/genPair->Pt(),charge);
424  // mapHisto_["ThetaResolutionGenVSMu"]->Fill(genMu,(-genPair->Theta()+recMu.Theta()),charge);
425  // mapHisto_["CotgThetaResolutionGenVSMu"]->Fill(genMu,(-cos(genPair->Theta())/sin(genPair->Theta())
426  // +cos(recMu.Theta())/sin(recMu.Theta())),charge);
427  // mapHisto_["EtaResolutionGenVSMu"]->Fill(genMu,(-genPair->Eta()+recMu.Eta()),charge);
428  // mapHisto_["PhiResolutionGenVSMu"]->Fill(genMu,MuScleFitUtils::deltaPhiNoFabs(recMu.Phi(), genPair->Phi()),charge);
429  // }
430  // }
431 }
432 
434  outputFile_->cd();
435 
436  // Resonances
437  // If no Z is required, use a smaller mass range.
438  double minMass = 0.;
439  double maxMass = 200.;
440  if (MuScleFitUtils::resfind[0] != 1) {
441  maxMass = 30.;
442  }
443  mapHisto_["GenMother"] = new HParticle(outputFile_, "GenMother", minMass, maxMass);
444  mapHisto_["SimResonance"] = new HParticle(outputFile_, "SimResonance", minMass, maxMass);
445  mapHisto_["RecoResonance"] = new HParticle(outputFile_, "RecoResonance", minMass, maxMass);
446 
447  // Resonance muons
448  mapHisto_["GenMotherMuons"] = new HParticle(outputFile_, "GenMotherMuons", minMass, 1.);
449  mapHisto_["SimResonanceMuons"] = new HParticle(outputFile_, "SimResonanceMuons", minMass, 1.);
450  mapHisto_["RecoResonanceMuons"] = new HParticle(outputFile_, "RecoResonanceMuons", minMass, 1.);
451 
452  // Deltas between resonance muons
453  mapHisto_["DeltaGenMotherMuons"] = new HDelta(outputFile_, "DeltaGenMotherMuons");
454  mapHisto_["DeltaSimResonanceMuons"] = new HDelta(outputFile_, "DeltaSimResonanceMuons");
455  mapHisto_["DeltaRecoResonanceMuons"] = new HDelta(outputFile_, "DeltaRecoResonanceMuons");
456 
457  // //Reconstructed muon kinematics
458  // //-----------------------------
459  // mapHisto_["hRecBestMu"] = new HParticle ("hRecBestMu");
460  // mapHisto_["hRecBestMu_Acc"] = new HParticle ("hRecBestMu_Acc");
461  // mapHisto_["hDeltaRecBestMu"] = new HDelta ("hDeltaRecBestMu");
462 
463  // mapHisto_["hRecBestRes"] = new HParticle ("hRecBestRes");
464  // mapHisto_["hRecBestRes_Acc"] = new HParticle ("hRecBestRes_Acc");
465  // mapHisto_["hRecBestResVSMu"] = new HMassVSPart ("hRecBestResVSMu");
466 
467  //Resolution VS muon kinematic
468  //----------------------------
469  mapHisto_["PtResolutionGenVSMu"] = new HResolutionVSPart(outputFile_, "PtResolutionGenVSMu");
470  mapHisto_["PtResolutionSimVSMu"] = new HResolutionVSPart(outputFile_, "PtResolutionSimVSMu");
471  mapHisto_["EtaResolutionGenVSMu"] = new HResolutionVSPart(outputFile_, "EtaResolutionGenVSMu");
472  mapHisto_["EtaResolutionSimVSMu"] = new HResolutionVSPart(outputFile_, "EtaResolutionSimVSMu");
473  mapHisto_["ThetaResolutionGenVSMu"] = new HResolutionVSPart(outputFile_, "ThetaResolutionGenVSMu");
474  mapHisto_["ThetaResolutionSimVSMu"] = new HResolutionVSPart(outputFile_, "ThetaResolutionSimVSMu");
475  mapHisto_["CotgThetaResolutionGenVSMu"] =
476  new HResolutionVSPart(outputFile_, "CotgThetaResolutionGenVSMu", -0.02, 0.02, -0.02, 0.02);
477  mapHisto_["CotgThetaResolutionSimVSMu"] = new HResolutionVSPart(outputFile_, "CotgThetaResolutionSimVSMu");
478  mapHisto_["PhiResolutionGenVSMu"] =
479  new HResolutionVSPart(outputFile_, "PhiResolutionGenVSMu", -0.002, 0.002, -0.002, 0.002);
480  mapHisto_["PhiResolutionSimVSMu"] = new HResolutionVSPart(outputFile_, "PhiResolutionSimVSMu");
481 
482  // Covariances between muons kinematic quantities
483  // ----------------------------------------------
484  double ptMax = ptMax_;
485 
486  // Mass resolution
487  // ---------------
488  mapHisto_["MassResolution"] = new HMassResolutionVSPart(outputFile_, "MassResolution");
489 
490  // mapHisto_["hResolRecoMassVSGenMassVSPt"] = new HResolutionVSPart
491 
492  // Mass resolution vs (pt, eta) of the muons from MC
493  massResolutionVsPtEta_ = new HCovarianceVSxy("Mass", "Mass", 100, 0., ptMax, 60, -3, 3);
494  // Mass resolution vs (pt, eta) of the muons from function
495  recoPtVsgenPt_ = new TH2D("recoPtVsgenPt", "recoPtVsgenPt", 100, 0, ptMax, 100, 0, ptMax);
496  recoPtVsgenPtEta12_ = new TH2D("recoPtVsgenPtEta12", "recoPtVsgenPtEta12", 100, 0, ptMax, 100, 0, ptMax);
497  deltaPtOverPt_ = new TH1D("deltaPtOverPt", "deltaPtOverPt", 100, -0.1, 0.1);
498  deltaPtOverPtForEta12_ = new TH1D("deltaPtOverPtForEta12", "deltaPtOverPtForEta12", 100, -0.1, 0.1);
499 
500  // Muons resolutions from resolution functions
501  // -------------------------------------------
502  int totBinsY = 60;
503  mapHisto_["hFunctionResolMass"] = new HFunctionResolution(outputFile_, "hFunctionResolMass", ptMax, totBinsY);
504  mapHisto_["hFunctionResolPt"] = new HFunctionResolution(outputFile_, "hFunctionResolPt", ptMax, totBinsY);
505  mapHisto_["hFunctionResolCotgTheta"] =
506  new HFunctionResolution(outputFile_, "hFunctionResolCotgTheta", ptMax, totBinsY);
507  mapHisto_["hFunctionResolPhi"] = new HFunctionResolution(outputFile_, "hFunctionResolPhi", ptMax, totBinsY);
508 
509  if (readCovariances_) {
510  // Covariances read from file. Used to compare the terms in the expression of mass error
511  mapHisto_["ReadCovariances"] = new HCovarianceVSParts(theCovariancesRootFileName_, "Covariance");
512 
513  // Variances
514  mapHisto_["MassResolutionPt"] = new HFunctionResolutionVarianceCheck(outputFile_, "functionResolMassPt", ptMax);
515  mapHisto_["MassResolutionCotgTheta"] =
516  new HFunctionResolutionVarianceCheck(outputFile_, "functionResolMassCotgTheta", ptMax);
517  mapHisto_["MassResolutionPhi"] = new HFunctionResolutionVarianceCheck(outputFile_, "functionResolMassPhi", ptMax);
518  // Covariances
519  mapHisto_["MassResolutionPt-CotgTheta"] =
520  new HFunctionResolution(outputFile_, "functionResolMassPt-CotgTheta", ptMax, totBinsY);
521  mapHisto_["MassResolutionPt-Phi"] =
522  new HFunctionResolution(outputFile_, "functionResolMassPt-Phi", ptMax, totBinsY);
523  mapHisto_["MassResolutionCotgTheta-Phi"] =
524  new HFunctionResolution(outputFile_, "functionResolMassCotgTheta-Phi", ptMax, totBinsY);
525  mapHisto_["MassResolutionPt1-Pt2"] =
526  new HFunctionResolution(outputFile_, "functionResolMassPt1-Pt2", ptMax, totBinsY);
527  mapHisto_["MassResolutionCotgTheta1-CotgTheta2"] =
528  new HFunctionResolution(outputFile_, "functionResolMassCotgTheta1-CotgTheta2", ptMax, totBinsY);
529  mapHisto_["MassResolutionPhi1-Phi2"] =
530  new HFunctionResolution(outputFile_, "functionResolMassPhi1-Phi2", ptMax, totBinsY);
531  mapHisto_["MassResolutionPt12-CotgTheta21"] =
532  new HFunctionResolution(outputFile_, "functionResolMassPt12-CotgTheta21", ptMax, totBinsY);
533  mapHisto_["MassResolutionPt12-Phi21"] =
534  new HFunctionResolution(outputFile_, "functionResolMassPt12-Phi21", ptMax, totBinsY);
535  mapHisto_["MassResolutionCotgTheta12-Phi21"] =
536  new HFunctionResolution(outputFile_, "functionResolMassCotgTheta12-Phi21", ptMax, totBinsY);
537 
538  mapHisto_["sigmaPtFromVariance"] = new HFunctionResolution(outputFile_, "sigmaPtFromVariance", ptMax, totBinsY);
539  mapHisto_["sigmaCotgThetaFromVariance"] =
540  new HFunctionResolution(outputFile_, "sigmaCotgThetaFromVariance", ptMax, totBinsY);
541  mapHisto_["sigmaPhiFromVariance"] = new HFunctionResolution(outputFile_, "sigmaPhiFromVariance", ptMax, totBinsY);
542 
543  // Derivatives
544  mapHisto_["DerivativePt"] = new HFunctionResolution(outputFile_, "derivativePt", ptMax);
545  mapHisto_["DerivativeCotgTheta"] = new HFunctionResolution(outputFile_, "derivativeCotgTheta", ptMax);
546  mapHisto_["DerivativePhi"] = new HFunctionResolution(outputFile_, "derivativePhi", ptMax);
547 
548  // Pt term from function
549  mapHisto_["MassResolutionPtFromFunction"] =
550  new HFunctionResolutionVarianceCheck(outputFile_, "functionResolMassPtFromFunction", ptMax);
551 
552  mapHisto_["FullMassResolution"] = new HFunctionResolution(outputFile_, "fullMassResolution", ptMax);
553  mapHisto_["MassRes1"] = new HFunctionResolution(outputFile_, "massRes1", ptMax);
554  mapHisto_["MassRes2"] = new HFunctionResolution(outputFile_, "massRes2", ptMax);
555  mapHisto_["MassRes3"] = new HFunctionResolution(outputFile_, "massRes3", ptMax);
556  mapHisto_["MassRes4"] = new HFunctionResolution(outputFile_, "massRes4", ptMax);
557  mapHisto_["MassRes5"] = new HFunctionResolution(outputFile_, "massRes5", ptMax);
558  mapHisto_["MassRes6"] = new HFunctionResolution(outputFile_, "massRes6", ptMax);
559  mapHisto_["MassResPtAndPt12"] = new HFunctionResolution(outputFile_, "massResPtAndPt12", ptMax);
560  } else {
561  mapHisto_["Covariances"] = new HCovarianceVSParts(outputFile_, "Covariance", ptMax);
562  }
563 }
564 
566  for (std::map<std::string, Histograms*>::const_iterator histo = mapHisto_.begin(); histo != mapHisto_.end();
567  histo++) {
568  (*histo).second->Write();
569  }
570  outputFile_->cd();
572  recoPtVsgenPt_->Write();
573  recoPtVsgenPtEta12_->Write();
574  deltaPtOverPt_->Write();
575  deltaPtOverPtForEta12_->Write();
576 }
577 
579  const reco::Particle::LorentzVector& recMu) {
580  //first is always mu-, second is always mu+
581  double deltaR =
582  sqrt(MuScleFitUtils::deltaPhi(recMu.Phi(), genMu.Phi()) * MuScleFitUtils::deltaPhi(recMu.Phi(), genMu.Phi()) +
583  ((recMu.Eta() - genMu.Eta()) * (recMu.Eta() - genMu.Eta())));
584  if (deltaR < 0.01)
585  return true;
586  else if (debug_ > 0)
587  std::cout << "Reco muon " << recMu << " with eta " << recMu.Eta() << " and phi " << recMu.Phi() << std::endl
588  << " DOES NOT MATCH with generated muon from resonance: " << std::endl
589  << genMu << " with eta " << genMu.Eta() << " and phi " << genMu.Phi() << std::endl;
590  return false;
591 }
592 
593 //define this as a plug-in
595 
596 #endif // RESOLUTIONANALYZER_CC
ResolutionAnalyzer::theRootFileName_
std::string theRootFileName_
Definition: ResolutionAnalyzer.h:107
ResolutionAnalyzer::ptMax_
double ptMax_
Definition: ResolutionAnalyzer.h:120
MuScleFitUtils::deltaPhiNoFabs
static double deltaPhiNoFabs(const double &phi1, const double &phi2)
Without fabs at the end, used to have a symmetric distribution for the resolution fits and variance c...
Definition: MuScleFitUtils.h:118
electrons_cff.bool
bool
Definition: electrons_cff.py:393
mps_fire.i
i
Definition: mps_fire.py:428
RootTreeHandler
Definition: RootTreeHandler.h:24
ResolutionAnalyzer::deltaPtOverPtForEta12_
TH1D * deltaPtOverPtForEta12_
Definition: ResolutionAnalyzer.h:126
ResolutionAnalyzer::fillHistoMap
void fillHistoMap()
Used to fill the map with the histograms needed.
Definition: ResolutionAnalyzer.cc:433
ResolutionAnalyzer::~ResolutionAnalyzer
~ResolutionAnalyzer() override
Definition: ResolutionAnalyzer.cc:52
MuScleFitUtils::deltaPhi
static double deltaPhi(const double &phi1, const double &phi2)
Definition: MuScleFitUtils.h:109
ResolutionAnalyzer::checkDeltaR
bool checkDeltaR(const reco::Particle::LorentzVector &genMu, const reco::Particle::LorentzVector &recMu)
Returns true if the two particles have DeltaR < cut.
Definition: ResolutionAnalyzer.cc:578
HDelta
Definition: Histograms.h:423
edm
HLT enums.
Definition: AlignableModifier.h:19
ResolutionAnalyzer::mapHisto_
std::map< std::string, Histograms * > mapHisto_
Definition: ResolutionAnalyzer.h:110
gather_cfg.cout
cout
Definition: gather_cfg.py:144
MuScleFitUtils::mMu2
static const double mMu2
Definition: MuScleFitUtils.h:139
HFunctionResolution
Definition: Histograms.h:1548
ResolutionAnalyzer::massResolutionVsPtEta_
HCovarianceVSxy * massResolutionVsPtEta_
Definition: ResolutionAnalyzer.h:122
ResolutionAnalyzer::resonance_
bool resonance_
Definition: ResolutionAnalyzer.h:114
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:85964
ResolutionAnalyzer::maxEvents_
int32_t maxEvents_
Definition: ResolutionAnalyzer.h:118
timingPdfMaker.histo
histo
Definition: timingPdfMaker.py:279
resolutionFunctionService
resolutionFunctionBase< double * > * resolutionFunctionService(const int identifier)
Service to build the resolution functor corresponding to the passed identifier.
Definition: Functions.cc:70
ResolutionAnalyzer::debug_
bool debug_
Definition: ResolutionAnalyzer.h:109
MuScleFitUtils::resfind
static std::vector< int > resfind
Definition: MuScleFitUtils.h:210
reco::Particle::LorentzVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
HCovarianceVSxy
Definition: Histograms.h:1906
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
HCovarianceVSxy::Write
void Write() override
Definition: Histograms.h:2018
resolutionFunctionVecService
resolutionFunctionBase< std::vector< double > > * resolutionFunctionVecService(const int identifier)
Service to build the resolution functor corresponding to the passed identifier when receiving a std::...
Definition: Functions.cc:92
HMassResolutionVSPart
Definition: Histograms.h:2265
MuonPairVector
std::vector< std::pair< lorentzVector, lorentzVector > > MuonPairVector
Definition: RootTreeHandler.h:13
ResolutionAnalyzer::outputFile_
TFile * outputFile_
Definition: ResolutionAnalyzer.h:111
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MuScleFitUtils::parResol
static std::vector< double > parResol
Definition: MuScleFitUtils.h:192
AlignmentTrackSelector_cfi.ptMax
ptMax
Definition: AlignmentTrackSelector_cfi.py:12
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
HLT_FULL_cff.eta2
eta2
Definition: HLT_FULL_cff.py:9605
HResolutionVSPart
A set of histograms for resolution.
Definition: Histograms.h:1199
HLT_FULL_cff.pt1
pt1
Definition: HLT_FULL_cff.py:9933
RootTreeHandler.h
MuScleFitUtils::massResolution
static double massResolution(const lorentzVector &mu1, const lorentzVector &mu2)
MuScleFitUtils::ResolFitType
static int ResolFitType
Definition: MuScleFitUtils.h:152
MuScleFitUtils::resolutionFunctionForVec
static resolutionFunctionBase< std::vector< double > > * resolutionFunctionForVec
Definition: MuScleFitUtils.h:154
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HLT_FULL_cff.eta1
eta1
Definition: HLT_FULL_cff.py:9604
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
edm::ParameterSet
Definition: ParameterSet.h:47
lorentzVector
reco::Particle::LorentzVector lorentzVector
Definition: GenMuonPair.h:9
B2GTnPMonitor_cfi.maxMass
maxMass
Definition: B2GTnPMonitor_cfi.py:55
ResolutionAnalyzer::readCovariances_
bool readCovariances_
Definition: ResolutionAnalyzer.h:115
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
ResolutionAnalyzer::ResolutionAnalyzer
ResolutionAnalyzer(const edm::ParameterSet &)
Definition: ResolutionAnalyzer.cc:19
HCovarianceVSParts
Definition: Histograms.h:2086
ResolutionAnalyzer::theCovariancesRootFileName_
std::string theCovariancesRootFileName_
Definition: ResolutionAnalyzer.h:108
ResolutionAnalyzer
Definition: ResolutionAnalyzer.h:65
ResolutionAnalyzer::recoPtVsgenPt_
TH2D * recoPtVsgenPt_
Definition: ResolutionAnalyzer.h:123
edm::EventSetup
Definition: EventSetup.h:57
HLT_FULL_cff.pt2
pt2
Definition: HLT_FULL_cff.py:9935
ResolutionAnalyzer::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: ResolutionAnalyzer.cc:64
ResolutionAnalyzer::deltaPtOverPt_
TH1D * deltaPtOverPt_
Definition: ResolutionAnalyzer.h:125
HFunctionResolutionVarianceCheck
Definition: Histograms.h:1739
ResolutionAnalyzer::recoPtVsgenPtEta12_
TH2D * recoPtVsgenPtEta12_
Definition: ResolutionAnalyzer.h:124
ResolutionAnalyzer.h
std
Definition: JetResolutionObject.h:76
HCovarianceVSxy::Fill
void Fill(const double &x, const double &y, const double &a, const double &b) override
Definition: Histograms.h:1987
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
B2GTnPMonitor_cfi.minMass
minMass
Definition: B2GTnPMonitor_cfi.py:54
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ResolutionAnalyzer::eventCounter_
int eventCounter_
Definition: ResolutionAnalyzer.h:113
ResolutionAnalyzer::treeFileName_
TString treeFileName_
Definition: ResolutionAnalyzer.h:117
RootTreeHandler::readTree
void readTree(const int maxEvents, const TString &fileName, MuonPairVector *savedPair, const int muonType, std::vector< std::pair< unsigned int, unsigned long long > > *evtRun, MuonPairVector *genPair=nullptr)
Definition: RootTreeHandler.h:86
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
Functions.h
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
edm::Event
Definition: Event.h:73
MuScleFitUtils::resolutionFunction
static resolutionFunctionBase< double * > * resolutionFunction
Definition: MuScleFitUtils.h:153
HParticle
Definition: Histograms.h:242
ResolutionAnalyzer::writeHistoMap
void writeHistoMap()
Writes the histograms in the map.
Definition: ResolutionAnalyzer.cc:565