CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/MuonAnalysis/MomentumScaleCalibration/plugins/MuScleFitMuonSelector.cc

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