CMS 3D CMS Logo

MuScleFitMuonSelector.cc
Go to the documentation of this file.
5 
6 const double MuScleFitMuonSelector::mMu2 = 0.011163612;
7 const unsigned int MuScleFitMuonSelector::motherPdgIdArray[] = {23, 100553, 100553, 553, 100443, 443};
8 
9 const reco::Candidate*
11  const reco::Candidate* tempMuon = status3Muon;
12  // bool lastCopy = ((reco::GenParticle*)tempMuon)->isLastCopy(); // isLastCopy() likely not enough robust
13  bool isPromptFinalState = ((reco::GenParticle*)tempMuon)->isPromptFinalState(); // pre-CMSSW_74X code: int status = tempStatus1Muon->status();
14  while(tempMuon == nullptr || tempMuon->numberOfDaughters()!=0){
15  if ( isPromptFinalState ) break; // pre-CMSSW_74X code: if (status == 1) break;
16  //std::vector<const reco::Candidate*> daughters;
17  for (unsigned int i=0; i<tempMuon->numberOfDaughters(); ++i){
18  if ( tempMuon->daughter(i)->pdgId()==tempMuon->pdgId() ){
19  tempMuon = tempMuon->daughter(i);
20  isPromptFinalState = ((reco::GenParticle*)tempMuon)->isPromptFinalState(); // pre-CMSSW_74X code: status = tempStatus1Muon->status();
21  break;
22  }else continue;
23  }//for loop
24  }//while loop
25 
26  return tempMuon;
27 }
28 
29 const reco::Candidate*
31  const reco::Candidate* tempMuon = status3Muon;
32  bool lastCopy = ((reco::GenParticle*)tempMuon)->isLastCopyBeforeFSR(); // pre-CMSSW_74X code: int status = tempStatus1Muon->status();
33  while(tempMuon == nullptr || tempMuon->numberOfDaughters()!=0){
34  if ( lastCopy ) break; // pre-CMSSW_74X code: if (status == 3) break;
35  //std::vector<const reco::Candidate*> daughters;
36  for (unsigned int i=0; i<tempMuon->numberOfDaughters(); ++i){
37  if ( tempMuon->daughter(i)->pdgId()==tempMuon->pdgId() ){
38  tempMuon = tempMuon->daughter(i);
39  lastCopy = ((reco::GenParticle*)tempMuon)->isLastCopyBeforeFSR(); // pre-CMSSW_74X code: status = tempStatus1Muon->status();
40  break;
41  }else continue;
42  }//for loop
43  }//while loop
44 
45  return tempMuon;
46 }
47 
48 
49 
51 {
52  reco::TrackRef iTrack = aMuon->innerTrack();
53  const reco::HitPattern& p = iTrack->hitPattern();
54 
55  reco::TrackRef gTrack = aMuon->globalTrack();
56  const reco::HitPattern& q = gTrack->hitPattern();
57 
58  return (//isMuonInAccept(aMuon) &&// no acceptance cuts!
59  iTrack->found() > 11 &&
60  gTrack->chi2()/gTrack->ndof() < 20.0 &&
61  q.numberOfValidMuonHits() > 0 &&
62  iTrack->chi2()/iTrack->ndof() < 4.0 &&
63  aMuon->muonID("TrackerMuonArbitrated") &&
64  aMuon->muonID("TMLastStationAngTight") &&
66  std::abs(iTrack->dxy()) < 3.0 && //should be done w.r.t. PV!
67  std::abs(iTrack->dz()) < 15.0 //should be done w.r.t. PV!
68  );
69 }
70 
72 {
73  reco::TrackRef iTrack = aMuon->innerTrack();
74  const reco::HitPattern& p = iTrack->hitPattern();
75 
76  return (//isMuonInAccept(aMuon) // no acceptance cuts!
77  iTrack->found() > 11 &&
78  iTrack->chi2()/iTrack->ndof() < 4.0 &&
79  aMuon->muonID("TrackerMuonArbitrated") &&
80  aMuon->muonID("TMLastStationAngTight") &&
82  std::abs(iTrack->dxy()) < 3.0 && //should be done w.r.t. PV!
83  std::abs(iTrack->dz()) < 15.0 //should be done w.r.t. PV!
84  );
85 }
86 
87 // Note that at this level we save all the information. Events for which no suitable muon pair is found
88 // are removed from the tree (together with the gen and sim information) at a later stage.
89 // It would be better to remove them directly at this point.
90 void MuScleFitMuonSelector::selectMuons(const edm::Event & event, std::vector<MuScleFitMuon> & muons,
91  std::vector<GenMuonPair> & genPair,
92  std::vector<std::pair<lorentzVector, lorentzVector> > & simPair,
94 {
96  try { event.getByLabel("onia2MuMuPatTrkTrk", collAll); }
97  catch (...) { std::cout << "J/psi not present in event!" << std::endl; }
98  std::vector<const pat::Muon*> collMuSel;
99 
100  //================onia cuts===========================/
101 
102  if (muonType_ <= -1 && PATmuons_) {
103  std::vector<const pat::CompositeCandidate*> collSelGG;
104  std::vector<const pat::CompositeCandidate*> collSelGT;
105  std::vector<const pat::CompositeCandidate*> collSelTT;
106  if (collAll.isValid()) {
107 
108  for (std::vector<pat::CompositeCandidate>::const_iterator it=collAll->begin();
109  it!=collAll->end(); ++it) {
110 
111  const pat::CompositeCandidate* cand = &(*it);
112  // cout << "Now checking candidate of type " << theJpsiCat << " with pt = " << cand->pt() << endl;
113  const pat::Muon* muon1 = dynamic_cast<const pat::Muon*>(cand->daughter("muon1"));
114  const pat::Muon* muon2 = dynamic_cast<const pat::Muon*>(cand->daughter("muon2"));
115 
116  if ((muon1->charge() * muon2->charge())>0)
117  continue;
118  // global + global?
119  if (muon1->isGlobalMuon() && muon2->isGlobalMuon()) {
120  if (selGlobalMuon(muon1) && selGlobalMuon(muon2) && cand->userFloat("vProb") > 0.001) {
121  collSelGG.push_back(cand);
122  continue;
123  }
124  }
125  // global + tracker? (x2)
126  if (muon1->isGlobalMuon() && muon2->isTrackerMuon()) {
127  if (selGlobalMuon(muon1) && selTrackerMuon(muon2) && cand->userFloat("vProb") > 0.001) {
128  collSelGT.push_back(cand);
129  continue;
130  }
131  }
132  if (muon2->isGlobalMuon() && muon1->isTrackerMuon()) {
133  if (selGlobalMuon(muon2) && selTrackerMuon(muon1) && cand->userFloat("vProb") > 0.001) {
134  collSelGT.push_back(cand);
135  continue;
136  }
137  }
138  // tracker + tracker?
139  if (muon1->isTrackerMuon() && muon2->isTrackerMuon()) {
140  if (selTrackerMuon(muon1) && selTrackerMuon(muon2) && cand->userFloat("vProb") > 0.001) {
141  collSelTT.push_back(cand);
142  continue;
143  }
144  }
145  }
146  }
147  // Split them in independent collections if using muonType_ == -2, -3 or -4. Take them all if muonType_ == -1.
148  std::vector<reco::Track> tracks;
149  if (!collSelGG.empty()){
150  //CHECK THAT THEY ARE ORDERED BY PT !!!!!!!!!!!!!!!!!!!!!!!
151  const pat::Muon* muon1 = dynamic_cast<const pat::Muon*>(collSelGG[0]->daughter("muon1"));
152  const pat::Muon* muon2 = dynamic_cast<const pat::Muon*>(collSelGG[0]->daughter("muon2"));
153  if (muonType_ == -1 || muonType_ == -2) {
154  tracks.push_back(*(muon1->innerTrack()));
155  tracks.push_back(*(muon2->innerTrack()));
156  collMuSel.push_back(muon1);
157  collMuSel.push_back(muon2);
158  }
159  }
160  else if (collSelGG.empty() && !collSelGT.empty()){
161  //CHECK THAT THEY ARE ORDERED BY PT !!!!!!!!!!!!!!!!!!!!!!!
162  const pat::Muon* muon1 = dynamic_cast<const pat::Muon*>(collSelGT[0]->daughter("muon1"));
163  const pat::Muon* muon2 = dynamic_cast<const pat::Muon*>(collSelGT[0]->daughter("muon2"));
164  if (muonType_ == -1 || muonType_ == -3) {
165  tracks.push_back(*(muon1->innerTrack()));
166  tracks.push_back(*(muon2->innerTrack()));
167  collMuSel.push_back(muon1);
168  collMuSel.push_back(muon2);
169  }
170  }
171  else if (collSelGG.empty() && collSelGT.empty() && !collSelTT.empty()){
172  //CHECK THAT THEY ARE ORDERED BY PT !!!!!!!!!!!!!!!!!!!!!!!
173  const pat::Muon* muon1 = dynamic_cast<const pat::Muon*>(collSelTT[0]->daughter("muon1"));
174  const pat::Muon* muon2 = dynamic_cast<const pat::Muon*>(collSelTT[0]->daughter("muon2"));
175  if (muonType_ == -1 || muonType_ == -4) {
176  tracks.push_back(*(muon1->innerTrack()));
177  tracks.push_back(*(muon2->innerTrack()));
178  collMuSel.push_back(muon1);
179  collMuSel.push_back(muon2);
180  }
181  }
182  if (tracks.size() != 2 && !tracks.empty()){
183  std::cout<<"ERROR strange number of muons selected by onia cuts!"<<std::endl;
184  abort();
185  }
186  muons = fillMuonCollection(tracks);
187  }
188  else if ((muonType_<4 && muonType_>=0) || muonType_>=10) { // Muons (glb,sta,trk)
189  std::vector<reco::Track> tracks;
190  if (PATmuons_ == true) {
192  event.getByLabel(muonLabel_, allMuons);
193  if (muonType_ == 0) {
194  // Take directly the muon
195  muons = fillMuonCollection(*allMuons);
196  }
197  else {
198  for (std::vector<pat::Muon>::const_iterator muon = allMuons->begin(); muon != allMuons->end(); ++muon) {
199  //std::cout<<"pat muon is global "<<muon->isGlobalMuon()<<std::endl;
200  takeSelectedMuonType(muon, tracks);
201  }
202  muons = fillMuonCollection(tracks);
203  }
204  }
205  else {
207  event.getByLabel(muonLabel_, allMuons);
208  if (muonType_ == 0) {
209  // Take directly the muon
210  muons = fillMuonCollection(*allMuons);
211  }
212  else {
213  for (std::vector<reco::Muon>::const_iterator muon = allMuons->begin(); muon != allMuons->end(); ++muon) {
214  takeSelectedMuonType(muon, tracks);
215  }
216  muons = fillMuonCollection(tracks);
217  }
218  }
219  }
220  else if (muonType_==4){ //CaloMuons
222  event.getByLabel(muonLabel_, caloMuons);
223  std::vector<reco::Track> tracks;
224  for (std::vector<reco::CaloMuon>::const_iterator muon = caloMuons->begin(); muon != caloMuons->end(); ++muon){
225  tracks.push_back(*(muon->track()));
226  }
227  muons = fillMuonCollection(tracks);
228  }
229 
230  else if (muonType_==5) { // Inner tracker tracks
232  event.getByLabel(muonLabel_, tracks);
233  muons = fillMuonCollection(*tracks);
234  }
235  plotter->fillRec(muons);
236 
237  // Generation and simulation part
238  if (speedup_ == false)
239  {
240  if (PATmuons_) {
241  // EM 2015.04.02 temporary fix to run on MINIAODSIM (which contains PAT muons) but not the "onia2MuMuPatTrkTrk" collection
242  // selectGeneratedMuons(collAll, collMuSel, genPair, plotter);
243  selectGenSimMuons(event, genPair, simPair, plotter);
244  }
245  else {
246  selectGenSimMuons(event, genPair, simPair, plotter);
247  }
248  }
249 }
250 
253  const std::vector<const pat::Muon*> & collMuSel,
254  std::vector<GenMuonPair> & genPair,
256 {
258 
259  //explicitly for JPsi but can be adapted!!!!!
260  for(std::vector<pat::CompositeCandidate>::const_iterator it=collAll->begin();
261  it!=collAll->end();++it) {
262  reco::GenParticleRef genJpsi = it->genParticleRef();
263  bool isMatched = (genJpsi.isAvailable() && genJpsi->pdgId() == 443);
264  if (isMatched){
265  genPatParticles->push_back(*(const_cast<reco::GenParticle*>(genJpsi.get())));
266  }
267  }
268 
269  if(collMuSel.size()==2) {
270  reco::GenParticleRef genMu1 = collMuSel[0]->genParticleRef();
271  reco::GenParticleRef genMu2 = collMuSel[1]->genParticleRef();
272  bool isMuMatched = (genMu1.isAvailable() && genMu2.isAvailable() &&
273  genMu1->pdgId()*genMu2->pdgId() == -169);
274  if (isMuMatched) {
275  genPatParticles->push_back(*(const_cast<reco::GenParticle*>(genMu1.get())));
276  genPatParticles->push_back(*(const_cast<reco::GenParticle*>(genMu2.get())));
277 
278  unsigned int motherId = 0;
279  if( genMu1->mother() != nullptr ) {
280  motherId = genMu1->mother()->pdgId();
281  }
282  if(genMu1->pdgId()==13)
283  genPair.push_back(GenMuonPair(genMu1.get()->p4(), genMu2.get()->p4(), motherId));
284  else
285  // genPair.push_back(std::make_pair(genMu2.get()->p4(),genMu1.get()->p4()) );
286  genPair.push_back(GenMuonPair(genMu2.get()->p4(), genMu1.get()->p4(), motherId));
287 
288  plotter->fillGen(const_cast <reco::GenParticleCollection*> (genPatParticles), true);
289 
290  if (debug_>0) std::cout << "Found genParticles in PAT" << std::endl;
291  }
292  else {
293  std::cout << "No recomuon selected so no access to generated info"<<std::endl;
294  // Fill it in any case, otherwise it will not be in sync with the event number
295  // genPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
296  genPair.push_back( GenMuonPair(lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.), 0 ) );
297  }
298  }
299  else{
300  std::cout << "No recomuon selected so no access to generated info"<<std::endl;
301  // Fill it in any case, otherwise it will not be in sync with the event number
302  // genPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
303  genPair.push_back( GenMuonPair(lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.), 0 ) );
304  }
305  if(debug_>0) {
306  std::cout << "genParticles:" << std::endl;
307  // std::cout << genPair.back().first << " , " << genPair.back().second << std::endl;
308  std::cout << genPair.back().mu1 << " , " << genPair.back().mu2 << std::endl;
309  }
310 }
311 
313  std::vector<GenMuonPair> & genPair,
314  std::vector<std::pair<lorentzVector,lorentzVector> > & simPair,
316 {
317  // Find and store in histograms the generated and simulated resonance and muons
318  // ----------------------------------------------------------------------------
321 
322  // Fill gen information only in the first loop
323  bool ifHepMC=false;
324 
325  event.getByLabel( genParticlesName_, evtMC );
326  event.getByLabel( genParticlesName_, genParticles );
327  if( evtMC.isValid() ) {
328  genPair.push_back( findGenMuFromRes(evtMC.product()) );
329  plotter->fillGen(evtMC.product(), sherpa_);
330  ifHepMC = true;
331  if (debug_>0) std::cout << "Found hepMC" << std::endl;
332  }
333  else if( genParticles.isValid() ) {
334  genPair.push_back( findGenMuFromRes(genParticles.product()) );
335  plotter->fillGen(genParticles.product());
336  if (debug_>0) std::cout << "Found genParticles" << std::endl;
337  }
338  else {
339  std::cout<<"ERROR "<<"non generation info and speedup true!!!!!!!!!!!!"<<std::endl;
340  // Fill it in any case, otherwise it will not be in sync with the event number
341  // genPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
342  genPair.push_back( GenMuonPair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.), 0 ) );
343  }
344  if(debug_>0) {
345  std::cout << "genParticles:" << std::endl;
346  std::cout << genPair.back().mu1 << " , " << genPair.back().mu2 << std::endl;
347  }
348  selectSimulatedMuons(event, ifHepMC, evtMC, simPair, plotter);
349 }
350 
352  const bool ifHepMC, edm::Handle<edm::HepMCProduct> evtMC,
353  std::vector<std::pair<lorentzVector,lorentzVector> > & simPair,
355 {
357  bool simTracksFound = false;
358  event.getByLabel(simTracksCollectionName_, simTracks);
359  if( simTracks.isValid() ) {
360  plotter->fillSim(simTracks);
361  if(ifHepMC) {
362  simPair.push_back( findSimMuFromRes(evtMC, simTracks) );
363  simTracksFound = true;
364  plotter->fillGenSim(evtMC,simTracks);
365  }
366  }
367  else {
368  std::cout << "SimTracks not existent" << std::endl;
369  }
370  if( !simTracksFound ) {
371  simPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
372  }
373 }
374 
376 {
377  const HepMC::GenEvent* Evt = evtMC->GetEvent();
378  GenMuonPair muFromRes;
379  //Loop on generated particles
380  for (HepMC::GenEvent::particle_const_iterator part=Evt->particles_begin();
381  part!=Evt->particles_end(); part++) {
382  if (std::abs((*part)->pdg_id())==13 && (*part)->status()==1) {
383  bool fromRes = false;
384  unsigned int motherPdgId = 0;
385  for (HepMC::GenVertex::particle_iterator mother = (*part)->production_vertex()->particles_begin(HepMC::ancestors);
386  mother != (*part)->production_vertex()->particles_end(HepMC::ancestors); ++mother) {
387  motherPdgId = (*mother)->pdg_id();
388 
389  // For sherpa the resonance is not saved. The muons from the resonance can be identified
390  // by having as mother a muon of status 3.
391  if( sherpa_ ) {
392  if( motherPdgId == 13 && (*mother)->status() == 3 ) fromRes = true;
393  }
394  else {
395  for( int ires = 0; ires < 6; ++ires ) {
396  if( motherPdgId == motherPdgIdArray[ires] && resfind_[ires] ) fromRes = true;
397  }
398  }
399  }
400  if(fromRes){
401  if((*part)->pdg_id()==13) {
402  // muFromRes.first = (*part)->momentum();
403  muFromRes.mu1 = (lorentzVector((*part)->momentum().px(),(*part)->momentum().py(),
404  (*part)->momentum().pz(),(*part)->momentum().e()));
405  }
406  else {
407  muFromRes.mu2 = (lorentzVector((*part)->momentum().px(),(*part)->momentum().py(),
408  (*part)->momentum().pz(),(*part)->momentum().e()));
409  }
410  muFromRes.motherId = motherPdgId;
411  }
412  }
413  }
414  return muFromRes;
415 }
416 
417 // std::pair<lorentzVector, lorentzVector> MuScleFitMuonSelector::findGenMuFromRes( const reco::GenParticleCollection* genParticles)
419 {
420  // std::pair<lorentzVector,lorentzVector> muFromRes;
421  GenMuonPair muFromRes;
422 
423  //Loop on generated particles
424  if( debug_>0 ) std::cout << "Starting loop on " << genParticles->size() << " genParticles" << std::endl;
425  for( reco::GenParticleCollection::const_iterator part=genParticles->begin(); part!=genParticles->end(); ++part ) {
426  if (debug_>0) std::cout<<"genParticle has pdgId = "<<std::abs(part->pdgId())<<" and status = "<<part->status()<<std::endl;
427  if (std::abs(part->pdgId())==13){// && part->status()==3) {
428  bool fromRes = false;
429  unsigned int motherPdgId = part->mother()->pdgId();
430  if( debug_>0 ) {
431  std::cout << "Found a muon with mother: " << motherPdgId << std::endl;
432  }
433  for( int ires = 0; ires < 6; ++ires ) {
434  // if( motherPdgId == motherPdgIdArray[ires] && resfind_[ires] ) fromRes = true; // changed by EM 2015.07.30
435  // begin of comment
436  // the list of resonances motherPdgIdArray does not contain the photon (PdgId = 21) while ~1% of the
437  // mu+mu- events in the range [50,120] GeV has a photon as the mother.
438  // It needs to be fixed without spoiling the logic of the selection of different resonances
439  // e.g. mixing up onia etc.
440  // end of comment
441  if( ( motherPdgId == motherPdgIdArray[ires] && resfind_[ires] ) || motherPdgId == 21 ) fromRes = true;
442  }
443  if(fromRes){
444  if (debug_>0) std::cout<<"fromRes = true, motherPdgId = "<<motherPdgId<<std::endl;
445  const reco::Candidate* status3Muon = &(*part);
446  const reco::Candidate* status1Muon = getStatus1Muon(status3Muon);
447  if(part->pdgId()==13) {
448  if (status1Muon->p4().pt()!=0) muFromRes.mu1 = MuScleFitMuon(status1Muon->p4(),-1);
449  else muFromRes.mu1 = MuScleFitMuon(status3Muon->p4(),-1);
450  if( debug_>0 ) std::cout << "Found a genMuon - : " << muFromRes.mu1 << std::endl;
451  // muFromRes.first = (lorentzVector(status1Muon->p4().px(),status1Muon->p4().py(),
452  // status1Muon->p4().pz(),status1Muon->p4().e()));
453  }
454  else {
455  if (status1Muon->p4().pt()!=0) muFromRes.mu2 = MuScleFitMuon(status1Muon->p4(),+1);
456  else muFromRes.mu2 = MuScleFitMuon(status3Muon->p4(),+1);
457  if( debug_>0 ) std::cout << "Found a genMuon + : " << muFromRes.mu2 << std::endl;
458  // muFromRes.second = (lorentzVector(status1Muon->p4().px(),status1Muon->p4().py(),
459  // status1Muon->p4().pz(),status1Muon->p4().e()));
460  }
461  muFromRes.motherId = motherPdgId;
462  }
463  }
464  }
465  return muFromRes;
466 }
467 
468 std::pair<lorentzVector, lorentzVector> MuScleFitMuonSelector::findSimMuFromRes( const edm::Handle<edm::HepMCProduct> & evtMC,
470 {
471  //Loop on simulated tracks
472  std::pair<lorentzVector, lorentzVector> simMuFromRes;
473  for( edm::SimTrackContainer::const_iterator simTrack=simTracks->begin(); simTrack!=simTracks->end(); ++simTrack ) {
474  //Chose muons
475  if (std::abs((*simTrack).type())==13) {
476  //If tracks from IP than find mother
477  if ((*simTrack).genpartIndex()>0) {
478  HepMC::GenParticle* gp = evtMC->GetEvent()->barcode_to_particle ((*simTrack).genpartIndex());
479  if( gp != nullptr ) {
480 
481  for (HepMC::GenVertex::particle_iterator mother = gp->production_vertex()->particles_begin(HepMC::ancestors);
482  mother!=gp->production_vertex()->particles_end(HepMC::ancestors); ++mother) {
483 
484  bool fromRes = false;
485  unsigned int motherPdgId = (*mother)->pdg_id();
486  for( int ires = 0; ires < 6; ++ires ) {
487  if( ( motherPdgId == motherPdgIdArray[ires] && resfind_[ires] ) || motherPdgId == 21 ) fromRes = true;
488  }
489  if( fromRes ) {
490  if(gp->pdg_id() == 13)
491  simMuFromRes.first = lorentzVector(simTrack->momentum().px(),simTrack->momentum().py(),
492  simTrack->momentum().pz(),simTrack->momentum().e());
493  else
494  simMuFromRes.second = lorentzVector(simTrack->momentum().px(),simTrack->momentum().py(),
495  simTrack->momentum().pz(),simTrack->momentum().e());
496  }
497  }
498  }
499  // else LogDebug("MuScleFitUtils") << "WARNING: no matching genParticle found for simTrack" << std::endl;
500  }
501  }
502  }
503  return simMuFromRes;
504 }
Analysis-level particle class.
bool isAvailable() const
Definition: Ref.h:577
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
void selectGeneratedMuons(const edm::Handle< pat::CompositeCandidateCollection > &collAll, const std::vector< const pat::Muon * > &collMuSel, std::vector< GenMuonPair > &genPair, MuScleFitPlotter *plotter)
For PATmuons the generator information is read directly from the PAT object.
std::pair< lorentzVector, lorentzVector > findSimMuFromRes(const edm::Handle< edm::HepMCProduct > &evtMC, const edm::Handle< edm::SimTrackContainer > &simTracks)
const std::string genParticlesName_
void selectMuons(const edm::Event &event, std::vector< MuScleFitMuon > &muons, std::vector< GenMuonPair > &genPair, std::vector< std::pair< lorentzVector, lorentzVector > > &simPair, MuScleFitPlotter *plotter)
Main method used to select muons of type specified by muonType_ from the collection specified by muon...
void fillGenSim(edm::Handle< edm::HepMCProduct > evtMC, edm::Handle< edm::SimTrackContainer > simTracks)
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
int ires[2]
std::vector< MuScleFitMuon > fillMuonCollection(const std::vector< T > &tracks)
Template function used to convert the muon collection to a vector of reco::LeafCandidate.
bool selTrackerMuon(const pat::Muon *aMuon)
Apply the Onia cuts to select trackerMuons.
bool muonID(const std::string &name) const
int charge() const final
electric charge
Definition: LeafCandidate.h:91
void selectSimulatedMuons(const edm::Event &event, const bool ifHepMC, edm::Handle< edm::HepMCProduct > evtMC, std::vector< std::pair< lorentzVector, lorentzVector > > &simPair, MuScleFitPlotter *plotter)
MuScleFitMuon mu1
Definition: GenMuonPair.h:54
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:499
#define nullptr
Int_t motherId
Definition: GenMuonPair.h:56
float userFloat(const std::string &key) const
Definition: PATObject.h:791
reco::Particle::LorentzVector lorentzVector
Definition: GenMuonPair.h:9
bool isTrackerMuon() const override
Definition: Muon.h:273
static const unsigned int motherPdgIdArray[6]
const edm::InputTag muonLabel_
static const double mMu2
reco::TrackRef innerTrack() const override
reference to Track reconstructed in the tracker only (reimplemented from reco::Muon) ...
Definition: Muon.h:73
void fillRec(std::vector< MuScleFitMuon > &muons)
Used when running on the root tree containing preselected muon pairs.
bool isGlobalMuon() const override
Definition: Muon.h:272
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
virtual int pdgId() const =0
PDG identifier.
simTrack
per collection params
void fillSim(edm::Handle< edm::SimTrackContainer > simTracks)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isMatched(TrackingRecHit const &hit)
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:245
bool isValid() const
Definition: HandleBase.h:74
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
T const * product() const
Definition: Handle.h:81
reco::TrackRef globalTrack() const override
reference to Track reconstructed in both tracked and muon detector (reimplemented from reco::Muon) ...
Definition: Muon.h:81
part
Definition: HCALResponse.h:20
GenMuonPair findGenMuFromRes(const reco::GenParticleCollection *genParticles)
void selectGenSimMuons(const edm::Event &event, std::vector< GenMuonPair > &genPair, std::vector< std::pair< lorentzVector, lorentzVector > > &simPair, MuScleFitPlotter *plotter)
const edm::InputTag simTracksCollectionName_
void fillGen(const reco::GenParticleCollection *genParticles, bool=false)
bool selGlobalMuon(const pat::Muon *aMuon)
Apply the Onia cuts to select globalMuons.
void takeSelectedMuonType(const T &muon, std::vector< reco::Track > &tracks)
Template function used to extract the selected muon type from the muon collection.
int numberOfValidMuonHits() const
Definition: HitPattern.h:833
const reco::Candidate * getStatus3Muon(const reco::Candidate *status3Muon)
virtual size_type numberOfDaughters() const =0
number of daughters
const std::vector< int > resfind_
Analysis-level muon class.
Definition: Muon.h:50
MuScleFitMuon mu2
Definition: GenMuonPair.h:55
Definition: event.py:1
const reco::Candidate * getStatus1Muon(const reco::Candidate *status3Muon)