CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFIsolationEstimator.h
Go to the documentation of this file.
1 //--------------------------------------------------------------------------------------------------
2 // $Id $
3 //
4 // PFIsolationEstimator
5 //
6 // Helper Class for calculating PFIsolation for Photons & Electron onthe fly. This class takes
7 // PF Particle collection and the reconstructed vertex collection as input.
8 //
9 // Authors: Vasundhara Chetluru
10 //--------------------------------------------------------------------------------------------------
11 
12 
17 
18 //#define STANDALONE // <---- this line
19 
20 #ifndef PFIsolationEstimator_H
21 #define PFIsolationEstimator_H
22 
23 #ifndef STANDALONE
29 #endif
30 #include <TROOT.h>
31 #include "TMVA/Factory.h"
32 #include "TMVA/Tools.h"
33 #include "TMVA/Reader.h"
34 #include "TH1.h"
35 #include "TH2.h"
36 
37 
44 
47 
49 
52 
53 using namespace std;
54 using namespace edm;
55 using namespace reco;
56 
57 
59  public:
62 
63  enum VetoType {
64  kElectron = -1, // MVA for non-triggering electrons
65  kPhoton = 1 // MVA for triggering electrons
66  };
67 
68 
69  void initializeElectronIsolation( Bool_t bApplyVeto);
70  void initializePhotonIsolation( Bool_t bApplyVeto);
71  void initializeElectronIsolationInRings( Bool_t bApplyVeto, int iNumberOfRings, float fRingSize );
72  void initializePhotonIsolationInRings( Bool_t bApplyVeto, int iNumberOfRings, float fRingSize );
73  void initializeRings(int iNumberOfRings, float fRingSize);
74  Bool_t isInitialized() const { return fisInitialized; }
75 
76 
77  float fGetIsolation(const reco::PFCandidate * pfCandidate,const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices );
78  vector<float > fGetIsolationInRings(const reco::PFCandidate * pfCandidate,const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices);
79 
80  float fGetIsolation(const reco::Photon* photon,const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices );
81  vector<float > fGetIsolationInRings(const reco::Photon* photon,const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices);
82 
83  float fGetIsolation(const reco::GsfElectron* electron,const reco::PFCandidateCollection* pfParticlesColl,const reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices );
84  vector<float > fGetIsolationInRings(const reco::GsfElectron* electron,const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices);
85 
86 
87  VertexRef chargedHadronVertex(edm::Handle< reco::VertexCollection > verticies, const reco::PFCandidate& pfcand );
88 
89  void setConeSize(float fValue = 0.4){ fConeSize = fValue;};
90 
91  void setParticleType(int iValue){iParticleType = iValue;};
92 
93  //Veto booleans
94  void setApplyVeto(Bool_t bValue = kTRUE){ bApplyVeto = bValue;};
95  void setApplyPFPUVeto(Bool_t bValue = kFALSE){bApplyPFPUVeto = bValue;};
96  void setApplyDzDxyVeto(Bool_t bValue = kTRUE){ bApplyDzDxyVeto = bValue;};
97  void setApplyMissHitPhVeto(Bool_t bValue = kFALSE){ bApplyMissHitPhVeto = bValue;};
98  void setDeltaRVetoBarrel(Bool_t bValue = kTRUE){ bDeltaRVetoBarrel = bValue;};
99  void setDeltaRVetoEndcap(Bool_t bValue = kTRUE){ bDeltaRVetoEndcap = bValue;};
100  void setRectangleVetoBarrel(Bool_t bValue = kTRUE){ bRectangleVetoBarrel = bValue;};
101  void setRectangleVetoEndcap(Bool_t bValue = kTRUE){ bRectangleVetoEndcap = bValue;};
102  //Use crystal size
103  void setUseCrystalSize(Bool_t bValue = kFALSE) { bUseCrystalSize = bValue;};
104 
105  //Veto Values
106  void setDeltaRVetoBarrelPhotons(float fValue = -1.0){fDeltaRVetoBarrelPhotons=fValue;};
107  void setDeltaRVetoBarrelNeutrals(float fValue = -1.0){fDeltaRVetoBarrelNeutrals=fValue;};
108  void setDeltaRVetoBarrelCharged(float fValue = -1.0){fDeltaRVetoBarrelCharged=fValue;};
109  void setDeltaRVetoEndcapPhotons(float fValue = -1.0){fDeltaRVetoEndcapPhotons=fValue;};
110  void setDeltaRVetoEndcapNeutrals(float fValue = -1.0){fDeltaRVetoEndcapNeutrals=fValue;};
111  void setDeltaRVetoEndcapCharged(float fValue = -1.0){fDeltaRVetoEndcapCharged=fValue;};
112  void setNumberOfCrystalEndcapPhotons(float fValue = -1){fNumberOfCrystalEndcapPhotons=fValue;};
113 
114  void setRectangleDeltaPhiVetoBarrelPhotons(float fValue = -1.0){fRectangleDeltaPhiVetoBarrelPhotons=fValue;};
115  void setRectangleDeltaPhiVetoBarrelNeutrals(float fValue = -1.0){fRectangleDeltaPhiVetoBarrelNeutrals=fValue;};
116  void setRectangleDeltaPhiVetoBarrelCharged(float fValue = -1.0){fRectangleDeltaPhiVetoBarrelCharged=fValue;};
117  void setRectangleDeltaPhiVetoEndcapPhotons(float fValue = -1.0){fRectangleDeltaPhiVetoEndcapPhotons=fValue;};
118  void setRectangleDeltaPhiVetoEndcapNeutrals(float fValue = -1.0){fRectangleDeltaPhiVetoEndcapNeutrals=fValue;};
119  void setRectangleDeltaPhiVetoEndcapCharged(float fValue = -1.0){fRectangleDeltaPhiVetoEndcapCharged=fValue;};
120 
121 
122  void setRectangleDeltaEtaVetoBarrelPhotons(float fValue = -1.0){fRectangleDeltaEtaVetoBarrelPhotons=fValue;};
123  void setRectangleDeltaEtaVetoBarrelNeutrals(float fValue = -1.0){fRectangleDeltaEtaVetoBarrelNeutrals=fValue;};
124  void setRectangleDeltaEtaVetoBarrelCharged(float fValue = -1.0){fRectangleDeltaEtaVetoBarrelCharged=fValue;};
125  void setRectangleDeltaEtaVetoEndcapPhotons(float fValue = -1.0){fRectangleDeltaEtaVetoEndcapPhotons=fValue;};
126  void setRectangleDeltaEtaVetoEndcapNeutrals(float fValue = -1.0){fRectangleDeltaEtaVetoEndcapNeutrals=fValue;};
127  void setRectangleDeltaEtaVetoEndcapCharged(float fValue = -1.0){fRectangleDeltaEtaVetoEndcapCharged=fValue;};
128 
129  //Veto implementation
130  float isPhotonParticleVetoed( const reco::PFCandidate* pfIsoCand );
131  float isNeutralParticleVetoed( const reco::PFCandidate* pfIsoCand );
132  float isChargedParticleVetoed(const reco::PFCandidate* pfIsoCand, edm::Handle< reco::VertexCollection > vertices);
133  float isChargedParticleVetoed( const reco::PFCandidate* pfIsoCand,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices );
134 
135 
136  float getIsolationPhoton(){ fIsolationPhoton = fIsolationInRingsPhoton[0]; return fIsolationPhoton; };
137  float getIsolationNeutral(){ fIsolationNeutral = fIsolationInRingsNeutral[0]; return fIsolationNeutral; };
138  float getIsolationCharged(){ fIsolationCharged = fIsolationInRingsCharged[0]; return fIsolationCharged; };
139  float getIsolationChargedAll(){ return fIsolationChargedAll; };
140 
141  vector<float > getIsolationInRingsPhoton(){ return fIsolationInRingsPhoton; };
142  vector<float > getIsolationInRingsNeutral(){ return fIsolationInRingsNeutral; };
143  vector<float > getIsolationInRingsCharged(){ return fIsolationInRingsCharged; };
144  vector<float > getIsolationInRingsChargedAll(){ return fIsolationInRingsChargedAll; };
145 
146 
147  void setNumbersOfRings(int iValue = 1){iNumberOfRings = iValue;};
148  void setRingSize(float fValue = 0.4){fRingSize = fValue;};
149 
150  int getNumbersOfRings(){return iNumberOfRings;};
151  float getRingSize(){return fRingSize; };
152 
153  int matchPFObject(const reco::Photon* photon, const reco::PFCandidateCollection* pfParticlesColl );
154  int matchPFObject(const reco::GsfElectron* photon, const reco::PFCandidateCollection* pfParticlesColl );
155 
156  private:
157 
158 
160 
162  float fIsolation;
167 
168  vector<float > fIsolationInRings;
169  vector<float > fIsolationInRingsPhoton;
170  vector<float > fIsolationInRingsNeutral;
171  vector<float > fIsolationInRingsCharged;
173 
175  float fConeSize;
176  Bool_t bApplyVeto;
181 
184 
187 
191 
195 
197 
201 
205 
209 
213 
216 
217  float fRingSize;
218 
219  float fDeltaR;
220  float fDeltaEta;
221  float fDeltaPhi;
222 
223  float fEta;
224  float fPhi;
225  float fEtaSC;
226  float fPhiSC;
227 
228  float fPt;
229  float fVx;
230  float fVy;
231  float fVz;
232 
235 
237 
238  void initialize( Bool_t bApplyVeto, int iParticleType);
239 };
240 
241 #endif
void setRectangleVetoEndcap(Bool_t bValue=kTRUE)
void setDeltaRVetoBarrel(Bool_t bValue=kTRUE)
vector< float > fIsolationInRingsCharged
void setNumberOfCrystalEndcapPhotons(float fValue=-1)
void setDeltaRVetoEndcap(Bool_t bValue=kTRUE)
void setRectangleDeltaPhiVetoBarrelPhotons(float fValue=-1.0)
void setApplyMissHitPhVeto(Bool_t bValue=kFALSE)
vector< float > fIsolationInRingsChargedAll
void setApplyDzDxyVeto(Bool_t bValue=kTRUE)
vector< float > fIsolationInRings
void setNumbersOfRings(int iValue=1)
void setDeltaRVetoBarrelNeutrals(float fValue=-1.0)
math::XYZVector vtxWRTCandidate
void setDeltaRVetoEndcapNeutrals(float fValue=-1.0)
void setDeltaRVetoBarrelCharged(float fValue=-1.0)
void setDeltaRVetoEndcapPhotons(float fValue=-1.0)
void setApplyPFPUVeto(Bool_t bValue=kFALSE)
vector< float > getIsolationInRingsCharged()
void setRectangleVetoBarrel(Bool_t bValue=kTRUE)
void setRectangleDeltaEtaVetoBarrelNeutrals(float fValue=-1.0)
vector< float > fIsolationInRingsNeutral
void setUseCrystalSize(Bool_t bValue=kFALSE)
void setParticleType(int iValue)
vector< float > getIsolationInRingsNeutral()
void setRectangleDeltaPhiVetoBarrelCharged(float fValue=-1.0)
void setRectangleDeltaPhiVetoEndcapNeutrals(float fValue=-1.0)
void setRectangleDeltaPhiVetoBarrelNeutrals(float fValue=-1.0)
void setRectangleDeltaEtaVetoBarrelCharged(float fValue=-1.0)
vector< float > getIsolationInRingsChargedAll()
void setDeltaRVetoBarrelPhotons(float fValue=-1.0)
void setRectangleDeltaPhiVetoEndcapCharged(float fValue=-1.0)
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
void setRectangleDeltaEtaVetoEndcapCharged(float fValue=-1.0)
void setRectangleDeltaEtaVetoEndcapPhotons(float fValue=-1.0)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
void setDeltaRVetoEndcapCharged(float fValue=-1.0)
Bool_t isInitialized() const
void setApplyVeto(Bool_t bValue=kTRUE)
void setRingSize(float fValue=0.4)
vector< float > getIsolationInRingsPhoton()
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
void setRectangleDeltaPhiVetoEndcapPhotons(float fValue=-1.0)
void setRectangleDeltaEtaVetoBarrelPhotons(float fValue=-1.0)
void setRectangleDeltaEtaVetoEndcapNeutrals(float fValue=-1.0)
void setConeSize(float fValue=0.4)
vector< float > fIsolationInRingsPhoton