CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes | Private Attributes
HiGammaJetSignalDef Class Reference

#include <HiPhotonType.h>

Public Member Functions

double getDeltaPhi (const reco::Candidate &track1, const reco::Candidate &track2)
 
double getDeltaR (const reco::Candidate &track1, const reco::Candidate &track2)
 
 HiGammaJetSignalDef ()
 
 HiGammaJetSignalDef (const reco::GenParticleCollection *sigPartic)
 
bool IsIsolated (const reco::GenParticle &pp)
 
bool IsIsolatedJP (const reco::GenParticle &pp)
 
bool IsIsolatedPP (const reco::GenParticle &pp)
 

Public Attributes

double PI
 

Private Attributes

const reco::GenParticleCollectionfSigParticles
 

Detailed Description

Definition at line 23 of file HiPhotonType.h.

Constructor & Destructor Documentation

HiGammaJetSignalDef::HiGammaJetSignalDef ( )

Definition at line 63 of file HiPhotonType.cc.

References PI.

63  :
64  fSigParticles(0)
65 { PI = 3.141592653589793238462643383279502884197169399375105820974945;
66 }
const reco::GenParticleCollection * fSigParticles
Definition: HiPhotonType.h:40
HiGammaJetSignalDef::HiGammaJetSignalDef ( const reco::GenParticleCollection sigPartic)

Definition at line 67 of file HiPhotonType.cc.

References fSigParticles, and PI.

67  :
68  fSigParticles(0)
69 {
70  using namespace std;
71 
72  fSigParticles = sigParticles;
73  PI = 3.141592653589793238462643383279502884197169399375105820974945;
74 
75 }
const reco::GenParticleCollection * fSigParticles
Definition: HiPhotonType.h:40

Member Function Documentation

double HiGammaJetSignalDef::getDeltaPhi ( const reco::Candidate track1,
const reco::Candidate track2 
)

Definition at line 230 of file HiPhotonType.cc.

References dPhi(), reco::Candidate::phi(), and PI.

231 {
232  double dPhi = track1.phi() - track2.phi();
233 
234  while(dPhi >= PI) dPhi -= (2.0*PI);
235  while(dPhi < (-1.0*PI)) dPhi += (2.0*PI);
236 
237  return dPhi;
238 }
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
virtual double phi() const =0
momentum azimuthal angle
double HiGammaJetSignalDef::getDeltaR ( const reco::Candidate track1,
const reco::Candidate track2 
)

Definition at line 219 of file HiPhotonType.cc.

References dPhi(), reco::Candidate::eta(), reco::Candidate::phi(), PI, and mathSSE::sqrt().

Referenced by IsIsolated(), IsIsolatedJP(), and IsIsolatedPP().

220 {
221  double dEta = track1.eta() - track2.eta();
222  double dPhi = track1.phi() - track2.phi();
223 
224  while(dPhi >= PI) dPhi -= (2.0*PI);
225  while(dPhi < (-1.0*PI)) dPhi += (2.0*PI);
226 
227  return sqrt(dEta*dEta+dPhi*dPhi);
228 }
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
T sqrt(T t)
Definition: SSEVec.h:18
virtual double eta() const =0
momentum pseudorapidity
virtual double phi() const =0
momentum azimuthal angle
bool HiGammaJetSignalDef::IsIsolated ( const reco::GenParticle pp)

Definition at line 77 of file HiPhotonType.cc.

References funct::abs(), reco::GenParticle::collisionId(), stringResolutionProvider_cfi::et, reco::LeafCandidate::et(), fSigParticles, getDeltaR(), mps_fire::i, createfilelist::int, AlCaHLTBitMon_ParallelJobs::p, reco::LeafCandidate::pdgId(), EnergyCorrector::pt, reco::LeafCandidate::pt(), and reco::LeafCandidate::status().

78 {
79  using namespace std;
80  using namespace edm;
81  using namespace reco;
82  // Check if a given particle is isolated.
83 
84  double EtCone = 0;
85  double EtPhoton = 0;
86  double PtMaxHadron = 0;
87  double cone = 0.5;
88  const int maxindex = (int)fSigParticles->size();
89  for(int i=0; i < maxindex ; ++i) {
90 
91  const GenParticle &p=(*fSigParticles)[i];
92 
93  if(p.status()!=1)
94  continue;
95  if (p.collisionId() != pp.collisionId())
96  continue;
97 
98  int apid= abs(p.pdgId());
99  if(apid>11 && apid<20)
100  continue; //get rid of muons and neutrinos
101 
102  if(getDeltaR(p,pp)>cone)
103  continue;
104 
105  double et=p.et();
106  double pt = p.pt();
107  EtCone+=et;
108  bool isHadron = false;
109  if ( fabs(pp.pdgId())>100 && fabs(pp.pdgId()) != 310)
110  isHadron = true;
111 
112  if(apid>100 && apid!=310 && pt>PtMaxHadron)
113  {
114  if ( (isHadron == true) && (fabs(pp.pt()-pt)<0.001) && (pp.pdgId()==p.pdgId()) )
115  continue;
116  else
117  PtMaxHadron=pt;
118  }
119 
120  }
121  EtPhoton = pp.et();
122 
123  // isolation cuts
124  if(EtCone-EtPhoton > 5+EtPhoton/20)
125  return kFALSE;
126  if(PtMaxHadron > 4.5+EtPhoton/40)
127  return kFALSE;
128 
129  return kTRUE;
130 
131 }
double getDeltaR(const reco::Candidate &track1, const reco::Candidate &track2)
virtual double pt() const final
transverse momentum
virtual int status() const final
status word
int collisionId() const
Definition: GenParticle.h:39
virtual double et() const final
transverse energy
virtual int pdgId() const final
PDG identifier.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
et
define resolution functions of each parameter
fixed size matrix
HLT enums.
const reco::GenParticleCollection * fSigParticles
Definition: HiPhotonType.h:40
bool HiGammaJetSignalDef::IsIsolatedJP ( const reco::GenParticle pp)

