CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/Validation/RecoEgamma/plugins/ElectronConversionRejectionValidator.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 //
00003 #include "FWCore/ServiceRegistry/interface/Service.h"
00004 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00005 //
00006 #include "Validation/RecoEgamma/plugins/ElectronConversionRejectionValidator.h"
00007 
00008 //
00009 #include "FWCore/Framework/interface/MakerMacros.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/Utilities/interface/Exception.h"
00012 //
00013 
00014 //
00015 #include "DataFormats/Common/interface/Handle.h"
00016 #include "DataFormats/TrackReco/interface/Track.h"
00017 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00018 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00019 #include "DataFormats/EgammaCandidates/interface/Conversion.h"
00020 #include "DataFormats/EgammaCandidates/interface/ConversionFwd.h"
00021 #include "DataFormats/EgammaCandidates/interface/Photon.h"
00022 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
00023 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00024 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00025 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00026 #include "DataFormats/Math/interface/deltaPhi.h"
00027 #include "RecoEgamma/EgammaTools/interface/ConversionTools.h"
00028 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00029 
00030 //
00031 //
00032 #include "TFile.h"
00033 #include "TH1.h"
00034 #include "TH2.h"
00035 #include "TTree.h"
00036 #include "TVector3.h"
00037 #include "TProfile.h"
00038 #include "TMath.h"
00039 //
00050 using namespace std;
00051 
00052 
00053 ElectronConversionRejectionValidator::ElectronConversionRejectionValidator( const edm::ParameterSet& pset )
00054   {
00055 
00056     fName_     = pset.getUntrackedParameter<std::string>("Name");
00057     verbosity_ = pset.getUntrackedParameter<int>("Verbosity");
00058     parameters_ = pset;
00059     
00060     gsfElectronCollectionProducer_ = pset.getParameter<std::string>("gsfElectronProducer");
00061     gsfElectronCollection_ = pset.getParameter<std::string>("gsfElectronCollection");
00062 
00063     conversionCollectionProducer_ = pset.getParameter<std::string>("convProducer");
00064     conversionCollection_ = pset.getParameter<std::string>("conversionCollection");
00065     // conversionTrackProducer_ = pset.getParameter<std::string>("trackProducer");
00066 
00067     isRunCentrally_=   pset.getParameter<bool>("isRunCentrally");
00068     
00069     elePtMin_ = pset.getParameter<double>("elePtMin");
00070     eleExpectedHitsInnerMax_ = pset.getParameter<int>("eleExpectedHitsInnerMax");
00071     eleD0Max_ = pset.getParameter<double>("eleD0Max");
00072     
00073   }
00074 
00075 
00076 
00077 
00078 
00079 ElectronConversionRejectionValidator::~ElectronConversionRejectionValidator() {}
00080 
00081 
00082 
00083 
00084 void  ElectronConversionRejectionValidator::beginJob() {
00085 
00086   nEvt_=0;
00087   nEntry_=0;
00088 
00089 
00090   dbe_ = 0;
00091   dbe_ = edm::Service<DQMStore>().operator->();
00092 
00093 
00094   double ptMin = parameters_.getParameter<double>("ptMin");
00095   double ptMax = parameters_.getParameter<double>("ptMax");
00096   int ptBin = parameters_.getParameter<int>("ptBin");
00097 
00098   double trackptMin = parameters_.getParameter<double>("trackptMin");
00099   double trackptMax = parameters_.getParameter<double>("trackptMax");
00100   int trackptBin = parameters_.getParameter<int>("trackptBin");  
00101 
00102 //   double resMin = parameters_.getParameter<double>("resMin");
00103 //   double resMax = parameters_.getParameter<double>("resMax");
00104 //   int resBin = parameters_.getParameter<int>("resBin");
00105 
00106   double etaMin = parameters_.getParameter<double>("etaMin");
00107   double etaMax = parameters_.getParameter<double>("etaMax");
00108   int etaBin = parameters_.getParameter<int>("etaBin");
00109 
00110   double phiMin = -TMath::Pi();
00111   double phiMax =  TMath::Pi();
00112   int    phiBin = parameters_.getParameter<int>("phiBin");
00113 
00114 
00115   double rhoMin = parameters_.getParameter<double>("rhoMin");
00116   double rhoMax = parameters_.getParameter<double>("rhoMax");
00117   int    rhoBin = parameters_.getParameter<int>("rhoBin");
00118 
00119   double zMin = parameters_.getParameter<double>("zMin");
00120   double zMax = parameters_.getParameter<double>("zMax");
00121   int    zBin = parameters_.getParameter<int>("zBin");
00122 
00123 //   double dPhiTracksMin = parameters_.getParameter<double>("dPhiTracksMin");
00124 //   double dPhiTracksMax = parameters_.getParameter<double>("dPhiTracksMax");
00125 //   int dPhiTracksBin = parameters_.getParameter<int>("dPhiTracksBin");
00126 
00127 //   double eoverpMin = parameters_.getParameter<double>("eoverpMin");
00128 //   double eoverpMax = parameters_.getParameter<double>("eoverpMax");
00129 //   int    eoverpBin = parameters_.getParameter<int>("eoverpBin");
00130 
00131 
00132   //  double dEtaTracksMin = parameters_.getParameter<double>("dEtaTracksMin");  // unused
00133   //  double dEtaTracksMax = parameters_.getParameter<double>("dEtaTracksMax"); // unused
00134   //  int    dEtaTracksBin = parameters_.getParameter<int>("dEtaTracksBin");  // unused
00135 
00136 //   double dCotTracksMin = parameters_.getParameter<double>("dCotTracksMin");
00137 //   double dCotTracksMax = parameters_.getParameter<double>("dCotTracksMax");
00138 //   int    dCotTracksBin = parameters_.getParameter<int>("dCotTracksBin");
00139 
00140 
00141   if (dbe_) {
00142 
00144     // SC from reco photons
00145 
00146     //TString simfolder = TString(
00147     //std::string simpath = dqmpath_ + "SimulationInfo";
00148     dbe_->setCurrentFolder(dqmpath_);
00149     //
00150     // simulation information about conversions
00151     // Histograms for efficiencies
00152     h_elePtAll_ = dbe_->book1D("elePtAll","# of Electrons",ptBin,ptMin,ptMax);
00153     h_eleEtaAll_ = dbe_->book1D("eleEtaAll","# of Electrons",etaBin,etaMin,etaMax);
00154     h_elePhiAll_ = dbe_->book1D("elePhiAll","# of Electrons",phiBin,phiMin,phiMax);
00155     
00156     h_elePtPass_ = dbe_->book1D("elePtPass","# of Electrons",ptBin,ptMin,ptMax);
00157     h_eleEtaPass_ = dbe_->book1D("eleEtaPass","# of Electrons",etaBin,etaMin,etaMax);
00158     h_elePhiPass_ = dbe_->book1D("elePhiPass","# of Electrons",phiBin,phiMin,phiMax);
00159 
00160     h_elePtFail_ = dbe_->book1D("elePtFail","# of Electrons",ptBin,ptMin,ptMax);
00161     h_eleEtaFail_ = dbe_->book1D("eleEtaFail","# of Electrons",etaBin,etaMin,etaMax);
00162     h_elePhiFail_ = dbe_->book1D("elePhiFail","# of Electrons",phiBin,phiMin,phiMax);    
00163     
00164     h_convPt_ = dbe_->book1D("convPt","# of Electrons",ptBin,ptMin,ptMax);
00165     h_convEta_ = dbe_->book1D("convEta","# of Electrons",etaBin,etaMin,etaMax);
00166     h_convPhi_ = dbe_->book1D("convPhi","# of Electrons",phiBin,phiMin,phiMax); 
00167     h_convRho_ = dbe_->book1D("convRho","# of Electrons",rhoBin,rhoMin,rhoMax);        
00168     h_convZ_ = dbe_->book1D("convZ","# of Electrons",zBin,zMin,zMax);        
00169     h_convProb_ = dbe_->book1D("convProb","# of Electrons",100,0.0,1.0);    
00170 
00171     h_convLeadTrackpt_ = dbe_->book1D("convLeadTrackpt","# of Electrons",trackptBin,trackptMin,trackptMax);
00172     h_convTrailTrackpt_ = dbe_->book1D("convTrailTrackpt","# of Electrons",trackptBin,trackptMin,trackptMax);
00173     h_convLog10TrailTrackpt_ = dbe_->book1D("convLog10TrailTrackpt","# of Electrons",ptBin,-2.0,3.0);
00174 
00175     h_convLeadTrackAlgo_ = dbe_->book1D("convLeadTrackAlgo","# of Electrons",31,-0.5,30.5);
00176     h_convTrailTrackAlgo_ = dbe_->book1D("convLeadTrackAlgo","# of Electrons",31,-0.5,30.5);
00177 
00178     
00179   } // if DQM
00180 
00181 
00182 
00183 }
00184 
00185 
00186 
00187  void  ElectronConversionRejectionValidator::beginRun (edm::Run const & r, edm::EventSetup const & theEventSetup) {
00188 
00189 
00190 }
00191 
00192 void  ElectronConversionRejectionValidator::endRun (edm::Run& r, edm::EventSetup const & theEventSetup) {
00193 
00194 
00195 }
00196 
00197 
00198 
00199 void ElectronConversionRejectionValidator::analyze( const edm::Event& e, const edm::EventSetup& esup ) {
00200 
00201   using namespace edm;
00202 
00203 
00204   nEvt_++;
00205   LogInfo("ElectronConversionRejectionValidator") << "ElectronConversionRejectionValidator Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n";
00206   //  std::cout << "ElectronConversionRejectionValidator Analyzing event number: "  << e.id() << " Global Counter " << nEvt_ <<"\n";
00207 
00208 
00209 
00211   Handle<reco::ConversionCollection> convHandle;
00212   e.getByLabel(conversionCollectionProducer_, conversionCollection_ , convHandle);
00213   if (!convHandle.isValid()) {
00214     edm::LogError("ElectronConversionRejectionValidator") << "Error! Can't get the Conversion collection "<< std::endl;
00215     return;
00216   }
00217 
00219   Handle<reco::GsfElectronCollection> gsfElectronHandle;
00220   e.getByLabel(gsfElectronCollectionProducer_, gsfElectronCollection_ , gsfElectronHandle);
00221   const reco::GsfElectronCollection &gsfElectronCollection = *(gsfElectronHandle.product());
00222   if (!gsfElectronHandle.isValid()) {
00223     edm::LogError("ElectronConversionRejectionValidator") << "Error! Can't get the Electron collection "<< std::endl;
00224     return;
00225   }
00226 
00227   // offline  Primary vertex
00228   edm::Handle<reco::VertexCollection> vertexHandle;
00229   e.getByLabel("offlinePrimaryVertices", vertexHandle);
00230   if (!vertexHandle.isValid()) {
00231       edm::LogError("ElectronConversionRejectionValidator") << "Error! Can't get the product primary Vertex Collection "<< "\n";
00232       return;
00233   }
00234   const reco::Vertex &thevtx = vertexHandle->at(0);
00235   
00236   edm::Handle<reco::BeamSpot> bsHandle;
00237   e.getByLabel("offlineBeamSpot", bsHandle);
00238   if (!bsHandle.isValid()) {
00239       edm::LogError("ElectronConversionRejectionValidator") << "Error! Can't get the product beamspot Collection "<< "\n";
00240       return;
00241   }
00242   const reco::BeamSpot &thebs = *bsHandle.product();
00243 
00244  
00245   //loop over electrons
00246   for (reco::GsfElectronCollection::const_iterator iele = gsfElectronCollection.begin(); iele!=gsfElectronCollection.end(); ++iele) {
00247     //apply basic pre-selection cuts to remove the conversions with obviously displaced tracks which will anyways be
00248     //removed from the analysis by the hit pattern or impact parameter requirements
00249     if (iele->pt() < elePtMin_) continue;
00250     if (iele->gsfTrack()->trackerExpectedHitsInner().numberOfHits() > eleExpectedHitsInnerMax_) continue;
00251     if ( std::abs(iele->gsfTrack()->dxy(thevtx.position())) > eleD0Max_ ) continue;
00252     
00253     //fill information for all electrons
00254     h_elePtAll_->Fill(iele->pt());
00255     h_eleEtaAll_->Fill(iele->eta());
00256     h_elePhiAll_->Fill(iele->phi());
00257     
00258     
00259     //find matching conversion if any
00260     reco::ConversionRef convref = ConversionTools::matchedConversion(*iele,convHandle,thebs.position());
00261     //fill information on passing electrons only if there is no matching conversion (electron passed the conversion rejection cut!)
00262     if (convref.isNull()) {  
00263       h_elePtPass_->Fill(iele->pt());
00264       h_eleEtaPass_->Fill(iele->eta());
00265       h_elePhiPass_->Fill(iele->phi());
00266     }
00267     else {
00268       //matching conversion, electron failed conversion rejection cut
00269       //fill information on electron and matching conversion
00270       //(Note that in case of multiple matching conversions passing the requirements, the conversion tools returns the one closest to the IP,
00271       //which is most likely to be the conversion of the primary photon in case there was one.)
00272       
00273       //fill electron info
00274       h_elePtFail_->Fill(iele->pt());
00275       h_eleEtaFail_->Fill(iele->eta());
00276       h_elePhiFail_->Fill(iele->phi());
00277       
00278       //fill conversion info
00279       math::XYZVectorF convmom = convref->refittedPairMomentum();
00280       h_convPt_->Fill(convmom.rho());
00281       h_convEta_->Fill(convmom.eta());
00282       h_convPhi_->Fill(convmom.phi());
00283       h_convRho_->Fill(convref->conversionVertex().position().rho());
00284       h_convZ_->Fill(convref->conversionVertex().position().z());
00285       h_convProb_->Fill(ChiSquaredProbability(convref->conversionVertex().chi2(),convref->conversionVertex().ndof()));
00286      
00287       //fill information about conversion tracks
00288       if (convref->tracks().size()<2) continue;
00289       
00290       RefToBase<reco::Track> tk1 = convref->tracks().front();
00291       RefToBase<reco::Track> tk2 = convref->tracks().back();
00292       
00293       RefToBase<reco::Track> tklead;
00294       RefToBase<reco::Track> tktrail;
00295       if (tk1->pt() >= tk2->pt()) {
00296         tklead = tk1;
00297         tktrail = tk2;
00298       }
00299       else {
00300         tklead = tk2;
00301         tktrail = tk1;
00302       }
00303       h_convLeadTrackpt_->Fill(tklead->pt());
00304       h_convTrailTrackpt_->Fill(tktrail->pt());
00305       h_convLog10TrailTrackpt_->Fill(log10(tktrail->pt()));
00306       h_convLeadTrackAlgo_->Fill(tklead->algo());
00307       h_convTrailTrackAlgo_->Fill(tktrail->algo());
00308       
00309     }
00310   }
00311 
00312 }
00313 
00314 
00315 
00316 
00317 
00318 void ElectronConversionRejectionValidator::endJob() {
00319 
00320 
00321   std::string outputFileName = parameters_.getParameter<std::string>("OutputFileName");
00322   if ( ! isRunCentrally_ ) {
00323     dbe_->save(outputFileName);
00324   }
00325 
00326   edm::LogInfo("ElectronConversionRejectionValidator") << "Analyzed " << nEvt_  << "\n";
00327   // std::cout  << "::endJob Analyzed " << nEvt_ << " events " << " with total " << nPho_ << " Photons " << "\n";
00328   //  std::cout  << "ElectronConversionRejectionValidator::endJob Analyzed " << nEvt_ << " events " << "\n";
00329 
00330   return ;
00331 }
00332 
00333