CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/HeavyFlavorAnalysis/Skimming/src/Tau3MuReco.cc

Go to the documentation of this file.
00001 #include <TMath.h>
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 #include "HeavyFlavorAnalysis/Skimming/interface/Combinatorics.h"
00006 
00007 #include "HeavyFlavorAnalysis/Skimming/interface/Tau3MuReco.h"
00008 
00009 Tau3MuReco::Tau3MuReco(const edm::ParameterSet& iConfig):m_kMatchingDeltaR(iConfig.getParameter<double>("RecoAnalysisMatchingDeltaR")),
00010                                                          m_kMatchingPt(iConfig.getParameter<double>("RecoAnalysisMatchingPt")),
00011                                                          m_kTauMassCut(iConfig.getParameter<double>("RecoAnalysisTauMassCut")),
00012                                                          m_kTauMass(iConfig.getParameter<double>("RecoAnalysisTauMass")),
00013                                                          m_kMuonMass(iConfig.getParameter<double>("RecoAnalysisMuonMass")),
00014                                                          m_kMuonSource(iConfig.getParameter<edm::InputTag>("MuonSourceTag")),
00015                                                          m_kTrackSource(iConfig.getParameter<edm::InputTag>("TrackSourceTag"))
00016 {
00017 }
00018 
00019 Tau3MuReco::~Tau3MuReco()
00020 {
00021 
00022 }
00023 
00024 bool Tau3MuReco::doTau3MuReco(const edm::Event& iEvent, const edm::EventSetup& iSetup, reco::MuonCollection* muonCollection, reco::TrackCollection* trackCollection)
00025 {
00026     m_MuonCollection = muonCollection;
00027     m_TrackCollection = trackCollection;
00028     
00029     edm::Handle<reco::MuonCollection> muons;
00030     edm::Handle<reco::TrackCollection> tracks;
00031         
00032     reco::MuonCollection::const_iterator muon;
00033 
00034     iEvent.getByLabel(m_kMuonSource, muons);
00035     iEvent.getByLabel(m_kTrackSource, tracks);
00036         
00037     for( muon = muons->begin(); muon != muons->end(); ++muon )
00038     {
00039         m_TrackCollection->push_back(*(muon->track().get()));
00040         m_MuonCollection->push_back(*muon);
00041     }
00042 
00043     if( m_TrackCollection->size() > 3 )
00044     {
00045         //find the right three ones coming from tau
00046         if(findCorrectPairing()==false)
00047         {  //maybe implement something like in ==3 (throw away ....)
00048             LogDebug("Tau3MuReco") << "Could not find correct combination!" << std::endl;
00049             return false;
00050         }
00051         
00052         return true;
00053     }
00054 
00055     if( m_TrackCollection->size() == 3 )
00056     {
00057         if(fabs(getInvariantMass(m_TrackCollection, m_kMuonMass)-m_kTauMass)<=m_kTauMassCut)
00058             return true;
00059         else//throw away one muon which don't match
00060             removeIncorrectMuon();
00061         
00062     }
00063 
00064     if( m_TrackCollection->size() == 2)
00065     {
00066         //search the third track
00067 
00068         //get 3rd muon canidate from tracks
00069         if(find3rdTrack(iEvent, iSetup, *(tracks.product()))==false)
00070         {
00071             LogDebug("Tau3MuReco") << "A 3rd Track can not be assigned!" << std::endl;
00072             return false;
00073         }
00074         else
00075             return true;
00076     }
00077 
00078     // cannot use this event, because less than 2 muons have been found
00079     
00080     LogDebug("Tau3MuReco") << "Not enough (" << m_TrackCollection->size() << ") muons found! Event skipped!" << std::endl;
00081     
00082     return false;
00083 }
00084 
00085 //private
00086 bool Tau3MuReco::check4MuonTrack(const reco::Track& track)
00087 {
00088     reco::TrackCollection::const_iterator iter;
00089 
00090     for(iter = m_TrackCollection->begin(); iter!=m_TrackCollection->end(); iter++)
00091     {
00092         //check if the track has the right charge
00093         //and check if dR is smaller than fMatchingDeltaR
00094         if((*iter).charge() == track.charge()
00095            && getDeltaR(*iter,track)<m_kMatchingDeltaR
00096            && fabs(((*iter).pt())-(track.pt()))<=(track.pt()*m_kMatchingPt))
00097         {
00098             LogDebug("Tau3MuReco") << "Found muon track in Tracks with DeltaR: " << getDeltaR(*iter,track)
00099                                        << " Pt: " << track.pt()
00100                                        << " Muons Pt is: " << (*iter).pt() << std::endl;
00101             return true;
00102         }
00103     }
00104     return false;
00105 }
00106 
00107 
00108 bool Tau3MuReco::find3rdTrack(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::TrackCollection& Tracks)
00109 {
00110     //check size of TrackVector (should be two!)
00111     if(m_TrackCollection->size()!=2)
00112         return false;
00113     
00114     //more then two tracks should be in the event
00115     if(Tracks.size()<=2)
00116         return false;
00117     
00118     double SumDeltaR = 0;
00119     
00120     double MuonDeltaR = getDeltaR(m_TrackCollection->at(0),m_TrackCollection->at(1));
00121     
00122     //Loop overall tracks
00123 
00124     LogDebug("Tau3MuReco") << "Number of tracks found: " << Tracks.size() << std::endl;
00125 
00126     std::multimap<double, const reco::Track> TrackMultiMap;
00127 
00128     unsigned short muonCounter = 0;
00129     
00130     reco::TrackCollection::const_iterator track;
00131 
00132     for(track=Tracks.begin(); track!=Tracks.end(); track++)
00133     {
00134         if(check4MuonTrack(*track))
00135         {
00136             muonCounter++;
00137             continue;
00138         }
00139         
00140         SumDeltaR = MuonDeltaR;
00141 
00142         SumDeltaR += getDeltaR(m_TrackCollection->at(1), *track);
00143         SumDeltaR += getDeltaR(*track,m_TrackCollection->at(0));
00144 
00145         std::pair<double, const reco::Track> actTrack(SumDeltaR,*track);
00146 
00147         TrackMultiMap.insert(actTrack);
00148     }
00149 
00150     //two tracks should be clearly identified as muons by check4MuonTrack
00151     //else event is not useable
00152     if(muonCounter<2)
00153     {
00154         LogDebug("Tau3MuReco") << "Not enough muons (" << muonCounter << ") assigned to a track! Event skipped!" << std::endl;
00155         return false;
00156     }
00157 
00158     std::multimap<double, const reco::Track>::iterator it = TrackMultiMap.begin();
00159 
00160     if(it==TrackMultiMap.end())
00161     {
00162         LogDebug("Tau3MuReco") << "Not enough tracks (0) left! Event skipped!" << std::endl;
00163         return false;
00164     }
00165 
00166     //get 2mu+track with minimal DeltaR Sum (MultiMaps are sorted)
00167     m_TrackCollection->push_back((*it).second);
00168 
00169     //and check charge of this track
00170     //and check invariant mass of this combination
00171     //and make a vertex fit
00172 
00173     char Charge = m_TrackCollection->at(0).charge() * m_TrackCollection->at(1).charge();
00174 
00175     unsigned int  count = 0;
00176 
00177     //Charge > 0 means the two muons have same charge, so the third track has to have the opposit charge
00178     while((Charge > 0 && ((*it).second).charge()==(m_TrackCollection->at(0)).charge())
00179           || fabs(getInvariantMass(m_TrackCollection)-m_kTauMass) > m_kTauMassCut
00180         )
00181     {
00182 
00183         count++;
00184 
00185         LogDebug("Tau3MuReco") << "Track canidate: " << count << std::endl;
00186         
00187         if(Charge > 0 && ((*it).second).charge()!=(m_TrackCollection->at(0)).charge())
00188         LogDebug("Tau3MuReco") << "\tWrong charge!" << std::endl;
00189         LogDebug("Tau3MuReco") << "\tInvariant Mass deviation! " << fabs(getInvariantMass(m_TrackCollection)-m_kTauMass) << std::endl;
00190         LogDebug("Tau3MuReco") << "\tTrack Pt: " << (*it).second.pt() << std::endl;
00191         LogDebug("Tau3MuReco") << "\tDelta R: " << (*it).first << std::endl;
00192         LogDebug("Tau3MuReco") << "\tChi2: " << ((*it).second).normalizedChi2() << std::endl;
00193 
00194         ++it;
00195 
00196         //was not the best canidate
00197         m_TrackCollection->pop_back();
00198 
00199         if(it==TrackMultiMap.end())
00200             return false;
00201 
00202         //get next to minimal (Delta R Sum) track
00203         m_TrackCollection->push_back((*it).second);
00204     }
00205 
00206     LogDebug("Tau3MuReco") << "Found corresponding 3rd track: " << std::endl;
00207     LogDebug("Tau3MuReco") << "Track canidate: " << count << std::endl;
00208     LogDebug("Tau3MuReco") << "\tInvariant Mass deviation! " << fabs(getInvariantMass(m_TrackCollection)-m_kTauMass) << std::endl;
00209     LogDebug("Tau3MuReco") << "\tDelta R: " << (*it).first << std::endl;
00210     LogDebug("Tau3MuReco") << "\tNormChi2: " << ((*it).second).normalizedChi2() << std::endl;
00211 
00212     //choose this track, because it is the best canidate
00213     return true;
00214 }
00215 
00216 bool Tau3MuReco::findCorrectPairing()
00217 {
00218     Combinatorics myCombinatorics(m_TrackCollection->size(), 3);
00219 
00220     std::vector < std::vector<UInt_t> > CombinationVec = myCombinatorics.GetCombinations();
00221 
00222     std::vector< std::vector<UInt_t> >::iterator it = CombinationVec.begin();
00223 
00224     char Charge = 0;
00225 
00226     reco::TrackCollection tempTrackCollection;
00227     reco::MuonCollection tempMuonCollection;
00228 
00229     do
00230     {
00231         if(it==CombinationVec.end())
00232             return false;
00233 
00234         Charge = 0;
00235         
00236         tempMuonCollection.clear();
00237         tempTrackCollection.clear();
00238         
00239         for(UInt_t i=0; i< (*it).size(); i++)
00240         {
00241             Charge += m_TrackCollection->at((*it).at(i)).charge();
00242             tempTrackCollection.push_back(m_TrackCollection->at((*it).at(i)));
00243             tempMuonCollection.push_back(m_MuonCollection->at((*it).at(i)));
00244         }
00245 
00246         LogDebug("Tau3MuReco") << "Charge is: " << (int)Charge << " Have to be -1 or 1!!!" << std::endl;
00247         LogDebug("Tau3MuReco") << "Invariant mass is: " << fabs(getInvariantMass(&tempTrackCollection)-m_kTauMass) << " Have to be smaller than " << m_kTauMassCut << std::endl;
00248 
00249         it++;
00250     }
00251     while(abs(Charge)!=1 || fabs(getInvariantMass(&tempTrackCollection)-m_kTauMass)>m_kTauMassCut);
00252 
00253     *m_MuonCollection = tempMuonCollection;
00254     *m_TrackCollection = tempTrackCollection;
00255 
00256     return true;
00257 }
00258 
00259 bool Tau3MuReco::removeIncorrectMuon()
00260 {
00261     
00262     double deltaR12 = getDeltaR(m_TrackCollection->at(0),m_TrackCollection->at(1));
00263     double deltaR23 = getDeltaR(m_TrackCollection->at(1),m_TrackCollection->at(2));
00264     double deltaR31 = getDeltaR(m_TrackCollection->at(2),m_TrackCollection->at(0));
00265     
00266     //if DeltaR12 is the smallest, than the 3rd one seems to be wrong
00267     //if DeltaR23 is the smallest, than the 2nd one seems to be wrong
00268     //if DeltaR31 is the smallest, than the 1st one seems to be wrong
00269     
00270     unsigned char temp;
00271     double junk;
00272     
00273     deltaR12 < deltaR23 ? temp=3 : temp=1;
00274     deltaR12 < deltaR23 ? junk=deltaR12 : junk=deltaR23;
00275 
00276     if(deltaR31 < junk)
00277         temp=2;
00278     
00279     m_TrackCollection->erase(m_TrackCollection->begin()+temp-1);
00280 
00281     return true;
00282 }
00283 
00284 double Tau3MuReco::getInvariantMass(const reco::TrackCollection* tracks, const double MuonMass)
00285 {
00286     unsigned int numOfParticles = tracks->size();
00287     
00288     double SumPx = 0;
00289     double SumPy = 0;
00290     double SumPz = 0;
00291     
00292     double SumE = 0;
00293     
00294     for(unsigned int i=0; i<numOfParticles; i++)
00295     {
00296         SumPx += tracks->at(i).px();
00297         SumPy += tracks->at(i).py();
00298         SumPz += tracks->at(i).pz();
00299         
00300         SumE += sqrt(pow(tracks->at(i).p(),2)+pow(MuonMass,2));
00301     }
00302     
00303     double invmass = sqrt(pow(SumE,2)-pow(SumPx,2)-pow(SumPy,2)-pow(SumPz,2));
00304             
00305     return invmass;
00306 }
00307 
00308 double Tau3MuReco::getDeltaR(const reco::Track& track1, const reco::Track& track2)
00309 {
00310     double dEta = track1.eta() - track2.eta();
00311     double dPhi = track1.phi() - track2.phi();
00312     
00313     while(dPhi >= TMath::Pi())       dPhi -= (2.0*TMath::Pi());
00314     while(dPhi < (-1.0*TMath::Pi())) dPhi += (2.0*TMath::Pi());
00315     
00316     return sqrt(pow(dEta,2)+pow(dPhi,2));
00317 }