CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/Validation/RecoMuon/src/GlobalMuonMatchAnalyzer.cc

Go to the documentation of this file.
00001 
00013 #include "Validation/RecoMuon/src/GlobalMuonMatchAnalyzer.h"
00014 
00015 // user include files
00016 #include "FWCore/Framework/interface/Frameworkfwd.h"
00017 #include "FWCore/Framework/interface/EDAnalyzer.h"
00018 #include "FWCore/Framework/interface/ESHandle.h"
00019 #include "FWCore/Framework/interface/Event.h"
00020 #include "FWCore/Framework/interface/MakerMacros.h"
00021 
00022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00023 #include "FWCore/Utilities/interface/InputTag.h"
00024 
00025 #include "SimTracker/Records/interface/TrackAssociatorRecord.h"
00026 #include "SimTracker/TrackAssociation/interface/TrackAssociatorBase.h"
00027 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00028 
00029 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00030 #include "DataFormats/Common/interface/View.h"
00031 #include "DataFormats/TrackReco/interface/Track.h"
00032 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00033 #include "DataFormats/MuonReco/interface/MuonTrackLinks.h"
00034 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00035 
00036 #include "DQMServices/Core/interface/DQMStore.h"
00037 #include "DQMServices/Core/interface/MonitorElement.h"
00038 #include "FWCore/ServiceRegistry/interface/Service.h"
00039 
00040 #include <TH2.h>
00041 
00042 
00043 GlobalMuonMatchAnalyzer::GlobalMuonMatchAnalyzer(const edm::ParameterSet& iConfig)
00044 
00045 {
00046   //now do what ever initialization is needed
00047   tkAssociatorName_ = iConfig.getUntrackedParameter<std::string>("tkAssociator");
00048   muAssociatorName_ = iConfig.getUntrackedParameter<std::string>("muAssociator");
00049 
00050   tpName_ = iConfig.getUntrackedParameter<edm::InputTag>("tpLabel");
00051   tkName_ = iConfig.getUntrackedParameter<edm::InputTag>("tkLabel");
00052   staName_ = iConfig.getUntrackedParameter<edm::InputTag>("muLabel");
00053   glbName_ = iConfig.getUntrackedParameter<edm::InputTag>("glbLabel");
00054 
00055   out = iConfig.getUntrackedParameter<std::string>("out");
00056   dbe_ = edm::Service<DQMStore>().operator->();
00057 }
00058 
00059 
00060 GlobalMuonMatchAnalyzer::~GlobalMuonMatchAnalyzer()
00061 {
00062  
00063    // do anything here that needs to be done at desctruction time
00064    // (e.g. close files, deallocate resources etc.)
00065 
00066 }
00067 
00068 
00069 //
00070 // member functions
00071 //
00072 
00073 // ------------ method called to for each event  ------------
00074 void
00075 GlobalMuonMatchAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00076 {
00077    using namespace edm;
00078    using namespace reco;
00079 
00080    Handle<TrackingParticleCollection> tpHandle;
00081    iEvent.getByLabel(tpName_,tpHandle);
00082    const TrackingParticleCollection tpColl = *(tpHandle.product());
00083 
00084    Handle<reco::MuonTrackLinksCollection> muHandle;
00085    iEvent.getByLabel(glbName_,muHandle);
00086    const reco::MuonTrackLinksCollection muColl = *(muHandle.product());
00087 
00088    Handle<View<Track> > staHandle;
00089    iEvent.getByLabel(staName_,staHandle);
00090    //   const reco::TrackCollection staColl = *(staHandle.product());
00091 
00092    Handle<View<Track> > glbHandle;
00093    iEvent.getByLabel(glbName_,glbHandle);
00094    //   const reco::TrackCollection glbColl = *(glbHandle.product());
00095 
00096    Handle<View<Track> > tkHandle;
00097    iEvent.getByLabel(tkName_,tkHandle);
00098    //   const reco::TrackCollection mtkColl = *(tkHandle.product());
00099 
00100    reco::RecoToSimCollection tkrecoToSimCollection;
00101    reco::SimToRecoCollection tksimToRecoCollection;
00102    tkrecoToSimCollection = tkAssociator_->associateRecoToSim(tkHandle,tpHandle,&iEvent);
00103    tksimToRecoCollection = tkAssociator_->associateSimToReco(tkHandle,tpHandle,&iEvent);
00104 
00105    reco::RecoToSimCollection starecoToSimCollection;
00106    reco::SimToRecoCollection stasimToRecoCollection;
00107    starecoToSimCollection = muAssociator_->associateRecoToSim(staHandle,tpHandle,&iEvent);
00108    stasimToRecoCollection = muAssociator_->associateSimToReco(staHandle,tpHandle,&iEvent);
00109 
00110    reco::RecoToSimCollection glbrecoToSimCollection;
00111    reco::SimToRecoCollection glbsimToRecoCollection;
00112    glbrecoToSimCollection = muAssociator_->associateRecoToSim(glbHandle,tpHandle,&iEvent);
00113    glbsimToRecoCollection = muAssociator_->associateSimToReco(glbHandle,tpHandle,&iEvent);
00114 
00115 
00116    for (TrackingParticleCollection::size_type i=0; i<tpColl.size(); ++i){
00117      TrackingParticleRef tp(tpHandle,i);
00118 
00119      std::vector<std::pair<RefToBase<Track>, double> > rvGlb;
00120      RefToBase<Track> rGlb;
00121      if(glbsimToRecoCollection.find(tp) != glbsimToRecoCollection.end()){
00122        rvGlb = glbsimToRecoCollection[tp];
00123        if(rvGlb.size() != 0) {
00124          rGlb = rvGlb.begin()->first;
00125        }
00126      }
00127 
00128      std::vector<std::pair<RefToBase<Track>, double> > rvSta;
00129      RefToBase<Track> rSta;
00130      if(stasimToRecoCollection.find(tp) != stasimToRecoCollection.end()){
00131        rvSta = stasimToRecoCollection[tp];
00132        if(rvSta.size() != 0) {
00133          rSta = rvSta.begin()->first;
00134        }
00135      }
00136 
00137      std::vector<std::pair<RefToBase<Track>, double> > rvTk;
00138      RefToBase<Track> rTk;
00139      if(tksimToRecoCollection.find(tp) != tksimToRecoCollection.end()){
00140        rvTk = tksimToRecoCollection[tp];
00141        if(rvTk.size() != 0) {
00142          rTk = rvTk.begin()->first;
00143        }
00144      }
00145      
00146      if( rvSta.size() != 0 && rvTk.size() != 0 ){
00147        //should have matched
00148        h_shouldMatch->Fill(rTk->eta(),rTk->pt());
00149      }
00150 
00151      for ( reco::MuonTrackLinksCollection::const_iterator links = muHandle->begin(); links != muHandle->end(); ++links ) {
00152        if( rGlb == RefToBase<Track>(links->globalTrack() ) ) {
00153          if( RefToBase<Track>(links->trackerTrack() ) == rTk && 
00154              RefToBase<Track>(links->standAloneTrack() ) == rSta ) {
00155            //goodMatch
00156            h_goodMatchSim->Fill(rGlb->eta(),rGlb->pt());
00157          } 
00158          if ( RefToBase<Track>(links->trackerTrack() ) == rTk &&
00159               RefToBase<Track>(links->standAloneTrack() ) != rSta ) {
00160            //tkOnlyMatch
00161            h_tkOnlySim->Fill(rGlb->eta(),rGlb->pt());
00162          } 
00163          if ( RefToBase<Track>(links->standAloneTrack() ) == rSta &&
00164               RefToBase<Track>(links->trackerTrack() ) != rTk ) {
00165            //staOnlyMatch
00166            h_staOnlySim->Fill(rGlb->eta(),rGlb->pt());
00167          }
00168        }
00169      }
00170 
00171    }
00172    
00174    
00175    for ( reco::MuonTrackLinksCollection::const_iterator links = muHandle->begin(); links != muHandle->end(); ++links ) {
00176      RefToBase<Track> glbRef = RefToBase<Track>(links->globalTrack() );
00177      RefToBase<Track> staRef = RefToBase<Track>(links->standAloneTrack() );
00178      RefToBase<Track> tkRef  = RefToBase<Track>(links->trackerTrack() );
00179      
00180      std::vector<std::pair<TrackingParticleRef, double> > tp1;
00181      TrackingParticleRef tp1r;
00182      if(glbrecoToSimCollection.find(glbRef) != glbrecoToSimCollection.end()){
00183        tp1 = glbrecoToSimCollection[glbRef];
00184        if(tp1.size() != 0) {
00185          tp1r = tp1.begin()->first;
00186        }
00187      }
00188      
00189      std::vector<std::pair<TrackingParticleRef, double> > tp2;
00190      TrackingParticleRef tp2r;
00191      if(starecoToSimCollection.find(staRef) != starecoToSimCollection.end()){
00192        tp2 = starecoToSimCollection[staRef];
00193        if(tp2.size() != 0) {
00194          tp2r = tp2.begin()->first;
00195        }
00196      }
00197      
00198      std::vector<std::pair<TrackingParticleRef, double> > tp3;
00199      TrackingParticleRef tp3r;
00200      if(tkrecoToSimCollection.find(tkRef) != tkrecoToSimCollection.end()){
00201        tp3 = tkrecoToSimCollection[tkRef];
00202        if(tp3.size() != 0) {
00203          tp3r = tp3.begin()->first;
00204        }
00205      }
00206      
00207      
00208      if(tp1.size() != 0) {
00209        //was reconstructed
00210        h_totReco->Fill(glbRef->eta(),glbRef->pt());
00211        if(tp2r == tp3r) { // && tp1r == tp3r) {
00212          //came from same TP
00213          h_goodMatch->Fill(glbRef->eta(),glbRef->pt());
00214        } else {
00215          //mis-match
00216          h_fakeMatch->Fill(glbRef->eta(),glbRef->pt());
00217        }
00218      }
00219      
00220    }   
00221 
00222    
00223 }
00224 
00225 
00226 // ------------ method called once each job just before starting event loop  ------------
00227 void 
00228 GlobalMuonMatchAnalyzer::beginJob()
00229 {
00230   dbe_->cd();
00231   std::string dirName="Matcher/";
00232   dbe_->setCurrentFolder("RecoMuonV/Matcher");
00233 
00234   h_shouldMatch = dbe_->book2D("h_shouldMatch","SIM associated to Tk and Sta",50,-2.5,2.5,100,0.,500.);
00235   h_goodMatchSim = dbe_->book2D("h_goodMatchSim","SIM associated to Glb Sta Tk",50,-2.5,2.5,100,0.,500.);
00236   h_tkOnlySim = dbe_->book2D("h_tkOnlySim","SIM associated to Glb Tk",50,-2.5,2.5,100,0.,500.);
00237   h_staOnlySim = dbe_->book2D("h_staOnlySim","SIM associated to Glb Sta",50,-2.5,2.5,100,0.,500.);
00238 
00239   h_totReco = dbe_->book2D("h_totReco","Total Glb Reconstructed",50,-2.5,2.5,100,0.,500.);
00240   h_goodMatch = dbe_->book2D("h_goodMatch","Sta and Tk from same SIM",50,-2.5,2.5,100, 0., 500.);
00241   h_fakeMatch = dbe_->book2D("h_fakeMatch","Sta and Tk not from same SIM",50,-2.5,2.5,100,0.,500.);
00242 
00243   h_effic = dbe_->book1D("h_effic","Efficiency vs #eta",50,-2.5,2.5);
00244   h_efficPt = dbe_->book1D("h_efficPt","Efficiency vs p_{T}",100,0.,100.);
00245 
00246   h_fake = dbe_->book1D("h_fake","Fake fraction vs #eta",50,-2.5,2.5);
00247   h_fakePt = dbe_->book1D("h_fakePt","Fake fraction vs p_{T}",100,0.,100.);
00248 }
00249 
00250 // ------------ method called once each job just after ending the event loop  ------------
00251 void 
00252 GlobalMuonMatchAnalyzer::endJob() {
00253   computeEfficiencyEta(h_effic,h_goodMatchSim,h_shouldMatch);
00254   computeEfficiencyPt(h_efficPt,h_goodMatchSim,h_shouldMatch);
00255 
00256   computeEfficiencyEta(h_fake,h_fakeMatch,h_totReco);
00257   computeEfficiencyPt(h_fakePt,h_fakeMatch,h_totReco);
00258 
00259   if( out.size() != 0 && dbe_ ) dbe_->save(out);
00260 }
00261 
00262 void GlobalMuonMatchAnalyzer::beginRun(const edm::Run&, const edm::EventSetup& setup)
00263 {
00264   // Tk Associator
00265   edm::ESHandle<TrackAssociatorBase> tkassociatorHandle;
00266   setup.get<TrackAssociatorRecord>().get(tkAssociatorName_,tkassociatorHandle);
00267   tkAssociator_ = tkassociatorHandle.product();
00268 
00269   // Mu Associator
00270   edm::ESHandle<TrackAssociatorBase> muassociatorHandle;
00271   setup.get<TrackAssociatorRecord>().get(muAssociatorName_,muassociatorHandle);
00272   muAssociator_ = muassociatorHandle.product();
00273 }
00274 
00275 
00276 
00277 void GlobalMuonMatchAnalyzer::computeEfficiencyEta(MonitorElement *effHist, MonitorElement *recoTH2, MonitorElement *simTH2){
00278   TH2F * h1 = recoTH2->getTH2F();
00279   TH1D* reco = h1->ProjectionX();
00280 
00281   TH2F * h2 =simTH2->getTH2F();
00282   TH1D* sim  = h2 ->ProjectionX();
00283 
00284     
00285   TH1F *hEff = (TH1F*) reco->Clone();  
00286   
00287   hEff->Divide(sim);
00288   
00289   hEff->SetName("tmp_"+TString(reco->GetName()));
00290   
00291   // Set the error accordingly to binomial statistics
00292   int nBinsEta = hEff->GetNbinsX();
00293   for(int bin = 1; bin <=  nBinsEta; bin++) {
00294     float nSimHit = sim->GetBinContent(bin);
00295     float eff = hEff->GetBinContent(bin);
00296     float error = 0;
00297     if(nSimHit != 0 && eff <= 1) {
00298       error = sqrt(eff*(1-eff)/nSimHit);
00299     }
00300     hEff->SetBinError(bin, error);
00301     effHist->setBinContent(bin,eff);
00302     effHist->setBinError(bin,error);
00303   }
00304   
00305 }
00306 
00307 void GlobalMuonMatchAnalyzer::computeEfficiencyPt(MonitorElement *effHist, MonitorElement *recoTH2, MonitorElement *simTH2){
00308   TH2F * h1 = recoTH2->getTH2F();
00309   TH1D* reco = h1->ProjectionY();
00310 
00311   TH2F * h2 = simTH2->getTH2F();
00312   TH1D* sim  = h2 ->ProjectionY();
00313 
00314     
00315   TH1F *hEff = (TH1F*) reco->Clone();  
00316   
00317   hEff->Divide(sim);
00318   
00319   hEff->SetName("tmp_"+TString(reco->GetName()));
00320   
00321   // Set the error accordingly to binomial statistics
00322   int nBinsPt = hEff->GetNbinsX();
00323   for(int bin = 1; bin <=  nBinsPt; bin++) {
00324     float nSimHit = sim->GetBinContent(bin);
00325     float eff = hEff->GetBinContent(bin);
00326     float error = 0;
00327     if(nSimHit != 0 && eff <= 1) {
00328       error = sqrt(eff*(1-eff)/nSimHit);
00329     }
00330     hEff->SetBinError(bin, error);
00331     effHist->setBinContent(bin,eff);
00332     effHist->setBinError(bin,error);
00333   }
00334   
00335 }
00336