CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/Alignment/HIPAlignmentAlgorithm/src/LhcTrackAnalyzer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    LhcTrackAnalyzer
00004 // Class:      LhcTrackAnalyzer
00005 // 
00016 //
00017 
00018 // updated to 25/2/2009 5.30 pm
00019 
00020 //
00021 //
00022 
00023 // system include files
00024 #include <memory>
00025 
00026 
00027 // user include files
00028 #include "Alignment/HIPAlignmentAlgorithm/interface/LhcTrackAnalyzer.h"
00029 #include "FWCore/Framework/interface/Frameworkfwd.h"
00030 #include "FWCore/Framework/interface/EDAnalyzer.h"
00031 
00032 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00033 #include "FWCore/Framework/interface/Event.h"
00034 #include <FWCore/Framework/interface/ESHandle.h>
00035 #include "FWCore/Framework/interface/MakerMacros.h"
00036 
00037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00038 
00039 
00040 #include "TH1F.h"
00041 #include "TH2F.h"
00042 #include "TFile.h"
00043 #include "TROOT.h"
00044 #include "TChain.h"
00045 #include "TNtuple.h"
00046 #include "FWCore/ServiceRegistry/interface/Service.h"
00047 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00048 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00049 #include <Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h>
00050 #include <DataFormats/GeometrySurface/interface/Surface.h>
00051 #include <DataFormats/GeometrySurface/interface/GloballyPositioned.h>
00052 #include <Geometry/CommonDetUnit/interface/GeomDet.h>
00053 #include "MagneticField/Engine/interface/MagneticField.h" 
00054 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h" 
00055 
00056 #include "DataFormats/TrackReco/interface/Track.h"
00057 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00058 #include "DataFormats/VertexReco/interface/Vertex.h"
00059 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00060 
00061 // Constructor
00062 
00063 LhcTrackAnalyzer::LhcTrackAnalyzer(const edm::ParameterSet& iConfig)
00064 
00065 {
00066   //now do what ever initialization is needed
00067   debug_    = iConfig.getParameter<bool>       ("Debug");  
00068   TrackCollectionTag_      = iConfig.getParameter<edm::InputTag>("TrackCollectionTag");  
00069   PVtxCollectionTag_      = iConfig.getParameter<edm::InputTag>("PVtxCollectionTag");  
00070   filename_ = iConfig.getParameter<std::string>("OutputFileName");
00071 }
00072    
00073 // Destructor
00074 LhcTrackAnalyzer::~LhcTrackAnalyzer()
00075 {
00076  
00077    // do anything here that needs to be done at desctruction time
00078    // (e.g. close files, deallocate resources etc.)
00079 
00080 }
00081 
00082 
00083 //
00084 // member functions
00085 //
00086 
00087 // ------------ method called to for each event  ------------
00088 void
00089 LhcTrackAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00090 {
00091   using namespace edm;
00092   using namespace reco;
00093   using namespace std;
00094 
00095   //=======================================================
00096   // Initialize Root-tuple variables
00097   //=======================================================
00098 
00099   SetVarToZero();
00100   
00101   //=======================================================
00102   // Retrieve the Track information
00103   //=======================================================
00104   
00105   Handle< TrackCollection>  trackCollectionHandle;
00106   iEvent.getByLabel(TrackCollectionTag_, trackCollectionHandle);
00107   Handle<VertexCollection>  vertexCollectionHandle;
00108   iEvent.getByLabel(PVtxCollectionTag_, vertexCollectionHandle);
00109   for(VertexCollection::const_iterator vtx = vertexCollectionHandle->begin();vtx!=vertexCollectionHandle->end(); ++vtx)
00110     {
00111       if(vtx==vertexCollectionHandle->begin()){
00112         if(vtx->isFake())goodvtx_=false;
00113         else goodvtx_=true;
00114       }
00115       else break;
00116     }
00117 
00118 
00119 
00120   goodbx_=true;
00121   // int bx = iEvent.bunchCrossing();
00122   //if (bx==51 || bx==2724) goodbx_=true;
00123   
00124 
00125  
00126   run_=iEvent.id().run();
00127   event_=iEvent.id().event();
00128  
00129   if(debug_)
00130     cout<<"LhcTrackAnalyzer::analyze() looping over "<< trackCollectionHandle->size()<< "tracks." << endl;    
00131   
00132   // unsigned int i = 0;   
00133   for(TrackCollection::const_iterator track = trackCollectionHandle->begin(); track!= trackCollectionHandle->end(); ++track)
00134     {
00135       if ( nTracks_ >= nMaxtracks_ ) {
00136         std::cout << " LhcTrackAnalyzer::analyze() : Warning - Run "<< run_<<" Event "<< event_<<"\tNumber of tracks: " <<  trackCollectionHandle->size()<< " , greater than " << nMaxtracks_ << std::endl;
00137           continue;
00138          
00139         }
00140           pt_[nTracks_]       = track->pt();
00141           eta_[nTracks_]      = track->eta();
00142           phi_[nTracks_]      = track->phi();
00143           chi2_[nTracks_]     = track->chi2();
00144           chi2ndof_[nTracks_] = track->normalizedChi2();
00145           charge_[nTracks_]   = track->charge();
00146           qoverp_[nTracks_]   = track->qoverp();
00147           dz_[nTracks_]       = track->dz();
00148           dxy_[nTracks_]      = track->dxy();
00149           xPCA_[nTracks_]     = track->vertex().x();
00150           yPCA_[nTracks_]     = track->vertex().y();
00151           zPCA_[nTracks_]     = track->vertex().z(); 
00152           validhits_[nTracks_][0]=track->numberOfValidHits();
00153           validhits_[nTracks_][1]=track->hitPattern().numberOfValidPixelBarrelHits();
00154           validhits_[nTracks_][2]=track->hitPattern().numberOfValidPixelEndcapHits();
00155           validhits_[nTracks_][3]=track->hitPattern().numberOfValidStripTIBHits();
00156           validhits_[nTracks_][4]=track->hitPattern().numberOfValidStripTIDHits();
00157           validhits_[nTracks_][5]=track->hitPattern().numberOfValidStripTOBHits();
00158           validhits_[nTracks_][6]=track->hitPattern().numberOfValidStripTECHits();
00159 
00160 
00161 
00162           int myalgo=-88;
00163           if(track->algo()==reco::TrackBase::undefAlgorithm)myalgo=0;
00164           if(track->algo()==reco::TrackBase::ctf)myalgo=1;
00165           if(track->algo()==reco::TrackBase::iter0)myalgo=4;
00166           if(track->algo()==reco::TrackBase::iter1)myalgo=5;
00167           if(track->algo()==reco::TrackBase::iter2)myalgo=6;
00168           if(track->algo()==reco::TrackBase::iter3)myalgo=7;
00169           if(track->algo()==reco::TrackBase::iter4)myalgo=8;
00170           if(track->algo()==reco::TrackBase::iter5)myalgo=9;
00171           if(track->algo()==reco::TrackBase::iter6)myalgo=10;
00172           if(track->algo()==reco::TrackBase::iter7)myalgo=11;
00173           trkAlgo_[nTracks_]  = myalgo;
00174 
00175           int myquality=-99;
00176           if(track->quality(reco::TrackBase::undefQuality))myquality=-1;
00177           if(track->quality(reco::TrackBase::loose))myquality=0;
00178           if(track->quality(reco::TrackBase::tight))myquality=1;
00179           if(track->quality(reco::TrackBase::highPurity))myquality=2;
00180           //if(track->quality(reco::TrackBase::confirmed))myquality=3;
00181           // if(track->quality(reco::TrackBase::goodIterative))myquality=4;
00182           // if(track->quality(reco::TrackBase::qualitySize))myquality=5;         
00183           trkQuality_[nTracks_]= myquality;
00184 
00185           if(track->quality(reco::TrackBase::highPurity))isHighPurity_[nTracks_]=1;
00186           else isHighPurity_[nTracks_]=0;
00187           nTracks_++;
00188 
00189 
00190         
00191 
00192      
00193     }//end loop on tracks
00194 
00195   for(int d=0;d<nTracks_;++d){
00196     if(abs(trkQuality_[d])>5)cout<<"MYQUALITY!!! " <<trkQuality_[d] <<" at track # "<<d<<"/"<< nTracks_<<endl;
00197   }
00198 
00199  
00200 
00201 
00202   rootTree_->Fill();
00203 } 
00204 
00205 
00206 // ------------ method called once each job before begining the event loop  ------------
00207 void LhcTrackAnalyzer::beginJob()
00208 {
00209   edm::LogInfo("beginJob") << "Begin Job" << std::endl;
00210   // Define TTree for output
00211   rootFile_ = new TFile(filename_.c_str(),"recreate");
00212   rootTree_ = new TTree("tree","Lhc Track tree");
00213   
00214   // Track Paramters 
00215   rootTree_->Branch("run",&run_,"run/I");
00216   rootTree_->Branch("event",&event_,"event/I");
00217   rootTree_->Branch("goodbx",&goodbx_,"goodbx/O");
00218   rootTree_->Branch("goodvtx",&goodvtx_,"goodvtx/O");
00219   rootTree_->Branch("nTracks",&nTracks_,"nTracks/I");
00220   rootTree_->Branch("pt",&pt_,"pt[nTracks]/D");
00221   rootTree_->Branch("eta",&eta_,"eta[nTracks]/D");
00222   rootTree_->Branch("phi",&phi_,"phi[nTracks]/D");
00223   rootTree_->Branch("chi2",&chi2_,"chi2[nTracks]/D");
00224   rootTree_->Branch("chi2ndof",&chi2ndof_,"chi2ndof[nTracks]/D");
00225   rootTree_->Branch("charge",&charge_,"charge[nTracks]/I");
00226   rootTree_->Branch("qoverp",&qoverp_,"qoverp[nTracks]/D");
00227   rootTree_->Branch("dz",&dz_,"dz[nTracks]/D");
00228   rootTree_->Branch("dxy",&dxy_,"dxy[nTracks]/D");
00229   rootTree_->Branch("xPCA",&xPCA_,"xPCA[nTracks]/D");
00230   rootTree_->Branch("yPCA",&yPCA_,"yPCA[nTracks]/D");
00231   rootTree_->Branch("zPCA",&zPCA_,"zPCA[nTracks]/D");
00232  rootTree_->Branch("isHighPurity",&isHighPurity_,"isHighPurity[nTracks]/I");
00233   rootTree_->Branch("trkQuality",&trkQuality_,"trkQuality[nTracks]/I");
00234   rootTree_->Branch("trkAlgo",&trkAlgo_,"trkAlgo[nTracks]/I");
00235   rootTree_->Branch("nValidHits",&validhits_,"nValidHits[nTracks][7]/I");
00236  
00237   
00238 }
00239 
00240 // ------------ method called once each job just after ending the event loop  ------------
00241 void LhcTrackAnalyzer::endJob() 
00242 {
00243    if ( rootFile_ ) {
00244      rootFile_->Write();
00245      rootFile_->Close();
00246    }
00247 
00248    
00249 
00250 }
00251 
00252 void LhcTrackAnalyzer::SetVarToZero() {
00253   run_=-1;
00254   event_=-99;
00255   nTracks_ = 0;
00256   for ( int i=0; i<nMaxtracks_; ++i ) {
00257     pt_[i]        = 0;
00258     eta_[i]       = 0;
00259     phi_[i]       = 0;
00260     chi2_[i]      = 0;
00261     chi2ndof_[i]  = 0;
00262     charge_[i]    = 0;
00263     qoverp_[i]    = 0;
00264     dz_[i]        = 0;
00265     dxy_[i]       = 0;
00266     xPCA_[i]        = 0;
00267     yPCA_[i]        = 0;
00268     zPCA_[i]        = 0;
00269     trkQuality_[i]  = 0;
00270     trkAlgo_[i]     = -1;
00271     isHighPurity_[i]=-3;
00272     for(int j=0;j<7;j++){
00273       validhits_[nTracks_][j]=-1*j;
00274     }
00275   }
00276 
00277 
00278 }
00279 
00280 //define this as a plug-in
00281 DEFINE_FWK_MODULE(LhcTrackAnalyzer);