CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/MuonAnalysis/MomentumScaleCalibration/src/MuScleFitMuonSelector.cc

Go to the documentation of this file.
00001 #include "MuonAnalysis/MomentumScaleCalibration/interface/MuScleFitMuonSelector.h"
00002 #include "DataFormats/MuonReco/interface/CaloMuon.h"
00003 
00004 const double MuScleFitMuonSelector::mMu2 = 0.011163612;
00005 const unsigned int MuScleFitMuonSelector::motherPdgIdArray[] = {23, 100553, 100553, 553, 100443, 443};
00006 
00007 
00008 bool MuScleFitMuonSelector::selGlobalMuon(const pat::Muon* aMuon)
00009 {
00010   reco::TrackRef iTrack = aMuon->innerTrack();
00011   const reco::HitPattern& p = iTrack->hitPattern();
00012 
00013   reco::TrackRef gTrack = aMuon->globalTrack();
00014   const reco::HitPattern& q = gTrack->hitPattern();
00015 
00016   return (//isMuonInAccept(aMuon) &&// no acceptance cuts!
00017     iTrack->found() > 11 &&
00018     gTrack->chi2()/gTrack->ndof() < 20.0 &&
00019     q.numberOfValidMuonHits() > 0 &&
00020     iTrack->chi2()/iTrack->ndof() < 4.0 &&
00021     aMuon->muonID("TrackerMuonArbitrated") &&
00022     aMuon->muonID("TMLastStationAngTight") &&
00023     p.pixelLayersWithMeasurement() > 1 &&
00024     fabs(iTrack->dxy()) < 3.0 &&  //should be done w.r.t. PV!
00025     fabs(iTrack->dz()) < 15.0 //should be done w.r.t. PV!
00026   );
00027 }
00028 
00029 bool MuScleFitMuonSelector::selTrackerMuon(const pat::Muon* aMuon)
00030 {
00031   reco::TrackRef iTrack = aMuon->innerTrack();
00032   const reco::HitPattern& p = iTrack->hitPattern();
00033 
00034   return (//isMuonInAccept(aMuon) // no acceptance cuts!
00035     iTrack->found() > 11 &&
00036     iTrack->chi2()/iTrack->ndof() < 4.0 &&
00037     aMuon->muonID("TrackerMuonArbitrated") &&
00038     aMuon->muonID("TMLastStationAngTight") &&
00039     p.pixelLayersWithMeasurement() > 1 &&
00040     fabs(iTrack->dxy()) < 3.0 && //should be done w.r.t. PV!
00041     fabs(iTrack->dz()) < 15.0 //should be done w.r.t. PV!
00042   );
00043 }
00044 
00045 // Note that at this level we save all the information. Events for which no suitable muon pair is found
00046 // are removed from the tree (together with the gen and sim information) at a later stage.
00047 // It would be better to remove them directly at this point.
00048 void MuScleFitMuonSelector::selectMuons(const edm::Event & event, std::vector<reco::LeafCandidate> & muons,
00049                                         std::vector<GenMuonPair> & genPair,
00050                                         std::vector<std::pair<lorentzVector,lorentzVector> > & simPair,
00051                                         MuScleFitPlotter * plotter)
00052 {
00053   edm::Handle<pat::CompositeCandidateCollection > collAll;
00054   try {event.getByLabel("onia2MuMuPatTrkTrk",collAll);}
00055   catch (...) {std::cout << "J/psi not present in event!" << std::endl;}
00056   std::vector<const pat::Muon*> collMuSel;
00057 
00058   //================onia cuts===========================/
00059 
00060   if( muonType_ <= -1 && PATmuons_) {
00061     std::vector<const pat::CompositeCandidate*> collSelGG;
00062     std::vector<const pat::CompositeCandidate*> collSelGT;
00063     std::vector<const pat::CompositeCandidate*> collSelTT;
00064     if (collAll.isValid()) {
00065 
00066       for(std::vector<pat::CompositeCandidate>::const_iterator it=collAll->begin();
00067           it!=collAll->end();++it) {
00068       
00069         const pat::CompositeCandidate* cand = &(*it);   
00070         // cout << "Now checking candidate of type " << theJpsiCat << " with pt = " << cand->pt() << endl;
00071         const pat::Muon* muon1 = dynamic_cast<const pat::Muon*>(cand->daughter("muon1"));
00072         const pat::Muon* muon2 = dynamic_cast<const pat::Muon*>(cand->daughter("muon2"));
00073       
00074         if((muon1->charge() * muon2->charge())>0)
00075           continue;
00076         // global + global?
00077         if (muon1->isGlobalMuon() && muon2->isGlobalMuon() ) {
00078           if (selGlobalMuon(muon1) && selGlobalMuon(muon2) && cand->userFloat("vProb") > 0.001 ) {
00079             collSelGG.push_back(cand);
00080             continue;
00081           }
00082         }
00083         // global + tracker? (x2)    
00084         if (muon1->isGlobalMuon() && muon2->isTrackerMuon() ) {
00085           if (selGlobalMuon(muon1) &&  selTrackerMuon(muon2) && cand->userFloat("vProb") > 0.001 ) {
00086             collSelGT.push_back(cand);
00087             continue;
00088           }
00089         }
00090         if (muon2->isGlobalMuon() && muon1->isTrackerMuon() ) {
00091           if (selGlobalMuon(muon2) && selTrackerMuon(muon1) && cand->userFloat("vProb") > 0.001) {
00092             collSelGT.push_back(cand);
00093             continue;
00094           }
00095         }
00096         // tracker + tracker?  
00097         if (muon1->isTrackerMuon() && muon2->isTrackerMuon() ) {
00098           if (selTrackerMuon(muon1) && selTrackerMuon(muon2) && cand->userFloat("vProb") > 0.001) {
00099             collSelTT.push_back(cand);
00100             continue;
00101           }
00102         }
00103       }
00104     }
00105     // Split them in independent collections if using muonType_ == -2, -3 or -4. Take them all if muonType_ == -1.
00106     std::vector<reco::Track> tracks;
00107     if(collSelGG.size()){
00108       //CHECK THAT THEY ARE ORDERED BY PT !!!!!!!!!!!!!!!!!!!!!!!
00109       const pat::Muon* muon1 = dynamic_cast<const pat::Muon*>(collSelGG[0]->daughter("muon1"));
00110       const pat::Muon* muon2 = dynamic_cast<const pat::Muon*>(collSelGG[0]->daughter("muon2"));
00111       if( muonType_ == -1 || muonType_ == -2 ) {
00112         tracks.push_back(*(muon1->innerTrack()));
00113         tracks.push_back(*(muon2->innerTrack()));  
00114         collMuSel.push_back(muon1);
00115         collMuSel.push_back(muon2);
00116       }
00117     }
00118     else if(!collSelGG.size() && collSelGT.size()){
00119       //CHECK THAT THEY ARE ORDERED BY PT !!!!!!!!!!!!!!!!!!!!!!!
00120       const pat::Muon* muon1 = dynamic_cast<const pat::Muon*>(collSelGT[0]->daughter("muon1"));
00121       const pat::Muon* muon2 = dynamic_cast<const pat::Muon*>(collSelGT[0]->daughter("muon2"));
00122       if( muonType_ == -1 || muonType_ == -3 ) {
00123         tracks.push_back(*(muon1->innerTrack()));
00124         tracks.push_back(*(muon2->innerTrack()));   
00125         collMuSel.push_back(muon1);
00126         collMuSel.push_back(muon2);
00127      }
00128     }
00129     else if(!collSelGG.size() && !collSelGT.size() && collSelTT.size()){
00130       //CHECK THAT THEY ARE ORDERED BY PT !!!!!!!!!!!!!!!!!!!!!!!
00131       const pat::Muon* muon1 = dynamic_cast<const pat::Muon*>(collSelTT[0]->daughter("muon1"));
00132       const pat::Muon* muon2 = dynamic_cast<const pat::Muon*>(collSelTT[0]->daughter("muon2"));
00133       if( muonType_ == -1 || muonType_ == -4 ) {
00134         tracks.push_back(*(muon1->innerTrack()));
00135         tracks.push_back(*(muon2->innerTrack()));   
00136         collMuSel.push_back(muon1);
00137         collMuSel.push_back(muon2);
00138       }
00139     }
00140     if (tracks.size() != 2 && tracks.size() != 0){
00141       std::cout<<"ERROR strange number of muons selected by onia cuts!"<<std::endl;
00142       abort();
00143     }
00144     muons = fillMuonCollection(tracks); 
00145   }
00146   else if( (muonType_<4 && muonType_>0) || muonType_>=10 ) { // Muons (glb,sta,trk)
00147     std::vector<reco::Track> tracks;
00148     if( PATmuons_ == true ) {
00149       edm::Handle<pat::MuonCollection> allMuons;
00150       event.getByLabel( muonLabel_, allMuons );
00151       for( std::vector<pat::Muon>::const_iterator muon = allMuons->begin(); muon != allMuons->end(); ++muon ) {
00152         //std::cout<<"pat muon is global "<<muon->isGlobalMuon()<<std::endl;
00153         takeSelectedMuonType(muon, tracks);
00154       }
00155       muons = fillMuonCollection(tracks);
00156     }
00157     else {
00158       edm::Handle<reco::MuonCollection> allMuons;
00159       event.getByLabel (muonLabel_, allMuons);
00160       if( muonType_ == 0 ) {
00161         // Take directly the muon
00162         muons = fillMuonCollection(*allMuons);
00163       }
00164       else {
00165         for( std::vector<reco::Muon>::const_iterator muon = allMuons->begin(); muon != allMuons->end(); ++muon ) {
00166           takeSelectedMuonType(muon, tracks);
00167         }
00168         muons = fillMuonCollection(tracks);
00169       }
00170     }
00171   }
00172   else if(muonType_==4){  //CaloMuons
00173     edm::Handle<reco::CaloMuonCollection> caloMuons;
00174     event.getByLabel (muonLabel_, caloMuons);
00175     std::vector<reco::Track> tracks;
00176     for (std::vector<reco::CaloMuon>::const_iterator muon = caloMuons->begin(); muon != caloMuons->end(); ++muon){
00177       tracks.push_back(*(muon->track()));
00178     }
00179     muons = fillMuonCollection(tracks);
00180   }
00181 
00182   else if (muonType_==5) { // Inner tracker tracks
00183     edm::Handle<reco::TrackCollection> tracks;
00184     event.getByLabel (muonLabel_, tracks);
00185     muons = fillMuonCollection(*tracks);
00186   }
00187   plotter->fillRec(muons);
00188 
00189   // Generation and simulation part
00190   if( speedup_ == false )
00191   {
00192     if( PATmuons_ ) {
00193       selectGeneratedMuons(collAll, collMuSel, genPair, plotter);
00194     }
00195     else {
00196       selectGenSimMuons(event, genPair, simPair, plotter);
00197     }
00198   }
00199 }
00200 
00202 void MuScleFitMuonSelector::selectGeneratedMuons(const edm::Handle<pat::CompositeCandidateCollection> & collAll,
00203                                                  const std::vector<const pat::Muon*> & collMuSel,
00204                                                  std::vector<GenMuonPair> & genPair,
00205                                                  MuScleFitPlotter * plotter)
00206 {
00207   reco::GenParticleCollection* genPatParticles = new reco::GenParticleCollection;
00208 
00209   //explicitly for JPsi but can be adapted!!!!!
00210   for(std::vector<pat::CompositeCandidate>::const_iterator it=collAll->begin();
00211       it!=collAll->end();++it) {
00212     reco::GenParticleRef genJpsi = it->genParticleRef();
00213     bool isMatched = (genJpsi.isAvailable() && genJpsi->pdgId() == 443);  
00214     if (isMatched){
00215       genPatParticles->push_back(*(const_cast<reco::GenParticle*>(genJpsi.get())));
00216     }
00217   }
00218 
00219   if(collMuSel.size()==2) {
00220     reco::GenParticleRef genMu1 = collMuSel[0]->genParticleRef();
00221     reco::GenParticleRef genMu2 = collMuSel[1]->genParticleRef();
00222     bool isMuMatched = (genMu1.isAvailable() && genMu2.isAvailable() && 
00223                         genMu1->pdgId()*genMu2->pdgId() == -169);
00224     if (isMuMatched) {
00225       genPatParticles->push_back(*(const_cast<reco::GenParticle*>(genMu1.get())));
00226       genPatParticles->push_back(*(const_cast<reco::GenParticle*>(genMu2.get())));
00227 
00228       unsigned int motherId = 0;
00229       if( genMu1->mother() != 0 ) {
00230          motherId = genMu1->mother()->pdgId();
00231       }
00232       if(genMu1->pdgId()==13)
00233         genPair.push_back(GenMuonPair(genMu1.get()->p4(), genMu2.get()->p4(), motherId));
00234       else
00235         // genPair.push_back(std::make_pair(genMu2.get()->p4(),genMu1.get()->p4()) );
00236         genPair.push_back(GenMuonPair(genMu2.get()->p4(), genMu1.get()->p4(), motherId));
00237 
00238       plotter->fillGen1(const_cast <reco::GenParticleCollection*> (genPatParticles), true);
00239 
00240       if (debug_>0) std::cout << "Found genParticles in PAT" << std::endl;
00241     }
00242     else {
00243       std::cout << "No recomuon selected so no access to generated info"<<std::endl;
00244       // Fill it in any case, otherwise it will not be in sync with the event number
00245       // genPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );    
00246       genPair.push_back( GenMuonPair(lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.), 0 ) );    
00247     }
00248   }
00249   else{
00250     std::cout << "No recomuon selected so no access to generated info"<<std::endl;
00251     // Fill it in any case, otherwise it will not be in sync with the event number
00252     // genPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
00253     genPair.push_back( GenMuonPair(lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.), 0 ) );    
00254   }
00255   if(debug_>0) {
00256     std::cout << "genParticles:" << std::endl;
00257     // std::cout << genPair.back().first << " , " << genPair.back().second << std::endl;
00258     std::cout << genPair.back().mu1 << " , " << genPair.back().mu2 << std::endl;
00259   }
00260 }
00261 
00262 void MuScleFitMuonSelector::selectGenSimMuons(const edm::Event & event,
00263                                               std::vector<GenMuonPair> & genPair,
00264                                               std::vector<std::pair<lorentzVector,lorentzVector> > & simPair,
00265                                               MuScleFitPlotter * plotter)
00266 {  
00267   // Find and store in histograms the generated and simulated resonance and muons
00268   // ----------------------------------------------------------------------------
00269   edm::Handle<edm::HepMCProduct> evtMC;
00270   edm::Handle<reco::GenParticleCollection> genParticles;
00271 
00272   // Fill gen information only in the first loop
00273   bool ifHepMC=false;
00274 
00275   event.getByLabel( genParticlesName_, evtMC );
00276   event.getByLabel( genParticlesName_, genParticles );
00277   if( evtMC.isValid() ) {
00278     genPair.push_back( findGenMuFromRes(evtMC.product()) );
00279     plotter->fillGen2(evtMC.product(), sherpa_);
00280     ifHepMC = true;
00281     if (debug_>0) std::cout << "Found hepMC" << std::endl;
00282   }
00283   else if( genParticles.isValid() ) {
00284     genPair.push_back( findGenMuFromRes(genParticles.product()) );
00285     plotter->fillGen1(genParticles.product());
00286     if (debug_>0) std::cout << "Found genParticles" << std::endl;
00287   }
00288   else {
00289     std::cout<<"ERROR "<<"non generation info and speedup true!!!!!!!!!!!!"<<std::endl;
00290     // Fill it in any case, otherwise it will not be in sync with the event number
00291     // genPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
00292     genPair.push_back( GenMuonPair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.), 0 ) );
00293   }
00294   if(debug_>0) {
00295     std::cout << "genParticles:" << std::endl;
00296     std::cout << genPair.back().mu1 << " , " << genPair.back().mu2 << std::endl;
00297   }
00298   selectSimulatedMuons(event, ifHepMC, evtMC, simPair, plotter);
00299 }
00300 
00301 void MuScleFitMuonSelector::selectSimulatedMuons(const edm::Event & event,
00302                                                  const bool ifHepMC, edm::Handle<edm::HepMCProduct> evtMC,
00303                                                  std::vector<std::pair<lorentzVector,lorentzVector> > & simPair,
00304                                                  MuScleFitPlotter * plotter)
00305 {
00306   edm::Handle<edm::SimTrackContainer> simTracks;
00307   bool simTracksFound = false;
00308   event.getByLabel(simTracksCollectionName_, simTracks);
00309   if( simTracks.isValid() ) {
00310     plotter->fillSim(simTracks);
00311     if(ifHepMC) {
00312       simPair.push_back( findSimMuFromRes(evtMC, simTracks) );
00313       simTracksFound = true;
00314       plotter->fillGenSim(evtMC,simTracks);
00315     }
00316   }
00317   else {
00318     std::cout << "SimTracks not existent" << std::endl;
00319   }
00320   if( !simTracksFound ) {
00321     simPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
00322   }
00323 }
00324 
00325 GenMuonPair MuScleFitMuonSelector::findGenMuFromRes( const edm::HepMCProduct* evtMC )
00326 {
00327   const HepMC::GenEvent* Evt = evtMC->GetEvent();
00328   GenMuonPair muFromRes;
00329   //Loop on generated particles
00330   for (HepMC::GenEvent::particle_const_iterator part=Evt->particles_begin();
00331        part!=Evt->particles_end(); part++) {
00332     if (fabs((*part)->pdg_id())==13 && (*part)->status()==1) {
00333       bool fromRes = false;
00334       unsigned int motherPdgId = 0;
00335       for (HepMC::GenVertex::particle_iterator mother = (*part)->production_vertex()->particles_begin(HepMC::ancestors);
00336            mother != (*part)->production_vertex()->particles_end(HepMC::ancestors); ++mother) {
00337         motherPdgId = (*mother)->pdg_id();
00338 
00339         // For sherpa the resonance is not saved. The muons from the resonance can be identified
00340         // by having as mother a muon of status 3.
00341         if( sherpa_ ) {
00342           if( motherPdgId == 13 && (*mother)->status() == 3 ) fromRes = true;
00343         }
00344         else {
00345           for( int ires = 0; ires < 6; ++ires ) {
00346             if( motherPdgId == motherPdgIdArray[ires] && resfind_[ires] ) fromRes = true;
00347           }
00348         }
00349       }
00350       if(fromRes){
00351         if((*part)->pdg_id()==13) {
00352           //   muFromRes.first = (*part)->momentum();
00353           muFromRes.mu1 = (lorentzVector((*part)->momentum().px(),(*part)->momentum().py(),
00354                                          (*part)->momentum().pz(),(*part)->momentum().e()));
00355         }
00356         else {
00357           muFromRes.mu2 = (lorentzVector((*part)->momentum().px(),(*part)->momentum().py(),
00358                                          (*part)->momentum().pz(),(*part)->momentum().e()));
00359         }
00360         muFromRes.motherId = motherPdgId;
00361       }
00362     }
00363   }
00364   return muFromRes;
00365 }
00366 
00367 // std::pair<lorentzVector, lorentzVector> MuScleFitMuonSelector::findGenMuFromRes( const reco::GenParticleCollection* genParticles)
00368 GenMuonPair MuScleFitMuonSelector::findGenMuFromRes( const reco::GenParticleCollection* genParticles)
00369 {
00370   // std::pair<lorentzVector,lorentzVector> muFromRes;
00371   GenMuonPair muFromRes;
00372 
00373   //Loop on generated particles
00374   if( debug_>0 ) std::cout << "Starting loop on " << genParticles->size() << " genParticles" << std::endl;
00375   for( reco::GenParticleCollection::const_iterator part=genParticles->begin(); part!=genParticles->end(); ++part ) {
00376     if (fabs(part->pdgId())==13 && part->status()==1) {
00377       bool fromRes = false;
00378       unsigned int motherPdgId = part->mother()->pdgId();
00379       if( debug_>0 ) {
00380         std::cout << "Found a muon with mother: " << motherPdgId << std::endl;
00381       }
00382       for( int ires = 0; ires < 6; ++ires ) {
00383         if( motherPdgId == motherPdgIdArray[ires] && resfind_[ires] ) fromRes = true;
00384       }
00385       if(fromRes){
00386         if(part->pdgId()==13) {
00387           muFromRes.mu1 = part->p4();
00388           if( debug_>0 ) std::cout << "Found a genMuon + : " << muFromRes.mu1 << std::endl;
00389           //      muFromRes.first = (lorentzVector(part->p4().px(),part->p4().py(),
00390           //                                       part->p4().pz(),part->p4().e()));
00391         }
00392         else {
00393           muFromRes.mu2 = part->p4();
00394           if( debug_>0 ) std::cout << "Found a genMuon - : " << muFromRes.mu2 << std::endl;
00395           //      muFromRes.second = (lorentzVector(part->p4().px(),part->p4().py(),
00396           //                                        part->p4().pz(),part->p4().e()));
00397         }
00398       }
00399     }
00400   }
00401   return muFromRes;
00402 }
00403 
00404 std::pair<lorentzVector, lorentzVector> MuScleFitMuonSelector::findSimMuFromRes( const edm::Handle<edm::HepMCProduct> & evtMC,
00405                                                                                  const edm::Handle<edm::SimTrackContainer> & simTracks )
00406 {
00407   //Loop on simulated tracks
00408   std::pair<lorentzVector, lorentzVector> simMuFromRes;
00409   for( edm::SimTrackContainer::const_iterator simTrack=simTracks->begin(); simTrack!=simTracks->end(); ++simTrack ) {
00410     //Chose muons
00411     if (fabs((*simTrack).type())==13) {
00412       //If tracks from IP than find mother
00413       if ((*simTrack).genpartIndex()>0) {
00414         HepMC::GenParticle* gp = evtMC->GetEvent()->barcode_to_particle ((*simTrack).genpartIndex());
00415         if( gp != 0 ) {
00416 
00417           for (HepMC::GenVertex::particle_iterator mother = gp->production_vertex()->particles_begin(HepMC::ancestors);
00418                mother!=gp->production_vertex()->particles_end(HepMC::ancestors); ++mother) {
00419 
00420             bool fromRes = false;
00421             unsigned int motherPdgId = (*mother)->pdg_id();
00422             for( int ires = 0; ires < 6; ++ires ) {
00423               if( motherPdgId == motherPdgIdArray[ires] && resfind_[ires] ) fromRes = true;
00424             }
00425             if( fromRes ) {
00426               if(gp->pdg_id() == 13)
00427                 simMuFromRes.first = lorentzVector(simTrack->momentum().px(),simTrack->momentum().py(),
00428                                                    simTrack->momentum().pz(),simTrack->momentum().e());
00429               else
00430                 simMuFromRes.second = lorentzVector(simTrack->momentum().px(),simTrack->momentum().py(),
00431                                                     simTrack->momentum().pz(),simTrack->momentum().e());
00432             }
00433           }
00434         }
00435         // else LogDebug("MuScleFitUtils") << "WARNING: no matching genParticle found for simTrack" << std::endl;
00436       }
00437     }
00438   }
00439   return simMuFromRes;
00440 }