CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/Alignment/OfflineValidation/plugins/PrimaryVertexValidation.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    PrimaryVertexValidation
00004 // Class:      PrimaryVertexValidation
00005 // 
00013 //
00014 // Original Author:  Marco Musich
00015 //         Created:  Tue Mar 02 10:39:34 CDT 2010
00016 //
00017 
00018 // system include files
00019 #include <memory>
00020 
00021 
00022 // user include files
00023 #include "Alignment/OfflineValidation/plugins/PrimaryVertexValidation.h"
00024 #include "FWCore/Framework/interface/Frameworkfwd.h"
00025 #include "FWCore/Framework/interface/EDAnalyzer.h"
00026 
00027 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include <FWCore/Framework/interface/ESHandle.h>
00030 #include "FWCore/Framework/interface/MakerMacros.h"
00031 
00032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00033 
00034 #include <SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h>
00035 #include <SimDataFormats/TrackingHit/interface/PSimHit.h>
00036 
00037 #include "TH1F.h"
00038 #include "TH2F.h"
00039 #include "TFile.h"
00040 #include "TROOT.h"
00041 #include "TChain.h"
00042 #include "TNtuple.h"
00043 #include "FWCore/ServiceRegistry/interface/Service.h"
00044 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00045 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00046 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
00047 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
00048 #include <Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h>
00049 #include <DataFormats/GeometrySurface/interface/Surface.h>
00050 #include <DataFormats/GeometrySurface/interface/GloballyPositioned.h>
00051 #include <Geometry/CommonDetUnit/interface/GeomDet.h>
00052 #include "MagneticField/Engine/interface/MagneticField.h" 
00053 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h" 
00054 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
00055 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00056 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00057 
00058 #include "DataFormats/TrackReco/interface/Track.h"
00059 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00060 #include "DataFormats/VertexReco/interface/Vertex.h"
00061 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00062 
00063  const int kBPIX = PixelSubdetector::PixelBarrel;
00064  const int kFPIX = PixelSubdetector::PixelEndcap;
00065 
00066 // Constructor
00067 
00068 PrimaryVertexValidation::PrimaryVertexValidation(const edm::ParameterSet& iConfig)
00069   : theConfig(iConfig), 
00070     Nevt_(0),
00071     theTrackFilter_(iConfig.getParameter<edm::ParameterSet>("TkFilterParameters")),
00072     rootFile_(0),
00073     rootTree_(0)
00074 {
00075   //now do what ever initialization is needed
00076   debug_    = iConfig.getParameter<bool>       ("Debug");  
00077   TrackCollectionTag_      = iConfig.getParameter<edm::InputTag>("TrackCollectionTag");  
00078   filename_ = iConfig.getParameter<std::string>("OutputFileName");
00079 
00080   SetVarToZero();
00081 }
00082    
00083 // Destructor
00084 PrimaryVertexValidation::~PrimaryVertexValidation()
00085 {
00086    // do anything here that needs to be done at desctruction time
00087    // (e.g. close files, deallocate resources etc.)
00088 }
00089 
00090 
00091 //
00092 // member functions
00093 //
00094 
00095 // ------------ method called to for each event  ------------
00096 void
00097 PrimaryVertexValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00098 {
00099 
00100   Nevt_++;  
00101 
00102   //=======================================================
00103   // Initialize Root-tuple variables
00104   //=======================================================
00105 
00106   SetVarToZero();
00107  
00108   //=======================================================
00109   // Retrieve the Magnetic Field information
00110   //=======================================================
00111 
00112   edm::ESHandle<MagneticField> theMGField;
00113   iSetup.get<IdealMagneticFieldRecord>().get( theMGField );
00114 
00115   //=======================================================
00116   // Retrieve the Tracking Geometry information
00117   //=======================================================
00118 
00119   edm::ESHandle<GlobalTrackingGeometry> theTrackingGeometry;
00120   iSetup.get<GlobalTrackingGeometryRecord>().get( theTrackingGeometry );
00121 
00122   //=======================================================
00123   // Retrieve the Track information
00124   //=======================================================
00125   
00126   edm::Handle<reco::TrackCollection>  trackCollectionHandle;
00127   iEvent.getByLabel(TrackCollectionTag_, trackCollectionHandle);
00128   
00129   //=======================================================
00130   // Retrieve offline vartex information (only for reco)
00131   //=======================================================
00132 
00133   /*
00134   double OfflineVertexX = 0.;
00135   double OfflineVertexY = 0.;
00136   double OfflineVertexZ = 0.;
00137   edm::Handle<reco::VertexCollection> vertices;
00138   try {
00139     iEvent.getByLabel("offlinePrimaryVertices", vertices);
00140   } catch (...) {
00141     std::cout << "No offlinePrimaryVertices found!" << std::endl;
00142   }
00143   if ( vertices.isValid() ) {
00144     OfflineVertexX = (*vertices)[0].x();
00145     OfflineVertexY = (*vertices)[0].y();
00146     OfflineVertexZ = (*vertices)[0].z();
00147   }
00148   */
00149   
00150   //=======================================================
00151   // Retrieve Beamspot information (only for reco)
00152   //=======================================================
00153 
00154   /*
00155     BeamSpot beamSpot;
00156     edm::Handle<BeamSpot> beamSpotHandle;
00157     iEvent.getByLabel("offlineBeamSpot", beamSpotHandle);
00158     
00159     if ( beamSpotHandle.isValid() )
00160     {
00161     beamSpot = *beamSpotHandle;
00162     
00163     } else
00164     {
00165     edm::LogInfo("PrimaryVertexValidation")
00166     << "No beam spot available from EventSetup \n";
00167     }
00168     
00169     double BSx0 = beamSpot.x0();
00170     double BSy0 = beamSpot.y0();
00171     double BSz0 = beamSpot.z0();
00172     
00173     if(debug_)
00174     std::cout<<"Beamspot x:"<<BSx0<<" y:"<<BSy0<<" z:"<<BSz0<std::<std::endl; 
00175     
00176     //double sigmaz = beamSpot.sigmaZ();
00177     //double dxdz = beamSpot.dxdz();
00178     //double BeamWidth = beamSpot.BeamWidth();
00179     */
00180   
00181   //=======================================================
00182   // Starts here ananlysis
00183   //=======================================================
00184   
00185   if(debug_)
00186     std::cout<<"PrimaryVertexValidation::analyze() looping over "<<trackCollectionHandle->size()<< "tracks." <<std::endl;       
00187   
00188   unsigned int i = 0;   
00189   for(reco::TrackCollection::const_iterator track = trackCollectionHandle->begin(); track!= trackCollectionHandle->end(); ++track, ++i)
00190     {
00191       if ( nTracks_ >= nMaxtracks_ ) {
00192         std::cout << " PrimaryVertexValidation::analyze() : Warning - Number of tracks: " << nTracks_ << " , greater than " << nMaxtracks_ << std::endl;
00193         continue;
00194       }
00195 
00196       pt_[nTracks_]       = track->pt();
00197       p_[nTracks_]        = track->p();
00198       nhits_[nTracks_]    = track->numberOfValidHits();
00199       eta_[nTracks_]      = track->eta();
00200       phi_[nTracks_]      = track->phi();
00201       chi2_[nTracks_]     = track->chi2();
00202       chi2ndof_[nTracks_] = track->normalizedChi2();
00203       charge_[nTracks_]   = track->charge();
00204       qoverp_[nTracks_]   = track->qoverp();
00205       dz_[nTracks_]       = track->dz();
00206       dxy_[nTracks_]      = track->dxy();
00207       xPCA_[nTracks_]     = track->vertex().x();
00208       yPCA_[nTracks_]     = track->vertex().y();
00209       zPCA_[nTracks_]     = track->vertex().z(); 
00210     
00211       //=======================================================
00212       // Retrieve rechit information
00213       //=======================================================  
00214 
00215       int nRecHit1D =0;
00216       int nRecHit2D =0;
00217       int nhitinTIB =0; 
00218       int nhitinTOB =0; 
00219       int nhitinTID =0; 
00220       int nhitinTEC =0; 
00221       int nhitinBPIX=0;
00222       int nhitinFPIX=0; 
00223       
00224       for (trackingRecHit_iterator iHit = track->recHitsBegin(); iHit != track->recHitsEnd(); ++iHit) {
00225         if((*iHit)->isValid()) {        
00226           
00227           if (this->isHit2D(**iHit)) {++nRecHit2D;}
00228           else {++nRecHit1D; }
00229    
00230           int type =(*iHit)->geographicalId().subdetId();
00231            
00232           if(type==int(StripSubdetector::TIB)){++nhitinTIB;}
00233           if(type==int(StripSubdetector::TOB)){++nhitinTOB;}
00234           if(type==int(StripSubdetector::TID)){++nhitinTID;}
00235           if(type==int(StripSubdetector::TEC)){++nhitinTEC;}
00236           if(type==int(                kBPIX)){++nhitinBPIX;}
00237           if(type==int(                kFPIX)){++nhitinFPIX;}
00238           
00239         }
00240       }      
00241 
00242       nhits1D_[nTracks_]     =nRecHit1D;
00243       nhits2D_[nTracks_]     =nRecHit2D;
00244       nhitsBPIX_[nTracks_]   =nhitinBPIX;
00245       nhitsFPIX_[nTracks_]   =nhitinFPIX;
00246       nhitsTIB_[nTracks_]    =nhitinTIB;
00247       nhitsTID_[nTracks_]    =nhitinTID;
00248       nhitsTOB_[nTracks_]    =nhitinTOB;
00249       nhitsTEC_[nTracks_]    =nhitinTEC;
00250 
00251       //=======================================================
00252       // Good tracks for vertexing selection
00253       //=======================================================  
00254 
00255       reco::TrackRef trackref(trackCollectionHandle,i);
00256       bool hasTheProbeFirstPixelLayerHit = false;
00257       hasTheProbeFirstPixelLayerHit = this->hasFirstLayerPixelHits(trackref);
00258       reco::TransientTrack theTTRef = reco::TransientTrack(trackref, &*theMGField, theTrackingGeometry );
00259       if (theTrackFilter_(theTTRef)&&hasTheProbeFirstPixelLayerHit){
00260         isGoodTrack_[nTracks_]=1;
00261       }
00262       
00263       //=======================================================
00264       // Fit unbiased vertex
00265       //=======================================================  
00266       
00267       std::vector<reco::TransientTrack> transientTracks;
00268       for(size_t j = 0; j < trackCollectionHandle->size(); j++)
00269         {
00270           reco::TrackRef tk(trackCollectionHandle, j);
00271           if( tk == trackref ) continue;
00272           bool hasTheTagFirstPixelLayerHit = false;
00273           hasTheTagFirstPixelLayerHit = this->hasFirstLayerPixelHits(tk);
00274           reco::TransientTrack theTT = reco::TransientTrack(tk, &*theMGField, theTrackingGeometry );
00275           if (theTrackFilter_(theTT)&&hasTheTagFirstPixelLayerHit){
00276             transientTracks.push_back(theTT);
00277           }
00278         }
00279       
00280       if(transientTracks.size() > 2){
00281         
00282         if(debug_)
00283           std::cout <<"PrimaryVertexValidation::analyze() :Transient Track Collection size: "<<transientTracks.size()<<std::endl;
00284         
00285         try{
00286 
00287           VertexFitter<5>* fitter = new AdaptiveVertexFitter;
00288           TransientVertex theFittedVertex = fitter->vertex(transientTracks);
00289           
00290           if(theFittedVertex.isValid ()){
00291             
00292             if(theFittedVertex.hasTrackWeight()){
00293               for(size_t rtracks= 0; rtracks < transientTracks.size(); rtracks++){
00294                 sumOfWeightsUnbiasedVertex_[nTracks_] += theFittedVertex.trackWeight(transientTracks[rtracks]);
00295               }
00296             }
00297 
00298             const math::XYZPoint myVertex(theFittedVertex.position().x(),theFittedVertex.position().y(),theFittedVertex.position().z());
00299             hasRecVertex_[nTracks_]    = 1;
00300             xUnbiasedVertex_[nTracks_] = theFittedVertex.position().x();
00301             yUnbiasedVertex_[nTracks_] = theFittedVertex.position().y();
00302             zUnbiasedVertex_[nTracks_] = theFittedVertex.position().z();
00303             chi2normUnbiasedVertex_[nTracks_] = theFittedVertex.normalisedChiSquared();
00304             chi2UnbiasedVertex_[nTracks_] = theFittedVertex.totalChiSquared();
00305             DOFUnbiasedVertex_[nTracks_] = theFittedVertex.degreesOfFreedom();   
00306             tracksUsedForVertexing_[nTracks_] = transientTracks.size();
00307             dxyFromMyVertex_[nTracks_] = track->dxy(myVertex);
00308             dzFromMyVertex_[nTracks_]  = track->dz(myVertex);
00309             dszFromMyVertex_[nTracks_] = track->dsz(myVertex);
00310             
00311             if(debug_){
00312               std::cout<<"PrimaryVertexValidation::analyze() :myVertex.x()= "<<myVertex.x()<<" myVertex.y()= "<<myVertex.y()<<" theFittedVertex.z()= "<<myVertex.z()<<std::endl;          
00313               std::cout<<"PrimaryVertexValidation::analyze() : track->dz(myVertex)= "<<track->dz(myVertex)<<std::endl;
00314               std::cout<<"PrimaryVertexValidation::analyze() : zPCA -myVertex.z() = "<<(track->vertex().z() -myVertex.z() )<<std::endl; 
00315             }// ends if debug_
00316           } // ends if the fitted vertex is Valid
00317           
00318         }  catch ( cms::Exception& er ) {
00319           LogTrace("PrimaryVertexValidation::analyze RECO")<<"caught std::exception "<<er.what()<<std::endl;
00320         }
00321       } //ends if transientTracks.size() > 2
00322       
00323       else {
00324         std::cout << "PrimaryVertexValidation::analyze() :Not enough tracks to make a vertex.  Returns no vertex info" << std::endl;
00325       }
00326 
00327       ++nTracks_;  
00328         
00329       if(debug_)
00330         std::cout<< "Track "<<i<<" : pT = "<<track->pt()<<std::endl;
00331       
00332     }// for loop on tracks
00333   
00334   rootTree_->Fill();
00335 } 
00336 
00337 // ------------ method called to discriminate 1D from 2D hits  ------------
00338 bool PrimaryVertexValidation::isHit2D(const TrackingRecHit &hit) const
00339 {
00340   if (hit.dimension() < 2) {
00341     return false; // some (muon...) stuff really has RecHit1D
00342   } else {
00343     const DetId detId(hit.geographicalId());
00344     if (detId.det() == DetId::Tracker) {
00345       if (detId.subdetId() == kBPIX || detId.subdetId() == kFPIX) {
00346         return true; // pixel is always 2D
00347       } else { // should be SiStrip now
00348         if (dynamic_cast<const SiStripRecHit2D*>(&hit)) return false; // normal hit
00349         else if (dynamic_cast<const SiStripMatchedRecHit2D*>(&hit)) return true; // matched is 2D
00350         else if (dynamic_cast<const ProjectedSiStripRecHit2D*>(&hit)) return false; // crazy hit...
00351         else {
00352           edm::LogError("UnkownType") << "@SUB=AlignmentTrackSelector::isHit2D"
00353                                       << "Tracker hit not in pixel and neither SiStripRecHit2D nor "
00354                                       << "SiStripMatchedRecHit2D nor ProjectedSiStripRecHit2D.";
00355           return false;
00356         }
00357       }
00358     } else { // not tracker??
00359       edm::LogWarning("DetectorMismatch") << "@SUB=AlignmentTrackSelector::isHit2D"
00360                                           << "Hit not in tracker with 'official' dimension >=2.";
00361       return true; // dimension() >= 2 so accept that...
00362     }
00363   }
00364   // never reached...
00365 }
00366 
00367 // ------------ method to check the presence of pixel hits  ------------
00368 bool PrimaryVertexValidation::hasFirstLayerPixelHits(const reco::TrackRef track)
00369 {
00370   bool accepted = false;
00371   // hit pattern of the track
00372   const reco::HitPattern& p = track->hitPattern();      
00373   for (int i=0; i<p.numberOfHits(); i++) {
00374     uint32_t pattern = p.getHitPattern(i);   
00375     if (p.pixelBarrelHitFilter(pattern) || p.pixelEndcapHitFilter(pattern) ) {
00376       if (p.getLayer(pattern) == 1) {
00377         if (p.validHitFilter(pattern)) {
00378           accepted = true;
00379         }
00380       }
00381     }
00382   }
00383   return accepted;
00384 } 
00385 
00386 // ------------ method called once each job before begining the event loop  ------------
00387 void PrimaryVertexValidation::beginJob()
00388 {
00389   edm::LogInfo("beginJob") << "Begin Job" << std::endl;
00390   // Define TTree for output
00391 
00392   Nevt_    = 0;
00393   
00394   rootFile_ = new TFile(filename_.c_str(),"recreate");
00395   rootTree_ = new TTree("tree","PV Validation tree");
00396   
00397   // Track Paramters 
00398   rootTree_->Branch("nTracks",&nTracks_,"nTracks/I");
00399   rootTree_->Branch("pt",&pt_,"pt[nTracks]/D");
00400   rootTree_->Branch("p",&p_,"p[nTracks]/D");
00401   rootTree_->Branch("nhits",&nhits_,"nhits[nTracks]/I");
00402   rootTree_->Branch("nhits1D",&nhits1D_,"nhits1D[nTracks]/I");
00403   rootTree_->Branch("nhits2D",&nhits2D_,"nhits2D[nTracks]/I");
00404   rootTree_->Branch("nhitsBPIX",&nhitsBPIX_,"nhitsBPIX[nTracks]/I");
00405   rootTree_->Branch("nhitsFPIX",&nhitsFPIX_,"nhitsFPIX[nTracks]/I");
00406   rootTree_->Branch("nhitsTIB",&nhitsTIB_,"nhitsTIB[nTracks]/I");
00407   rootTree_->Branch("nhitsTID",&nhitsTID_,"nhitsTID[nTracks]/I");
00408   rootTree_->Branch("nhitsTOB",&nhitsTOB_,"nhitsTOB[nTracks]/I");
00409   rootTree_->Branch("nhitsTEC",&nhitsTEC_,"nhitsTEC[nTracks]/I");
00410   rootTree_->Branch("eta",&eta_,"eta[nTracks]/D");
00411   rootTree_->Branch("phi",&phi_,"phi[nTracks]/D");
00412   rootTree_->Branch("chi2",&chi2_,"chi2[nTracks]/D");
00413   rootTree_->Branch("chi2ndof",&chi2ndof_,"chi2ndof[nTracks]/D");
00414   rootTree_->Branch("charge",&charge_,"charge[nTracks]/I");
00415   rootTree_->Branch("qoverp",&qoverp_,"qoverp[nTracks]/D");
00416   rootTree_->Branch("dz",&dz_,"dz[nTracks]/D");
00417   rootTree_->Branch("dxy",&dxy_,"dxy[nTracks]/D");
00418   rootTree_->Branch("xPCA",&xPCA_,"xPCA[nTracks]/D");
00419   rootTree_->Branch("yPCA",&yPCA_,"yPCA[nTracks]/D");
00420   rootTree_->Branch("zPCA",&zPCA_,"zPCA[nTracks]/D");
00421   rootTree_->Branch("xUnbiasedVertex",&xUnbiasedVertex_,"xUnbiasedVertex[nTracks]/D");
00422   rootTree_->Branch("yUnbiasedVertex",&yUnbiasedVertex_,"yUnbiasedVertex[nTracks]/D");
00423   rootTree_->Branch("zUnbiasedVertex",&zUnbiasedVertex_,"zUnbiasedVertex[nTracks]/D");
00424   rootTree_->Branch("chi2normUnbiasedVertex",&chi2normUnbiasedVertex_,"chi2normUnbiasedVertex[nTracks]/F");
00425   rootTree_->Branch("chi2UnbiasedVertex",&chi2UnbiasedVertex_,"chi2UnbiasedVertex[nTracks]/F");
00426   rootTree_->Branch("DOFUnbiasedVertex",&DOFUnbiasedVertex_," DOFUnbiasedVertex[nTracks]/F");
00427   rootTree_->Branch("sumOfWeightsUnbiasedVertex",&sumOfWeightsUnbiasedVertex_,"sumOfWeightsUnbiasedVertex[nTracks]/F");
00428   rootTree_->Branch("tracksUsedForVertexing",&tracksUsedForVertexing_,"tracksUsedForVertexing[nTracks]/I");
00429   rootTree_->Branch("dxyFromMyVertex",&dxyFromMyVertex_,"dxyFromMyVertex[nTracks]/D");
00430   rootTree_->Branch("dzFromMyVertex",&dzFromMyVertex_,"dzFromMyVertex[nTracks]/D");
00431   rootTree_->Branch("dszFromMyVertex",&dszFromMyVertex_,"dszFromMyVertex[nTracks]/D");
00432   rootTree_->Branch("hasRecVertex",&hasRecVertex_,"hasRecVertex[nTracks]/I");
00433   rootTree_->Branch("isGoodTrack",&isGoodTrack_,"isGoodTrack[nTracks]/I");
00434 }
00435 
00436 // ------------ method called once each job just after ending the event loop  ------------
00437 void PrimaryVertexValidation::endJob() 
00438 {
00439 
00440   std::cout<<"######################################"<<std::endl;
00441   std::cout<<"Number of analyzed events: "<<Nevt_<<std::endl;
00442   std::cout<<"######################################"<<std::endl;
00443   
00444    if ( rootFile_ ) {
00445      rootFile_->Write();
00446      rootFile_->Close();
00447    }
00448 }
00449 
00450 void PrimaryVertexValidation::SetVarToZero() {
00451   
00452   nTracks_ = 0;
00453   for ( int i=0; i<nMaxtracks_; ++i ) {
00454     pt_[i]        = 0;
00455     p_[i]         = 0;
00456     nhits_[i]     = 0;
00457     nhits1D_[i]   = 0;
00458     nhits2D_[i]   = 0;
00459     nhitsBPIX_[i]  = 0;
00460     nhitsFPIX_[i]  = 0;
00461     nhitsTIB_[i]   = 0;
00462     nhitsTID_[i]   = 0;
00463     nhitsTOB_[i]   = 0;
00464     nhitsTEC_[i]   = 0;
00465     eta_[i]       = 0;
00466     phi_[i]       = 0;
00467     chi2_[i]      = 0;
00468     chi2ndof_[i]  = 0;
00469     charge_[i]    = 0;
00470     qoverp_[i]    = 0;
00471     dz_[i]        = 0;
00472     dxy_[i]       = 0;
00473     xPCA_[i]      = 0;
00474     yPCA_[i]      = 0;
00475     zPCA_[i]      = 0;
00476     xUnbiasedVertex_[i] =0;    
00477     yUnbiasedVertex_[i] =0;
00478     zUnbiasedVertex_[i] =0;
00479     chi2normUnbiasedVertex_[i]=0;
00480     chi2UnbiasedVertex_[i]=0;
00481     DOFUnbiasedVertex_[i]=0;
00482     sumOfWeightsUnbiasedVertex_[i]=0;
00483     tracksUsedForVertexing_[i]=0;
00484     dxyFromMyVertex_[i]=0;
00485     dzFromMyVertex_[i]=0;
00486     dszFromMyVertex_[i]=0;
00487     hasRecVertex_[i] = 0;
00488     isGoodTrack_[i]  = 0;
00489   } 
00490 }
00491 
00492 //define this as a plug-in
00493 DEFINE_FWK_MODULE(PrimaryVertexValidation);