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