CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/DQMOffline/Muon/src/SegmentTrackAnalyzer.cc

Go to the documentation of this file.
00001 
00002 /*
00003  *  See header file for a description of this class.
00004  *
00005  *  $Date: 2009/12/22 17:43:38 $
00006  *  $Revision: 1.16 $
00007  *  \author G. Mila - INFN Torino
00008  */
00009 
00010 #include "DQMOffline/Muon/src/SegmentTrackAnalyzer.h"
00011 
00012 #include "DataFormats/Common/interface/Handle.h"
00013 #include "DataFormats/TrackReco/interface/Track.h"
00014 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00015 #include "DataFormats/MuonDetId/interface/MuonSubdetId.h"
00016 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00017 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00018 #include "RecoMuon/TransientTrackingRecHit/interface/MuonTransientTrackingRecHit.h"
00019 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00020 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
00021 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00022 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00023 
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025 
00026 #include <string>
00027 using namespace std;
00028 using namespace edm;
00029 
00030 
00031 
00032 SegmentTrackAnalyzer::SegmentTrackAnalyzer(const edm::ParameterSet& pSet, MuonServiceProxy *theService):MuonAnalyzerBase(theService) {
00033 
00034   parameters = pSet;
00035 
00036   const ParameterSet SegmentsTrackAssociatorParameters = parameters.getParameter<ParameterSet>("SegmentsTrackAssociatorParameters");
00037   theSegmentsAssociator = new SegmentsTrackAssociator(SegmentsTrackAssociatorParameters);
00038 
00039 }
00040 
00041 
00042 SegmentTrackAnalyzer::~SegmentTrackAnalyzer() { }
00043 
00044 
00045 void SegmentTrackAnalyzer::beginJob(DQMStore * dbe) {
00046 
00047 
00048   metname = "segmTrackAnalyzer";
00049   string trackCollection = parameters.getParameter<edm::InputTag>("MuTrackCollection").label() + parameters.getParameter<edm::InputTag>("MuTrackCollection").instance();
00050   LogTrace(metname)<<"[SegmentTrackAnalyzer] Parameters initialization";
00051   dbe->setCurrentFolder("Muons/SegmentTrackAnalyzer");
00052   
00053   // histograms initalization
00054   hitsNotUsed = dbe->book1D("HitsNotUsedForGlobalTracking_"+trackCollection, "recHits not used for GLB    ["+trackCollection+"]", 50, -0.5, 49.5);
00055   hitsNotUsedPercentual = dbe->book1D("HitsNotUsedForGlobalTrackingDvHitUsed_"+trackCollection, "(recHits_{notUsedForGLB}) / (recHits_{GLB})    ["+trackCollection+"]", 100, 0, 1.);
00056 
00057   TrackSegm = dbe->book2D("trackSegments_"+trackCollection, "Number of segments associated to the track    ["+trackCollection+"]", 3, 0.5, 3.5, 8, 0, 8);
00058   TrackSegm->setBinLabel(1,"DT+CSC",1);
00059   TrackSegm->setBinLabel(2,"DT",1);
00060   TrackSegm->setBinLabel(3,"CSC",1);
00061   
00062   hitStaProvenance = dbe->book1D("trackHitStaProvenance_"+trackCollection, "Number of recHits_{STAinTrack}    ["+trackCollection+"]", 7, 0.5, 7.5);
00063   hitStaProvenance->setBinLabel(1,"DT");
00064   hitStaProvenance->setBinLabel(2,"CSC");
00065   hitStaProvenance->setBinLabel(3,"RPC");
00066   hitStaProvenance->setBinLabel(4,"DT+CSC");
00067   hitStaProvenance->setBinLabel(5,"DT+RPC");
00068   hitStaProvenance->setBinLabel(6,"CSC+RPC");
00069   hitStaProvenance->setBinLabel(7,"DT+CSC+RPC");
00070 
00071 
00072   if(trackCollection!="standAloneMuons"){
00073     hitTkrProvenance = dbe->book1D("trackHitTkrProvenance_"+trackCollection, "Number of recHits_{TKinTrack}    ["+trackCollection+"]", 6, 0.5, 6.5);
00074     hitTkrProvenance->setBinLabel(1,"PixBarrel");
00075     hitTkrProvenance->setBinLabel(2,"PixEndCap");
00076     hitTkrProvenance->setBinLabel(3,"TIB");
00077     hitTkrProvenance->setBinLabel(4,"TID");
00078     hitTkrProvenance->setBinLabel(5,"TOB");
00079     hitTkrProvenance->setBinLabel(6,"TEC");
00080   }
00081 
00082   int etaBin = parameters.getParameter<int>("etaBin");
00083   double etaMin = parameters.getParameter<double>("etaMin");
00084   double etaMax = parameters.getParameter<double>("etaMax");
00085   trackHitPercentualVsEta = dbe->book2D("trackHitDivtrackSegmHitVsEta_"+trackCollection, "(recHits_{Track} / recHits_{associatedSegm}) vs #eta    [" +trackCollection+"]", etaBin, etaMin, etaMax, 20, 0, 1);
00086   dtTrackHitPercentualVsEta = dbe->book2D("dtTrackHitDivtrackSegmHitVsEta_"+trackCollection, "(recHits_{DTinTrack} / recHits_{associatedSegm}) vs #eta    [" +trackCollection+"]", etaBin, etaMin, etaMax, 20, 0, 1);
00087   cscTrackHitPercentualVsEta = dbe->book2D("cscTrackHitDivtrackSegmHitVsEta_"+trackCollection, "(recHits_{CSCinTrack} / recHits_{associatedSegm}) vs #eta    [" +trackCollection+"]", etaBin, etaMin, etaMax, 20, 0, 1);
00088 
00089   int phiBin = parameters.getParameter<int>("phiBin");
00090   double phiMin = parameters.getParameter<double>("phiMin");
00091   double phiMax = parameters.getParameter<double>("phiMax");
00092   trackHitPercentualVsPhi = dbe->book2D("trackHitDivtrackSegmHitVsPhi_"+trackCollection, "(recHits_{Track} / recHits_{associatedSegm}) vs #phi    [" +trackCollection+"]", phiBin, phiMin, phiMax, 20, 0, 1);
00093   trackHitPercentualVsPhi->setAxisTitle("rad",2);
00094   dtTrackHitPercentualVsPhi = dbe->book2D("dtTrackHitDivtrackSegmHitVsPhi_"+trackCollection, "(recHits_{DTinTrack} / recHits_{associatedSegm}) vs #phi    [" +trackCollection+"]", phiBin, phiMin, phiMax, 20, 0, 1);
00095   dtTrackHitPercentualVsPhi->setAxisTitle("rad",2);
00096   cscTrackHitPercentualVsPhi = dbe->book2D("cscTrackHitDivtrackSegmHitVsPhi_"+trackCollection, "(recHits_{CSCinTrack} / recHits_{associatedSegm}) vs #phi    [" +trackCollection+"]", phiBin, phiMin, phiMax, 20, 0, 1);
00097   cscTrackHitPercentualVsPhi->setAxisTitle("rad",2);
00098 
00099   int ptBin = parameters.getParameter<int>("ptBin");
00100   double ptMin = parameters.getParameter<double>("ptMin");
00101   double ptMax = parameters.getParameter<double>("ptMax");
00102   trackHitPercentualVsPt = dbe->book2D("trackHitDivtrackSegmHitVsPt_"+trackCollection, "(recHits_{Track} / recHits_{associatedSegm}) vs 1/p_{t}    [" +trackCollection+"]", ptBin, ptMin, ptMax, 20, 0, 1);
00103   trackHitPercentualVsPt->setAxisTitle("GeV",2);
00104   dtTrackHitPercentualVsPt = dbe->book2D("dtTrackHitDivtrackSegmHitVsPt_"+trackCollection, "(recHits_{DTinTrack} / recHits_{associatedSegm}) vs 1/p_{t}    [" +trackCollection+"]", ptBin, ptMin, ptMax, 20, 0, 1);
00105   dtTrackHitPercentualVsPt->setAxisTitle("GeV",2);
00106   cscTrackHitPercentualVsPt = dbe->book2D("cscTrackHitDivtrackSegmHitVsPt_"+trackCollection, "(recHits_{CSCinTrack} / recHits_{associatedSegm}) vs 1/p_{t}    [" +trackCollection+"]", ptBin, ptMin, ptMax, 20, 0, 1);
00107   cscTrackHitPercentualVsPt->setAxisTitle("GeV",2);
00108 
00109 }
00110 
00111 
00112 void SegmentTrackAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::Track& recoTrack){
00113 
00114   LogTrace(metname)<<"[SegmentTrackAnalyzer] Filling the histos";
00115   
00116   MuonTransientTrackingRecHit::MuonRecHitContainer segments = theSegmentsAssociator->associate(iEvent, iSetup, recoTrack );
00117  
00118   LogTrace(metname)<<"[SegmentTrackAnalyzer] # of segments associated to the track: "<<(segments).size();
00119 
00120   // hit counters
00121   int hitsFromDt=0;
00122   int hitsFromCsc=0;
00123   int hitsFromRpc=0;
00124   int hitsFromTk=0;
00125   int hitsFromTrack=0;
00126   int hitsFromSegmDt=0;
00127   int hitsFromSegmCsc=0;
00128   // segment counters
00129   int segmFromDt=0;
00130   int segmFromCsc=0;
00131 
00132   for (MuonTransientTrackingRecHit::MuonRecHitContainer::const_iterator segment=segments.begin();
00133        segment!=segments.end(); segment++) {
00134    
00135     DetId id = (*segment)->geographicalId();
00136     
00137     // hits from DT segments
00138     if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::DT ) {
00139       ++segmFromDt;
00140       const DTRecSegment4D *seg4D = dynamic_cast<const DTRecSegment4D*>((*segment)->hit());
00141       if((*seg4D).hasPhi())
00142         hitsFromSegmDt+=(*seg4D).phiSegment()->specificRecHits().size();
00143       if((*seg4D).hasZed())
00144         hitsFromSegmDt+=(*seg4D).zSegment()->specificRecHits().size();
00145       
00146     }
00147     
00148     // hits from CSC segments
00149     if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::CSC ) {
00150       hitsFromSegmCsc+=(*segment)->recHits().size();
00151       segmFromCsc++;
00152     }
00153 
00154   }
00155 
00156 
00157   // hits from track
00158   for(trackingRecHit_iterator recHit =  recoTrack.recHitsBegin(); recHit != recoTrack.recHitsEnd(); ++recHit){
00159 
00160     hitsFromTrack++;
00161      DetId id = (*recHit)->geographicalId();
00162      // hits from DT
00163      if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::DT ) 
00164        hitsFromDt++;   
00165      // hits from CSC
00166       if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::CSC ) 
00167        hitsFromCsc++;
00168      // hits from RPC
00169      if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::RPC ) 
00170        hitsFromRpc++;
00171      // hits from Tracker
00172      if (id.det() == DetId::Tracker){
00173        hitsFromTk++;
00174        if(id.subdetId() == PixelSubdetector::PixelBarrel )
00175          hitTkrProvenance->Fill(1);
00176        if(id.subdetId() == PixelSubdetector::PixelEndcap )
00177          hitTkrProvenance->Fill(2);
00178        if(id.subdetId() == SiStripDetId::TIB )
00179          hitTkrProvenance->Fill(3);
00180        if(id.subdetId() == SiStripDetId::TID )
00181          hitTkrProvenance->Fill(4);
00182        if(id.subdetId() == SiStripDetId::TOB )
00183          hitTkrProvenance->Fill(5);
00184        if(id.subdetId() == SiStripDetId::TEC )
00185          hitTkrProvenance->Fill(6);
00186      }
00187 
00188   }
00189 
00190   // fill the histos
00191   hitsNotUsed->Fill(hitsFromSegmDt+hitsFromSegmCsc+hitsFromRpc+hitsFromTk-hitsFromTrack);
00192   hitsNotUsedPercentual->Fill(double(hitsFromSegmDt+hitsFromSegmCsc+hitsFromRpc+hitsFromTk-hitsFromTrack)/hitsFromTrack);
00193 
00194   if(hitsFromDt!=0 && hitsFromCsc!=0)
00195     TrackSegm->Fill(1,segmFromDt+segmFromCsc);
00196   if(hitsFromDt!=0 && hitsFromCsc==0)
00197     TrackSegm->Fill(2,segmFromDt);
00198   if(hitsFromDt==0 && hitsFromCsc!=0)
00199     TrackSegm->Fill(3,segmFromCsc);
00200 
00201   if(hitsFromDt!=0 && hitsFromCsc==0 && hitsFromRpc==0) hitStaProvenance->Fill(1);
00202   if(hitsFromCsc!=0 && hitsFromDt==0 && hitsFromRpc==0) hitStaProvenance->Fill(2);
00203   if(hitsFromRpc!=0 && hitsFromDt==0 && hitsFromCsc==0) hitStaProvenance->Fill(3);
00204   if(hitsFromDt!=0 && hitsFromCsc!=0 && hitsFromRpc==0) hitStaProvenance->Fill(4);
00205   if(hitsFromDt!=0 && hitsFromRpc!=0 && hitsFromCsc==0) hitStaProvenance->Fill(5);
00206   if(hitsFromCsc!=0 && hitsFromRpc!=0 && hitsFromDt==0) hitStaProvenance->Fill(6);
00207   if(hitsFromDt!=0 && hitsFromCsc!=0 && hitsFromRpc!=0) hitStaProvenance->Fill(7);
00208 
00209   if(hitsFromSegmDt+hitsFromSegmCsc !=0){
00210     trackHitPercentualVsEta->Fill(recoTrack.eta(), double(hitsFromDt+hitsFromCsc)/(hitsFromSegmDt+hitsFromSegmCsc));
00211     trackHitPercentualVsPhi->Fill(recoTrack.phi(), double(hitsFromDt+hitsFromCsc)/(hitsFromSegmDt+hitsFromSegmCsc));
00212     trackHitPercentualVsPt->Fill(recoTrack.pt(), double(hitsFromDt+hitsFromCsc)/(hitsFromSegmDt+hitsFromSegmCsc));
00213   }
00214 
00215   if(hitsFromSegmDt!=0){
00216     dtTrackHitPercentualVsEta->Fill(recoTrack.eta(), double(hitsFromDt)/hitsFromSegmDt);
00217     dtTrackHitPercentualVsPhi->Fill(recoTrack.phi(), double(hitsFromDt)/hitsFromSegmDt);
00218     dtTrackHitPercentualVsPt->Fill(recoTrack.pt(), double(hitsFromDt)/hitsFromSegmDt);
00219   }
00220 
00221   if(hitsFromSegmCsc!=0){
00222     cscTrackHitPercentualVsEta->Fill(recoTrack.eta(), double(hitsFromCsc)/hitsFromSegmCsc);
00223     cscTrackHitPercentualVsPhi->Fill(recoTrack.phi(), double(hitsFromCsc)/hitsFromSegmCsc);
00224     cscTrackHitPercentualVsPt->Fill(recoTrack.pt(), double(hitsFromCsc)/hitsFromSegmCsc);
00225   }
00226 
00227 } 
00228