Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:16

0001 //--------------------------------------------------------------------------------------------------
0002 // $Id $
0003 //
0004 // PFIsolationEstimator
0005 //
0006 // Helper Class for calculating PFIsolation for Photons & Electron onthe fly. This class takes
0007 //        PF Particle collection and the reconstructed vertex collection as input.
0008 //
0009 // Authors: Vasundhara Chetluru
0010 //--------------------------------------------------------------------------------------------------
0011 
0012 /// --> NOTE if you want to use this class as standalone without the CMSSW part
0013 ///  you need to uncomment the below line and compile normally with scramv1 b
0014 ///  Then you need just to load it in your root macro the lib with the correct path, eg:
0015 ///  gSystem->Load("/data/benedet/CMSSW_5_2_2/lib/slc5_amd64_gcc462/pluginEGammaEGammaAnalysisTools.so");
0016 
0017 //#define STANDALONE   // <---- this line
0018 
0019 #ifndef PFIsolationEstimator_H
0020 #define PFIsolationEstimator_H
0021 
0022 #ifndef STANDALONE
0023 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0024 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0025 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0026 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h"
0027 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0028 #endif
0029 #include <TROOT.h>
0030 #include "TMVA/Factory.h"
0031 #include "TMVA/Tools.h"
0032 #include "TMVA/Reader.h"
0033 #include "TH1.h"
0034 #include "TH2.h"
0035 
0036 #include "DataFormats/VertexReco/interface/Vertex.h"
0037 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0038 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0039 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0040 #include "DataFormats/ParticleFlowCandidate/interface/PileUpPFCandidate.h"
0041 #include "DataFormats/ParticleFlowCandidate/interface/PileUpPFCandidateFwd.h"
0042 
0043 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0044 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
0045 
0046 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0047 
0048 #include "DataFormats/VertexReco/interface/Vertex.h"
0049 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0050 
0051 class PFIsolationEstimator {
0052 public:
0053   PFIsolationEstimator();
0054   ~PFIsolationEstimator();
0055 
0056   enum VetoType {
0057     kElectron = -1,  // MVA for non-triggering electrons
0058     kPhoton = 1      // MVA for triggering electrons
0059   };
0060 
0061   void initializeElectronIsolation(Bool_t bApplyVeto);
0062   void initializePhotonIsolation(Bool_t bApplyVeto);
0063   void initializeElectronIsolationInRings(Bool_t bApplyVeto, int iNumberOfRings, float fRingSize);
0064   void initializePhotonIsolationInRings(Bool_t bApplyVeto, int iNumberOfRings, float fRingSize);
0065   void initializeRings(int iNumberOfRings, float fRingSize);
0066   Bool_t isInitialized() const { return fisInitialized; }
0067 
0068   float fGetIsolation(const reco::PFCandidate* pfCandidate,
0069                       const reco::PFCandidateCollection* pfParticlesColl,
0070                       reco::VertexRef vtx,
0071                       edm::Handle<reco::VertexCollection> vertices);
0072   std::vector<float> fGetIsolationInRings(const reco::PFCandidate* pfCandidate,
0073                                           const reco::PFCandidateCollection* pfParticlesColl,
0074                                           reco::VertexRef vtx,
0075                                           edm::Handle<reco::VertexCollection> vertices);
0076 
0077   float fGetIsolation(const reco::Photon* photon,
0078                       const reco::PFCandidateCollection* pfParticlesColl,
0079                       reco::VertexRef vtx,
0080                       edm::Handle<reco::VertexCollection> vertices);
0081   std::vector<float> fGetIsolationInRings(const reco::Photon* photon,
0082                                           const reco::PFCandidateCollection* pfParticlesColl,
0083                                           reco::VertexRef vtx,
0084                                           edm::Handle<reco::VertexCollection> vertices);
0085 
0086   float fGetIsolation(const reco::GsfElectron* electron,
0087                       const reco::PFCandidateCollection* pfParticlesColl,
0088                       const reco::VertexRef vtx,
0089                       edm::Handle<reco::VertexCollection> vertices);
0090   std::vector<float> fGetIsolationInRings(const reco::GsfElectron* electron,
0091                                           const reco::PFCandidateCollection* pfParticlesColl,
0092                                           reco::VertexRef vtx,
0093                                           edm::Handle<reco::VertexCollection> vertices);
0094 
0095   reco::VertexRef chargedHadronVertex(edm::Handle<reco::VertexCollection> verticies, const reco::PFCandidate& pfcand);
0096 
0097   void setConeSize(float fValue = 0.4) { fConeSize = fValue; };
0098 
0099   void setParticleType(int iValue) { iParticleType = iValue; };
0100 
0101   //Veto booleans
0102   void setApplyVeto(Bool_t bValue = kTRUE) { bApplyVeto = bValue; };
0103   void setApplyPFPUVeto(Bool_t bValue = kFALSE) { bApplyPFPUVeto = bValue; };
0104   void setApplyDzDxyVeto(Bool_t bValue = kTRUE) { bApplyDzDxyVeto = bValue; };
0105   void setApplyMissHitPhVeto(Bool_t bValue = kFALSE) { bApplyMissHitPhVeto = bValue; };
0106   void setDeltaRVetoBarrel(Bool_t bValue = kTRUE) { bDeltaRVetoBarrel = bValue; };
0107   void setDeltaRVetoEndcap(Bool_t bValue = kTRUE) { bDeltaRVetoEndcap = bValue; };
0108   void setRectangleVetoBarrel(Bool_t bValue = kTRUE) { bRectangleVetoBarrel = bValue; };
0109   void setRectangleVetoEndcap(Bool_t bValue = kTRUE) { bRectangleVetoEndcap = bValue; };
0110   //Use crystal size
0111   void setUseCrystalSize(Bool_t bValue = kFALSE) { bUseCrystalSize = bValue; };
0112 
0113   //Veto Values
0114   void setDeltaRVetoBarrelPhotons(float fValue = -1.0) { fDeltaRVetoBarrelPhotons = fValue; };
0115   void setDeltaRVetoBarrelNeutrals(float fValue = -1.0) { fDeltaRVetoBarrelNeutrals = fValue; };
0116   void setDeltaRVetoBarrelCharged(float fValue = -1.0) { fDeltaRVetoBarrelCharged = fValue; };
0117   void setDeltaRVetoEndcapPhotons(float fValue = -1.0) { fDeltaRVetoEndcapPhotons = fValue; };
0118   void setDeltaRVetoEndcapNeutrals(float fValue = -1.0) { fDeltaRVetoEndcapNeutrals = fValue; };
0119   void setDeltaRVetoEndcapCharged(float fValue = -1.0) { fDeltaRVetoEndcapCharged = fValue; };
0120   void setNumberOfCrystalEndcapPhotons(float fValue = -1) { fNumberOfCrystalEndcapPhotons = fValue; };
0121 
0122   void setRectangleDeltaPhiVetoBarrelPhotons(float fValue = -1.0) { fRectangleDeltaPhiVetoBarrelPhotons = fValue; };
0123   void setRectangleDeltaPhiVetoBarrelNeutrals(float fValue = -1.0) { fRectangleDeltaPhiVetoBarrelNeutrals = fValue; };
0124   void setRectangleDeltaPhiVetoBarrelCharged(float fValue = -1.0) { fRectangleDeltaPhiVetoBarrelCharged = fValue; };
0125   void setRectangleDeltaPhiVetoEndcapPhotons(float fValue = -1.0) { fRectangleDeltaPhiVetoEndcapPhotons = fValue; };
0126   void setRectangleDeltaPhiVetoEndcapNeutrals(float fValue = -1.0) { fRectangleDeltaPhiVetoEndcapNeutrals = fValue; };
0127   void setRectangleDeltaPhiVetoEndcapCharged(float fValue = -1.0) { fRectangleDeltaPhiVetoEndcapCharged = fValue; };
0128 
0129   void setRectangleDeltaEtaVetoBarrelPhotons(float fValue = -1.0) { fRectangleDeltaEtaVetoBarrelPhotons = fValue; };
0130   void setRectangleDeltaEtaVetoBarrelNeutrals(float fValue = -1.0) { fRectangleDeltaEtaVetoBarrelNeutrals = fValue; };
0131   void setRectangleDeltaEtaVetoBarrelCharged(float fValue = -1.0) { fRectangleDeltaEtaVetoBarrelCharged = fValue; };
0132   void setRectangleDeltaEtaVetoEndcapPhotons(float fValue = -1.0) { fRectangleDeltaEtaVetoEndcapPhotons = fValue; };
0133   void setRectangleDeltaEtaVetoEndcapNeutrals(float fValue = -1.0) { fRectangleDeltaEtaVetoEndcapNeutrals = fValue; };
0134   void setRectangleDeltaEtaVetoEndcapCharged(float fValue = -1.0) { fRectangleDeltaEtaVetoEndcapCharged = fValue; };
0135 
0136   //Veto implementation
0137   float isPhotonParticleVetoed(const reco::PFCandidate* pfIsoCand);
0138   float isNeutralParticleVetoed(const reco::PFCandidate* pfIsoCand);
0139   float isChargedParticleVetoed(const reco::PFCandidate* pfIsoCand, edm::Handle<reco::VertexCollection> vertices);
0140   float isChargedParticleVetoed(const reco::PFCandidate* pfIsoCand,
0141                                 reco::VertexRef vtx,
0142                                 edm::Handle<reco::VertexCollection> vertices);
0143 
0144   float getIsolationPhoton() {
0145     fIsolationPhoton = fIsolationInRingsPhoton[0];
0146     return fIsolationPhoton;
0147   };
0148   float getIsolationNeutral() {
0149     fIsolationNeutral = fIsolationInRingsNeutral[0];
0150     return fIsolationNeutral;
0151   };
0152   float getIsolationCharged() {
0153     fIsolationCharged = fIsolationInRingsCharged[0];
0154     return fIsolationCharged;
0155   };
0156   float getIsolationChargedAll() { return fIsolationChargedAll; };
0157 
0158   std::vector<float> getIsolationInRingsPhoton() { return fIsolationInRingsPhoton; };
0159   std::vector<float> getIsolationInRingsNeutral() { return fIsolationInRingsNeutral; };
0160   std::vector<float> getIsolationInRingsCharged() { return fIsolationInRingsCharged; };
0161   std::vector<float> getIsolationInRingsChargedAll() { return fIsolationInRingsChargedAll; };
0162 
0163   void setNumbersOfRings(int iValue = 1) { iNumberOfRings = iValue; };
0164   void setRingSize(float fValue = 0.4) { fRingSize = fValue; };
0165 
0166   int getNumbersOfRings() { return iNumberOfRings; };
0167   float getRingSize() { return fRingSize; };
0168 
0169   int matchPFObject(const reco::Photon* photon, const reco::PFCandidateCollection* pfParticlesColl);
0170   int matchPFObject(const reco::GsfElectron* photon, const reco::PFCandidateCollection* pfParticlesColl);
0171 
0172 private:
0173   int iParticleType;
0174 
0175   Bool_t fisInitialized;
0176   float fIsolation;
0177   float fIsolationPhoton;
0178   float fIsolationNeutral;
0179   float fIsolationCharged;
0180   float fIsolationChargedAll;
0181 
0182   std::vector<float> fIsolationInRings;
0183   std::vector<float> fIsolationInRingsPhoton;
0184   std::vector<float> fIsolationInRingsNeutral;
0185   std::vector<float> fIsolationInRingsCharged;
0186   std::vector<float> fIsolationInRingsChargedAll;
0187 
0188   Bool_t checkClosestZVertex;
0189   float fConeSize;
0190   Bool_t bApplyVeto;
0191   Bool_t bApplyDzDxyVeto;
0192   Bool_t bApplyPFPUVeto;
0193   Bool_t bApplyMissHitPhVeto;
0194   Bool_t bUseCrystalSize;
0195 
0196   Bool_t bDeltaRVetoBarrel;
0197   Bool_t bDeltaRVetoEndcap;
0198 
0199   Bool_t bRectangleVetoBarrel;
0200   Bool_t bRectangleVetoEndcap;
0201 
0202   float fDeltaRVetoBarrelPhotons;
0203   float fDeltaRVetoBarrelNeutrals;
0204   float fDeltaRVetoBarrelCharged;
0205 
0206   float fDeltaRVetoEndcapPhotons;
0207   float fDeltaRVetoEndcapNeutrals;
0208   float fDeltaRVetoEndcapCharged;
0209 
0210   float fNumberOfCrystalEndcapPhotons;
0211 
0212   float fRectangleDeltaPhiVetoBarrelPhotons;
0213   float fRectangleDeltaPhiVetoBarrelNeutrals;
0214   float fRectangleDeltaPhiVetoBarrelCharged;
0215 
0216   float fRectangleDeltaPhiVetoEndcapPhotons;
0217   float fRectangleDeltaPhiVetoEndcapNeutrals;
0218   float fRectangleDeltaPhiVetoEndcapCharged;
0219 
0220   float fRectangleDeltaEtaVetoBarrelPhotons;
0221   float fRectangleDeltaEtaVetoBarrelNeutrals;
0222   float fRectangleDeltaEtaVetoBarrelCharged;
0223 
0224   float fRectangleDeltaEtaVetoEndcapPhotons;
0225   float fRectangleDeltaEtaVetoEndcapNeutrals;
0226   float fRectangleDeltaEtaVetoEndcapCharged;
0227 
0228   int iNumberOfRings;
0229   int iMissHits;
0230 
0231   float fRingSize;
0232 
0233   float fDeltaR;
0234   float fDeltaEta;
0235   float fDeltaPhi;
0236 
0237   float fEta;
0238   float fPhi;
0239   float fEtaSC;
0240   float fPhiSC;
0241 
0242   float fPt;
0243   float fVx;
0244   float fVy;
0245   float fVz;
0246 
0247   reco::SuperClusterRef refSC;
0248   bool pivotInBarrel;
0249 
0250   math::XYZVector vtxWRTCandidate;
0251 
0252   void initialize(Bool_t bApplyVeto, int iParticleType);
0253 };
0254 
0255 #endif