CMS 3D CMS Logo

V0Validator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: V0Validator
4 // Class: V0Validator
5 //
13 //
14 // Original Author: Brian Drell
15 // Created: Wed Feb 18 17:21:04 MST 2009
16 //
17 //
18 
21 
22 typedef std::vector<TrackingVertex> TrackingVertexCollection;
26 
28  : theDQMRootFileName(iConfig.getUntrackedParameter<std::string>("DQMRootFileName")),
29  dirName(iConfig.getUntrackedParameter<std::string>("dirName")),
30  recoRecoToSimCollectionToken_(
31  consumes<reco::RecoToSimCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackAssociatorMap"))),
32  recoSimToRecoCollectionToken_(
33  consumes<reco::SimToRecoCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackAssociatorMap"))),
34  trackingVertexCollection_Token_(
35  consumes<TrackingVertexCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackingVertexCollection"))),
36  vec_recoVertex_Token_(
37  consumes<std::vector<reco::Vertex> >(iConfig.getUntrackedParameter<edm::InputTag>("vertexCollection"))),
38  recoVertexCompositeCandidateCollection_k0s_Token_(consumes<reco::VertexCompositeCandidateCollection>(
39  iConfig.getUntrackedParameter<edm::InputTag>("kShortCollection"))),
40  recoVertexCompositeCandidateCollection_lambda_Token_(consumes<reco::VertexCompositeCandidateCollection>(
41  iConfig.getUntrackedParameter<edm::InputTag>("lambdaCollection"))) {}
42 
44 
46  double minKsMass = 0.49767 - 0.07;
47  double maxKsMass = 0.49767 + 0.07;
48  double minLamMass = 1.1156 - 0.05;
49  double maxLamMass = 1.1156 + 0.05;
50  int ksMassNbins = 100;
51  double ksMassXmin = minKsMass;
52  double ksMassXmax = maxKsMass;
53  int lamMassNbins = 100;
54  double lamMassXmin = minLamMass;
55  double lamMassXmax = maxLamMass;
56 
57  ibooker.cd();
58  std::string subDirName = V0Validator::dirName + "/K0";
59  ibooker.setCurrentFolder(subDirName);
60 
62  ibooker.book1D("K0sEffVsR_num", "K^{0}_{S} Efficiency vs #rho", 80, 0., 40.);
64  ibooker.book1D("K0sEffVsEta_num", "K^{0}_{S} Efficiency vs #eta", 40, -2.5, 2.5);
66  ibooker.book1D("K0sEffVsPt_num", "K^{0}_{S} Efficiency vs p_{T}", 70, 0., 20.);
67 
69  ibooker.book1D("K0sTkEffVsR_num", "K^{0}_{S} Tracking Efficiency vs #rho", 80, 0., 40.);
71  ibooker.book1D("K0sTkEffVsEta_num", "K^{0}_{S} Tracking Efficiency vs #eta", 40, -2.5, 2.5);
73  ibooker.book1D("K0sTkEffVsPt_num", "K^{0}_{S} Tracking Efficiency vs p_{T}", 70, 0., 20.);
74 
76  ibooker.book1D("K0sEffVsR_denom", "K^{0}_{S} Efficiency vs #rho", 80, 0., 40.);
78  ibooker.book1D("K0sEffVsEta_denom", "K^{0}_{S} Efficiency vs #eta", 40, -2.5, 2.5);
80  ibooker.book1D("K0sEffVsPt_denom", "K^{0}_{S} Efficiency vs p_{T}", 70, 0., 20.);
81 
83  ibooker.book1D("K0sFakeVsR_num", "K^{0}_{S} Fake Rate vs #rho", 80, 0., 40.);
85  ibooker.book1D("K0sFakeVsEta_num", "K^{0}_{S} Fake Rate vs #eta", 40, -2.5, 2.5);
87  ibooker.book1D("K0sFakeVsPt_num", "K^{0}_{S} Fake Rate vs p_{T}", 70, 0., 20.);
89  ibooker.book1D("K0sTkFakeVsR_num", "K^{0}_{S} Tracking Fake Rate vs #rho", 80, 0., 80.);
91  ibooker.book1D("K0sTkFakeVsEta_num", "K^{0}_{S} Tracking Fake Rate vs #eta", 40, -2.5, 2.5);
93  ibooker.book1D("K0sTkFakeVsPt_num", "K^{0}_{S} Tracking Fake Rate vs p_{T}", 70, 0., 20.);
94 
96  ibooker.book1D("K0sFakeVsR_denom", "K^{0}_{S} Fake Rate vs #rho", 80, 0., 40.);
98  ibooker.book1D("K0sFakeVsEta_denom", "K^{0}_{S} Fake Rate vs #eta", 40, -2.5, 2.5);
100  ibooker.book1D("K0sFakeVsPt_denom", "K^{0}_{S} Fake Rate vs p_{T}", 70, 0., 20.);
101  nCandidates_[V0Validator::KSHORT] = ibooker.book1D("nK0s", "Number of K^{0}_{S} found per event", 60, 0., 60.);
103  ibooker.book1D("ksMassFake", "Mass of fake K0S", ksMassNbins, minKsMass, maxKsMass);
105  ibooker.book1D("ksMassGood", "Mass of good reco K0S", ksMassNbins, minKsMass, maxKsMass);
107  ibooker.book1D("ksMassAll", "Invariant mass of all K0S", ksMassNbins, ksMassXmin, ksMassXmax);
109  ibooker.book1D("radDistFakeKs", "Production radius of daughter particle of Ks fake", 100, 0., 15.);
110  candidateStatus_[V0Validator::KSHORT] = ibooker.book1D("ksCandStatus", "Fake type by cand status", 10, 0., 10.);
111 
112  // Lambda Plots follow
113 
114  subDirName = V0Validator::dirName + "/Lambda";
115  ibooker.setCurrentFolder(subDirName);
116 
118  ibooker.book1D("LamEffVsR_num", "#Lambda^{0} Efficiency vs #rho", 80, 0., 40.);
120  ibooker.book1D("LamEffVsEta_num", "#Lambda^{0} Efficiency vs #eta", 40, -2.5, 2.5);
122  ibooker.book1D("LamEffVsPt_num", "#Lambda^{0} Efficiency vs p_{T}", 70, 0., 20.);
123 
125  ibooker.book1D("LamTkEffVsR_num", "#Lambda^{0} TrackingEfficiency vs #rho", 80, 0., 40.);
127  ibooker.book1D("LamTkEffVsEta_num", "#Lambda^{0} Tracking Efficiency vs #eta", 40, -2.5, 2.5);
129  ibooker.book1D("LamTkEffVsPt_num", "#Lambda^{0} Tracking Efficiency vs p_{T}", 70, 0., 20.);
130 
132  ibooker.book1D("LamEffVsR_denom", "#Lambda^{0} Efficiency vs #rho", 80, 0., 40.);
134  ibooker.book1D("LamEffVsEta_denom", "#Lambda^{0} Efficiency vs #eta", 40, -2.5, 2.5);
136  ibooker.book1D("LamEffVsPt_denom", "#Lambda^{0} Efficiency vs p_{T}", 70, 0., 20.);
137 
139  ibooker.book1D("LamFakeVsR_num", "#Lambda^{0} Fake Rate vs #rho", 80, 0., 40.);
141  ibooker.book1D("LamFakeVsEta_num", "#Lambda^{0} Fake Rate vs #eta", 40, -2.5, 2.5);
143  ibooker.book1D("LamFakeVsPt_num", "#Lambda^{0} Fake Rate vs p_{T}", 70, 0., 20.);
145  ibooker.book1D("LamTkFakeVsR_num", "#Lambda^{0} Tracking Fake Rate vs #rho", 80, 0., 40.);
147  ibooker.book1D("LamTkFakeVsEta_num", "#Lambda^{0} Tracking Fake Rate vs #eta", 40, -2.5, 2.5);
149  ibooker.book1D("LamTkFakeVsPt_num", "#Lambda^{0} Tracking Fake Rate vs p_{T}", 70, 0., 20.);
150 
152  ibooker.book1D("LamFakeVsR_denom", "#Lambda^{0} Fake Rate vs #rho", 80, 0., 40.);
154  ibooker.book1D("LamFakeVsEta_denom", "#Lambda^{0} Fake Rate vs #eta", 40, -2.5, 2.5);
156  ibooker.book1D("LamFakeVsPt_denom", "#Lambda^{0} Fake Rate vs p_{T}", 70, 0., 20.);
157 
158  nCandidates_[V0Validator::LAMBDA] = ibooker.book1D("nLam", "Number of #Lambda^{0} found per event", 60, 0., 60.);
160  ibooker.book1D("lamMassFake", "Mass of fake Lambda", lamMassNbins, minLamMass, maxLamMass);
162  ibooker.book1D("lamMassGood", "Mass of good Lambda", lamMassNbins, minLamMass, maxLamMass);
163 
165  ibooker.book1D("lamMassAll", "Invariant mass of all #Lambda^{0}", lamMassNbins, lamMassXmin, lamMassXmax);
167  ibooker.book1D("radDistFakeLam", "Production radius of daughter particle of Lam fake", 100, 0., 15.);
168 
169  candidateStatus_[V0Validator::LAMBDA] = ibooker.book1D("ksCandStatus", "Fake type by cand status", 10, 0., 10.);
170 }
171 
173  const reco::RecoToSimCollection& recotosimCollection,
174  V0Type v0_type,
175  int particle_pdgid,
176  int misreconstructed_particle_pdgid) {
177  using namespace edm;
178 
179  int numCandidateFound = 0;
180  int realCandidateFound = 0;
181  double mass = 0.;
182  float CandidatepT = 0.;
183  float CandidateEta = 0.;
184  float CandidateR = 0.;
185  int CandidateStatus = 0;
186  const unsigned int NUM_DAUGHTERS = 2;
187  if (!collection.empty()) {
188  for (reco::VertexCompositeCandidateCollection::const_iterator iCandidate = collection.begin();
189  iCandidate != collection.end();
190  iCandidate++) {
191  // Fill values to be histogrammed
192  mass = iCandidate->mass();
193  CandidatepT = (sqrt(iCandidate->momentum().perp2()));
194  CandidateEta = iCandidate->momentum().eta();
195  CandidateR = (sqrt(iCandidate->vertex().perp2()));
196  candidateMassAll[v0_type]->Fill(mass);
197  CandidateStatus = 0;
198 
199  std::array<reco::TrackRef, NUM_DAUGHTERS> theDaughterTracks = {
200  {(*(dynamic_cast<const reco::RecoChargedCandidate*>(iCandidate->daughter(0)))).track(),
201  (*(dynamic_cast<const reco::RecoChargedCandidate*>(iCandidate->daughter(1)))).track()}};
202 
203  TrackingParticleRef tpref;
204  TrackingParticleRef firstDauTP;
205  TrackingVertexRef candidateVtx;
206 
207  std::array<double, NUM_DAUGHTERS> radDist;
208  // Loop through candidate's daugher tracks
209  for (View<reco::Track>::size_type i = 0; i < theDaughterTracks.size(); ++i) {
210  radDist = {{-1., -1.}};
211  // Found track from theDaughterTracks
212  RefToBase<reco::Track> track(theDaughterTracks.at(i));
213 
214  if (recotosimCollection.find(track) != recotosimCollection.end()) {
215  const std::vector<std::pair<TrackingParticleRef, double> >& tp = recotosimCollection[track];
216  if (!tp.empty()) {
217  tpref = tp.begin()->first;
218 
219  TrackingVertexRef parentVertex = tpref->parentVertex();
220  if (parentVertex.isNonnull()) {
221  radDist[i] = parentVertex->position().R();
222  if (candidateVtx.isNonnull()) {
223  if (candidateVtx->position() == parentVertex->position()) {
224  if (parentVertex->nDaughterTracks() == 2) {
225  if (parentVertex->nSourceTracks() == 0) {
226  // No source tracks found for candidate's
227  // vertex: it shouldn't happen, but does for
228  // evtGen events
229  CandidateStatus = 6;
230  }
231 
232  for (TrackingVertex::tp_iterator iTP = parentVertex->sourceTracks_begin();
233  iTP != parentVertex->sourceTracks_end();
234  iTP++) {
235  if (abs((*iTP)->pdgId()) == particle_pdgid) {
236  CandidateStatus = 1;
237  realCandidateFound++;
238  numCandidateFound += 1.;
239  goodCandidateMass[v0_type]->Fill(mass);
240  } else {
241  CandidateStatus = 2;
242  if (abs((*iTP)->pdgId()) == misreconstructed_particle_pdgid) {
243  CandidateStatus = 7;
244  }
245  }
246  }
247  } else {
248  // Found a bad match because the mother has too
249  // many daughters
250  CandidateStatus = 3;
251  }
252  } else {
253  // Found a bad match because the parent vertices
254  // from the two tracks are different
255  CandidateStatus = 4;
256  }
257  } else {
258  // if candidateVtx is null, fill it with parentVertex
259  // to compare to the parentVertex from the second
260  // track
261  candidateVtx = parentVertex;
262  firstDauTP = tpref;
263  }
264  } // parent vertex is null
265  } // check on associated tp size zero
266  } else {
267  CandidateStatus = 5;
268  }
269  } // Loop on candidate's daughter tracks
270 
271  // fill the fake rate histograms
272  if (CandidateStatus > 1) {
273  candidateFakeVsR_num_[v0_type]->Fill(CandidateR);
274  candidateFakeVsEta_num_[v0_type]->Fill(CandidateEta);
275  candidateFakeVsPt_num_[v0_type]->Fill(CandidatepT);
276  candidateStatus_[v0_type]->Fill((float)CandidateStatus);
277  fakeCandidateMass_[v0_type]->Fill(mass);
278  for (auto distance : radDist) {
279  if (distance > 0)
280  candidateFakeDauRadDist_[v0_type]->Fill(distance);
281  }
282  }
283  if (CandidateStatus == 5) {
284  candidateTkFakeVsR_num_[v0_type]->Fill(CandidateR);
285  candidateTkFakeVsEta_num_[v0_type]->Fill(CandidateEta);
286  candidateTkFakeVsPt_num_[v0_type]->Fill(CandidatepT);
287  }
288  candidateFakeVsR_denom_[v0_type]->Fill(CandidateR);
289  candidateFakeVsEta_denom_[v0_type]->Fill(CandidateEta);
290  candidateFakeVsPt_denom_[v0_type]->Fill(CandidatepT);
291  } // Loop on candidates
292  } // check on presence of candidate's collection in the event
293  nCandidates_[v0_type]->Fill((float)numCandidateFound);
294 }
295 
297  V0Type v0_type,
298  int parent_particle_id,
299  int first_daughter_id, /* give only positive charge */
300  int second_daughter_id, /* give only positive charge */
302  const reco::SimToRecoCollection& simtorecoCollection) {
303  /* We store the TrackRef of the tracks that have been used to
304  * produce the V0 under consideration here. This is used later to
305  * check if a specific V0 has been really reconstructed or not. The
306  * ordering is based on the key_index of the reference, since it
307  * indeed does not matter that much. */
308 
309  std::set<V0Couple> reconstructed_V0_couples;
310  if (!collection.empty()) {
311  for (reco::VertexCompositeCandidateCollection::const_iterator iCandidate = collection.begin();
312  iCandidate != collection.end();
313  iCandidate++) {
314  reconstructed_V0_couples.insert(
315  V0Couple((dynamic_cast<const reco::RecoChargedCandidate*>(iCandidate->daughter(0)))->track(),
316  (dynamic_cast<const reco::RecoChargedCandidate*>(iCandidate->daughter(1)))->track()));
317  }
318  }
319 
320  /* PSEUDO CODE
321  for v in gen_vertices
322  if v.eventId().BX() !=0 continue
323  if v.nDaughterTracks != 2 continue
324  for source in v.sourceTracks_begin
325  if source is parent_particle_id
326  for daughter in v.daughterTracks_begin
327  if daughter in region_and_kine_cuts
328  decay_found
329  */
330  unsigned int candidateEff[2] = {0, 0};
331  for (auto const& gen_vertex : gen_vertices) {
332  if (gen_vertex.eventId().bunchCrossing() != 0)
333  continue; // Consider only in-time events
334  if (gen_vertex.nDaughterTracks() != 2)
335  continue; // Keep only V0 vertices
336  for (TrackingVertex::tp_iterator source = gen_vertex.sourceTracks_begin(); source != gen_vertex.sourceTracks_end();
337  ++source) {
338  if (std::abs((*source)->pdgId()) == parent_particle_id) {
339  if ((std::abs((gen_vertex.daughterTracks().at(0))->pdgId()) == first_daughter_id &&
340  std::abs((gen_vertex.daughterTracks().at(1))->pdgId()) == second_daughter_id) ||
341  (std::abs((gen_vertex.daughterTracks().at(0))->pdgId()) == second_daughter_id &&
342  std::abs((gen_vertex.daughterTracks().at(1))->pdgId()) == first_daughter_id)) {
343  if ((std::abs((gen_vertex.daughterTracks().at(0))->momentum().eta()) < 2.4 &&
344  gen_vertex.daughterTracks().at(0)->pt() > 0.9) &&
345  (std::abs((gen_vertex.daughterTracks().at(1))->momentum().eta()) < 2.4 &&
346  gen_vertex.daughterTracks().at(1)->pt() > 0.9)) {
347  // found desired generated Candidate
348  float candidateGenpT = sqrt((*source)->momentum().perp2());
349  float candidateGenEta = (*source)->momentum().eta();
350  float candidateGenR = sqrt((*source)->vertex().perp2());
351  candidateEffVsPt_denom_[v0_type]->Fill(candidateGenpT);
352  candidateEffVsEta_denom_[v0_type]->Fill(candidateGenEta);
353  candidateEffVsR_denom_[v0_type]->Fill(candidateGenR);
354 
355  std::array<reco::TrackRef, 2> reco_daughter;
356 
357  for (unsigned int daughter = 0; daughter < 2; ++daughter) {
358  if (simtorecoCollection.find(gen_vertex.daughterTracks()[daughter]) != simtorecoCollection.end()) {
359  if (!simtorecoCollection[gen_vertex.daughterTracks()[daughter]].empty()) {
360  candidateEff[daughter] = 1; // Found a daughter track
361  reco_daughter[daughter] =
362  simtorecoCollection[gen_vertex.daughterTracks()[daughter]].begin()->first.castTo<reco::TrackRef>();
363  }
364  } else {
365  candidateEff[daughter] = 2; // First daughter not found
366  }
367  }
368  if ((candidateEff[0] == 1 && candidateEff[1] == 1) && (reco_daughter[0].key() != reco_daughter[1].key()) &&
369  (reconstructed_V0_couples.find(V0Couple(reco_daughter[0], reco_daughter[1])) !=
370  reconstructed_V0_couples.end())) {
371  candidateEffVsPt_num_[v0_type]->Fill(candidateGenpT);
372  candidateEffVsEta_num_[v0_type]->Fill(candidateGenEta);
373  candidateEffVsR_num_[v0_type]->Fill(candidateGenR);
374  }
375  } // Check that daughters are inside the desired kinematic region
376  } // Check decay products of the current generatex vertex
377  } // Check pdgId of the source of the current generated vertex
378  } // Loop over all sources of the current generated vertex
379  } // Loop over all generated vertices
380 }
381 
383  using std::cout;
384  using std::endl;
385  using namespace edm;
386  using namespace std;
387 
388  // Make matching collections
389  Handle<reco::RecoToSimCollection> recotosimCollectionH;
390  iEvent.getByToken(recoRecoToSimCollectionToken_, recotosimCollectionH);
391 
392  Handle<reco::SimToRecoCollection> simtorecoCollectionH;
393  iEvent.getByToken(recoSimToRecoCollectionToken_, simtorecoCollectionH);
394 
395  // Get Monte Carlo information
397  iEvent.getByToken(trackingVertexCollection_Token_, TVCollectionH);
398 
399  // Select the primary vertex, create a new reco::Vertex to hold it
400  edm::Handle<std::vector<reco::Vertex> > primaryVtxCollectionH;
401  iEvent.getByToken(vec_recoVertex_Token_, primaryVtxCollectionH);
402 
403  std::vector<reco::Vertex>::const_iterator iVtxPH = primaryVtxCollectionH->begin();
404  for (std::vector<reco::Vertex>::const_iterator iVtx = primaryVtxCollectionH->begin();
405  iVtx < primaryVtxCollectionH->end();
406  iVtx++) {
407  if (primaryVtxCollectionH->size() > 1) {
408  if (iVtx->tracksSize() > iVtxPH->tracksSize()) {
409  iVtxPH = iVtx;
410  }
411  } else
412  iVtxPH = iVtx;
413  }
414 
415  // get the V0s;
420 
421  // Do fake rate and efficiency calculation
422 
423  // Get gen vertex collection out of the event, as done in the Vertex
424  // validation package!!!
425  if (k0sCollection.isValid()) {
426  doFakeRates(*k0sCollection.product(), *recotosimCollectionH.product(), V0Type::KSHORT, 310, 3122);
427  doEfficiencies(*TVCollectionH.product(),
428  V0Type::KSHORT,
429  310,
430  211,
431  211,
432  *k0sCollection.product(),
433  *simtorecoCollectionH.product());
434  }
435  if (lambdaCollection.isValid()) {
436  doFakeRates(*lambdaCollection.product(), *recotosimCollectionH.product(), V0Type::LAMBDA, 3122, 310);
437  doEfficiencies(*TVCollectionH.product(),
438  V0Type::LAMBDA,
439  3122,
440  211,
441  2212,
442  *lambdaCollection.product(),
443  *simtorecoCollectionH.product());
444  }
445 }
446 
447 // define this as a plug-in
448 // DEFINE_FWK_MODULE(V0Validator);
unsigned int size_type
Definition: View.h:90
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: V0Validator.cc:45
std::vector< VertexCompositeCandidate > VertexCompositeCandidateCollection
collection of Candidate objects
edm::EDGetTokenT< reco::SimToRecoCollection > recoSimToRecoCollectionToken_
Definition: V0Validator.h:145
V0Validator(const edm::ParameterSet &)
Definition: V0Validator.cc:27
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > recoVertexCompositeCandidateCollection_k0s_Token_
Definition: V0Validator.h:148
std::array< MonitorElement *, 2 > candidateStatus_
Definition: V0Validator.h:136
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
edm::RefVector< edm::HepMCProduct, HepMC::GenParticle > GenParticleRefVector
Definition: V0Validator.cc:25
std::array< MonitorElement *, 2 > candidateEffVsEta_num_
Definition: V0Validator.h:116
std::array< MonitorElement *, 2 > candidateFakeVsEta_num_
Definition: V0Validator.h:122
T const * product() const
Definition: Handle.h:70
std::array< MonitorElement *, 2 > candidateFakeVsR_denom_
Definition: V0Validator.h:128
std::array< MonitorElement *, 2 > goodCandidateMass
Definition: V0Validator.h:140
std::array< MonitorElement *, 2 > candidateTkEffVsPt_num_
Definition: V0Validator.h:120
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
std::array< MonitorElement *, 2 > fakeCandidateMass_
Definition: V0Validator.h:137
edm::AssociationMap< edm::OneToManyWithQualityGeneric< CaloParticleCollection, reco::CaloClusterCollection, std::pair< float, float > > > SimToRecoCollection
std::array< MonitorElement *, 2 > candidateEffVsPt_num_
Definition: V0Validator.h:117
std::array< MonitorElement *, 2 > candidateEffVsR_denom_
Definition: V0Validator.h:131
std::array< MonitorElement *, 2 > candidateFakeVsR_num_
Definition: V0Validator.h:121
bool empty() const
return true if empty
const_iterator find(const key_type &k) const
find element with specified reference key
std::array< MonitorElement *, 2 > candidateFakeDauRadDist_
Definition: V0Validator.h:138
std::array< MonitorElement *, 2 > candidateTkFakeVsEta_num_
Definition: V0Validator.h:125
std::array< MonitorElement *, 2 > candidateTkEffVsEta_num_
Definition: V0Validator.h:119
const_iterator end() const
last iterator over the map (read only)
edm::Ref< TrackingVertexCollection > TrackingVertexRef
Definition: V0Validator.cc:23
int iEvent
Definition: GenABIO.cc:224
std::array< MonitorElement *, 2 > candidateTkEffVsR_num_
Definition: V0Validator.h:118
T sqrt(T t)
Definition: SSEVec.h:19
~V0Validator() override
Definition: V0Validator.cc:43
std::array< MonitorElement *, 2 > candidateFakeVsPt_num_
Definition: V0Validator.h:123
std::array< MonitorElement *, 2 > candidateTkFakeVsR_num_
Definition: V0Validator.h:124
edm::EDGetTokenT< TrackingVertexCollection > trackingVertexCollection_Token_
Definition: V0Validator.h:146
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: V0Validator.cc:382
std::array< MonitorElement *, 2 > candidateTkFakeVsPt_num_
Definition: V0Validator.h:126
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< std::vector< reco::Vertex > > vec_recoVertex_Token_
Definition: V0Validator.h:147
std::array< MonitorElement *, 2 > nCandidates_
Definition: V0Validator.h:135
void doFakeRates(const reco::VertexCompositeCandidateCollection &collection, const reco::RecoToSimCollection &recotosimCollection, V0Type t, int particle_pdgid, int misreconstructed_particle_pdgid)
Definition: V0Validator.cc:172
std::array< MonitorElement *, 2 > candidateFakeVsPt_denom_
Definition: V0Validator.h:130
std::vector< TrackingVertex > TrackingVertexCollection
edm::RefVector< edm::HepMCProduct, HepMC::GenVertex > GenVertexRefVector
Definition: V0Validator.cc:24
std::array< MonitorElement *, 2 > candidateFakeVsEta_denom_
Definition: V0Validator.h:129
bool isValid() const
Definition: HandleBase.h:70
const_iterator begin() const
first iterator over the map (read only)
fixed size matrix
HLT enums.
std::array< MonitorElement *, 2 > candidateEffVsPt_denom_
Definition: V0Validator.h:133
std::array< MonitorElement *, 2 > candidateMassAll
Definition: V0Validator.h:139
edm::EDGetTokenT< reco::RecoToSimCollection > recoRecoToSimCollectionToken_
Definition: V0Validator.h:144
edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > recoVertexCompositeCandidateCollection_lambda_Token_
Definition: V0Validator.h:148
std::array< MonitorElement *, 2 > candidateEffVsR_num_
Definition: V0Validator.h:115
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
std::array< MonitorElement *, 2 > candidateEffVsEta_denom_
Definition: V0Validator.h:132
std::string dirName
Definition: V0Validator.h:143
void doEfficiencies(const TrackingVertexCollection &gen_vertices, V0Type t, int parent_particle_id, int first_daughter_id, int second_daughter_id, const reco::VertexCompositeCandidateCollection &collection, const reco::SimToRecoCollection &simtorecoCollection)
Definition: V0Validator.cc:296
static std::string const source
Definition: EdmProvDump.cc:49
Definition: Run.h:45
edm::AssociationMap< edm::OneToManyWithQualityGeneric< reco::CaloClusterCollection, CaloParticleCollection, float > > RecoToSimCollection
std::vector< TrackingVertex > TrackingVertexCollection
Definition: V0Validator.cc:22