Definition at line 179 of file HiPhotonType.cc.

References reco::GenParticle::collisionId(), stringResolutionProvider_cfi::et, reco::LeafCandidate::et(), fSigParticles, getDeltaR(), mps_fire::i, createfilelist::int, AlCaHLTBitMon_ParallelJobs::p, and reco::LeafCandidate::status().

180 {
181  using namespace std;
182  using namespace edm;
183  using namespace reco;
184 
185  // Check if a given particle is isolated.
186 
187  double EtCone = 0;
188  double EtPhoton = 0;
189  double cone = 0.4;
190 
191  const int maxindex = (int)fSigParticles->size();
192  for(int i=0; i < maxindex ; ++i) {
193 
194  const GenParticle &p=(*fSigParticles)[i];
195 
196  if(p.status()!=1)
197  continue;
198  if (p.collisionId() != pp.collisionId())
199  continue;
200 
201  if(getDeltaR(p,pp)>cone)
202  continue;
203 
204  double et=p.et();
205  // double pt = p.pt();
206  EtCone+=et;
207 
208  }
209  EtPhoton = pp.et();
210 
211  // isolation cuts
212 
213  if(EtCone-EtPhoton > 5)
214  return kFALSE;
215  return kTRUE;
216 
217 }
double getDeltaR(const reco::Candidate &track1, const reco::Candidate &track2)
virtual int status() const final
status word
int collisionId() const
Definition: GenParticle.h:39
virtual double et() const final
transverse energy
et
define resolution functions of each parameter
fixed size matrix
HLT enums.
const reco::GenParticleCollection * fSigParticles
Definition: HiPhotonType.h:40
bool HiGammaJetSignalDef::IsIsolatedPP ( const reco::GenParticle pp)

Definition at line 133 of file HiPhotonType.cc.

References reco::GenParticle::collisionId(), stringResolutionProvider_cfi::et, reco::LeafCandidate::et(), fSigParticles, getDeltaR(), mps_fire::i, createfilelist::int, AlCaHLTBitMon_ParallelJobs::p, and reco::LeafCandidate::status().

134 {
135  using namespace std;
136  using namespace edm;
137  using namespace reco;
138 
139  // Check if a given particle is isolated.
140 
141  double EtCone = 0;
142  double EtPhoton = 0;
143  double cone = 0.4;
144 
145  const int maxindex = (int)fSigParticles->size();
146  for(int i=0; i < maxindex ; ++i) {
147 
148  const GenParticle &p=(*fSigParticles)[i];
149 
150  if(p.status()!=1)
151  continue;
152  if (p.collisionId() != pp.collisionId())
153  continue;
154 
155  //int apid= abs(p.pdgId());
156  // if(apid>11 && apid<20)
157  // continue; //get rid of muons and neutrinos
158 
159  if(getDeltaR(p,pp)>cone)
160  continue;
161  double et=p.et();
162  // double pt = p.pt();
163  EtCone+=et;
164 
165  }
166 
167  EtPhoton = pp.et();
168 
169  // isolation cuts
170 
171  if(EtCone-EtPhoton > 5)
172  return kFALSE;
173  //if(PtMaxHadron > 4.5+EtPhoton/40)
174  // return kFALSE;
175 
176  return kTRUE;
177 
178 }
double getDeltaR(const reco::Candidate &track1, const reco::Candidate &track2)
virtual int status() const final
status word
int collisionId() const
Definition: GenParticle.h:39
virtual double et() const final
transverse energy
et
define resolution functions of each parameter
fixed size matrix
HLT enums.
const reco::GenParticleCollection * fSigParticles
Definition: HiPhotonType.h:40

Member Data Documentation

const reco::GenParticleCollection* HiGammaJetSignalDef::fSigParticles
private

Definition at line 40 of file HiPhotonType.h.

Referenced by HiGammaJetSignalDef(), IsIsolated(), IsIsolatedJP(), and IsIsolatedPP().

double HiGammaJetSignalDef::PI

Definition at line 37 of file HiPhotonType.h.

Referenced by getDeltaPhi(), getDeltaR(), and HiGammaJetSignalDef().