CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoBTag/SoftLepton/plugins/SoftPFElectronProducer.cc

Go to the documentation of this file.
00001 #include "RecoBTag/SoftLepton/plugins/SoftPFElectronProducer.h"
00002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00003 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00004 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00005 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00006 #include "DataFormats/Math/interface/deltaR.h"
00007 #include "DataFormats/Math/interface/Vector3D.h"
00008 #include "DataFormats/Math/interface/Point3D.h"
00009 #include "FWCore/Utilities/interface/EDMException.h"
00010 #include "FWCore/Framework/interface/Frameworkfwd.h"
00011 
00012 #include <cmath>
00013 
00014 SoftPFElectronProducer::SoftPFElectronProducer (const edm::ParameterSet& conf)
00015 {
00016     gsfElectronTag_ = conf.getParameter<edm::InputTag>("Electrons");
00017 
00018     barrelPtCuts_                    = conf.getParameter<std::vector<double> >("BarrelPtCuts");
00019     barreldRGsfTrackElectronCuts_    = conf.getParameter<std::vector<double> >("BarreldRGsfTrackElectronCuts");
00020     barrelEemPinRatioCuts_           = conf.getParameter<std::vector<double> >("BarrelEemPinRatioCuts");
00021     barrelMVACuts_                   = conf.getParameter<std::vector<double> >("BarrelMVACuts");
00022     barrelInversedRFirstLastHitCuts_ = conf.getParameter<std::vector<double> >("BarrelInversedRFirstLastHitCuts");
00023     barrelRadiusFirstHitCuts_        = conf.getParameter<std::vector<double> >("BarrelRadiusFirstHitCuts");
00024     barrelZFirstHitCuts_             = conf.getParameter<std::vector<double> >("BarrelZFirstHitCuts");
00025 
00026     forwardPtCuts_                   = conf.getParameter<std::vector<double> >("ForwardPtCuts");
00027     forwardInverseFBremCuts_         = conf.getParameter<std::vector<double> >("ForwardInverseFBremCuts");
00028     forwarddRGsfTrackElectronCuts_   = conf.getParameter<std::vector<double> >("ForwarddRGsfTrackElectronCuts");
00029     forwardRadiusFirstHitCuts_       = conf.getParameter<std::vector<double> >("ForwardRadiusFirstHitCuts");
00030     forwardZFirstHitCuts_            = conf.getParameter<std::vector<double> >("ForwardZFirstHitCuts");
00031     forwardMVACuts_                  = conf.getParameter<std::vector<double> >("ForwardMVACuts");
00032 
00033     produces<reco::GsfElectronCollection>();
00034 }
00035 
00036 SoftPFElectronProducer::~SoftPFElectronProducer()
00037 {
00038 
00039 }
00040 
00041 void SoftPFElectronProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00042 {
00043   std::auto_ptr<reco::GsfElectronCollection> output(new reco::GsfElectronCollection);
00044 
00045   edm::Handle<reco::GsfElectronCollection> gsfCandidates;
00046   iEvent.getByLabel(gsfElectronTag_, gsfCandidates);
00047 
00048   // this will throw if the input collection is not present
00049   const reco::GsfElectronCollection & gsfCollection = *(gsfCandidates.product());
00050 
00051   for (unsigned int i = 0 ; i != gsfCollection.size(); ++i)
00052   {
00053     if (not isClean(gsfCollection[i]))
00054       continue;
00055 
00056     output->push_back(gsfCollection[i]);
00057   }
00058 
00059   iEvent.put(output);
00060 }
00061 
00062 bool SoftPFElectronProducer::isClean(const reco::GsfElectron& gsfcandidate)
00063 {
00064   const double ecalEnergy            = gsfcandidate.superCluster()->energy();
00065   const double pIn                   = gsfcandidate.gsfTrack()->innerMomentum().R();
00066   const double pOut                  = gsfcandidate.gsfTrack()->outerMomentum().R();
00067   const double EemPinRatio           = (ecalEnergy - pIn)/(ecalEnergy + pIn);
00068   //const double EemPoutRatio          = (ecalEnergy - pOut)/(ecalEnergy + pOut);
00069   const double pt                    = gsfcandidate.pt();
00070   const double fbrem                 = (pIn - pOut)/pIn;
00071   const double dRGsfTrackElectron    = deltaR(gsfcandidate.gsfTrack()->eta(), gsfcandidate.gsfTrack()->phi(), gsfcandidate.eta(), gsfcandidate.phi());
00072   const double mva                   = gsfcandidate.mva();
00073   const math::XYZPoint firstHit      = gsfcandidate.gsfTrack()->innerPosition();
00074   const math::XYZPoint lastHit       = gsfcandidate.gsfTrack()->outerPosition();
00075   const double inversedRFirstLastHit = 1.0/deltaR(firstHit.eta(), firstHit.phi(), lastHit.eta(), lastHit.phi());
00076   const double radiusFirstHit        = firstHit.Rho();
00077   const double zFirstHit             = firstHit.Z();
00078   /*std::cout << "This particle has " << pt << " "
00079                                     << 1.0/fbrem << " "
00080                                     << dRGsfTrackElectron << " "
00081                                     << EemPoutRatio << " "
00082                                     << mva << std::endl;*/
00083   if(fabs(gsfcandidate.eta()) < 1.5)
00084   {
00085     // use barrel cuts
00086     if( barrelPtCuts_.front()                    > pt                    || barrelPtCuts_.back()                    < pt)                    return false;
00087     if( barreldRGsfTrackElectronCuts_.front()    > dRGsfTrackElectron    || barreldRGsfTrackElectronCuts_.back()    < dRGsfTrackElectron)    return false;
00088     if( barrelEemPinRatioCuts_.front()           > EemPinRatio           || barrelEemPinRatioCuts_.back()           < EemPinRatio)           return false;
00089     if( barrelMVACuts_.front()                   > mva                   || barrelMVACuts_.back()                   < mva)                   return false;
00090     if( barrelInversedRFirstLastHitCuts_.front() > inversedRFirstLastHit || barrelInversedRFirstLastHitCuts_.back() < inversedRFirstLastHit) return false;
00091     if( barrelRadiusFirstHitCuts_.front()        > radiusFirstHit        || barrelRadiusFirstHitCuts_.back()        < radiusFirstHit)        return false;
00092     if( barrelZFirstHitCuts_.front()             > zFirstHit             || barrelZFirstHitCuts_.back()             < zFirstHit)             return false;
00093   }
00094   else
00095   {
00096     // use endcap cuts
00097     if( forwardPtCuts_.front()                 > pt                 || forwardPtCuts_.back()                 < pt)                 return false;
00098     if( forwardInverseFBremCuts_.front()       > 1.0/fbrem          || forwardInverseFBremCuts_.back()       < 1.0/fbrem)          return false;
00099     if( forwarddRGsfTrackElectronCuts_.front() > dRGsfTrackElectron || forwarddRGsfTrackElectronCuts_.back() < dRGsfTrackElectron) return false;
00100     if( forwardRadiusFirstHitCuts_.front()     > radiusFirstHit     || forwardRadiusFirstHitCuts_.back()     < radiusFirstHit)     return false;
00101     if( forwardZFirstHitCuts_.front()          > zFirstHit          || forwardZFirstHitCuts_.back()          < zFirstHit)          return false;
00102     if( forwardMVACuts_.front()                > mva                || forwardMVACuts_.back()                < mva)                return false;
00103   }
00104 
00105   return true;
00106 }
00107