CMS 3D CMS Logo

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