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() [1/2]

HiGammaJetSignalDef::HiGammaJetSignalDef ( )

Definition at line 50 of file HiPhotonType.cc.

50  : fSigParticles(nullptr) {
51  PI = 3.141592653589793238462643383279502884197169399375105820974945;
52 }

◆ HiGammaJetSignalDef() [2/2]

HiGammaJetSignalDef::HiGammaJetSignalDef ( const reco::GenParticleCollection sigPartic)

Definition at line 53 of file HiPhotonType.cc.

53  : fSigParticles(nullptr) {
54  using namespace std;
55 
56  fSigParticles = sigParticles;
57  PI = 3.141592653589793238462643383279502884197169399375105820974945;
58 }

References fSigParticles.

Member Function Documentation

◆ getDeltaPhi()

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

Definition at line 201 of file HiPhotonType.cc.

201  {
202  double dPhi = track1.phi() - track2.phi();
203 
204  while (dPhi >= PI)
205  dPhi -= (2.0 * PI);
206  while (dPhi < (-1.0 * PI))
207  dPhi += (2.0 * PI);
208 
209  return dPhi;
210 }

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

◆ getDeltaR()

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

Definition at line 189 of file HiPhotonType.cc.

189  {
190  double dEta = track1.eta() - track2.eta();
191  double dPhi = track1.phi() - track2.phi();
192 
193  while (dPhi >= PI)
194  dPhi -= (2.0 * PI);
195  while (dPhi < (-1.0 * PI))
196  dPhi += (2.0 * PI);
197 
198  return sqrt(dEta * dEta + dPhi * dPhi);
199 }

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

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

◆ IsIsolated()

bool HiGammaJetSignalDef::IsIsolated ( const reco::GenParticle pp)

Definition at line 60 of file HiPhotonType.cc.

60  {
61  using namespace std;
62  using namespace edm;
63  using namespace reco;
64  // Check if a given particle is isolated.
65 
66  double EtCone = 0;
67  double EtPhoton = 0;
68  double PtMaxHadron = 0;
69  double cone = 0.5;
70  const int maxindex = (int)fSigParticles->size();
71  for (int i = 0; i < maxindex; ++i) {
72  const GenParticle &p = (*fSigParticles)[i];
73 
74  if (p.status() != 1)
75  continue;
76  if (p.collisionId() != pp.collisionId())
77  continue;
78 
79  int apid = abs(p.pdgId());
80  if (apid > 11 && apid < 20)
81  continue; //get rid of muons and neutrinos
82 
83  if (getDeltaR(p, pp) > cone)
84  continue;
85 
86  double et = p.et();
87  double pt = p.pt();
88  EtCone += et;
89  bool isHadron = false;
90  if (fabs(pp.pdgId()) > 100 && fabs(pp.pdgId()) != 310)
91  isHadron = true;
92 
93  if (apid > 100 && apid != 310 && pt > PtMaxHadron) {
94  if ((isHadron == true) && (fabs(pp.pt() - pt) < 0.001) && (pp.pdgId() == p.pdgId()))
95  continue;
96  else
97  PtMaxHadron = pt;
98  }
99  }
100  EtPhoton = pp.et();
101 
102  // isolation cuts
103  if (EtCone - EtPhoton > 5 + EtPhoton / 20)
104  return kFALSE;
105  if (PtMaxHadron > 4.5 + EtPhoton / 40)
106  return kFALSE;
107 
108  return kTRUE;
109 }

References funct::abs(), EgHLTOffHistBins_cfi::et, fSigParticles, getDeltaR(), mps_fire::i, createfilelist::int, AlCaHLTBitMon_ParallelJobs::p, createTree::pp, and DiDispStaMuonMonitor_cfi::pt.

◆ IsIsolatedJP()

bool HiGammaJetSignalDef::IsIsolatedJP ( const reco::GenParticle pp)

Definition at line 153 of file HiPhotonType.cc.

153  {
154  using namespace std;
155  using namespace edm;
156  using namespace reco;
157 
158  // Check if a given particle is isolated.
159 
160  double EtCone = 0;
161  double EtPhoton = 0;
162  double cone = 0.4;
163 
164  const int maxindex = (int)fSigParticles->size();
165  for (int i = 0; i < maxindex; ++i) {
166  const GenParticle &p = (*fSigParticles)[i];
167 
168  if (p.status() != 1)
169  continue;
170  if (p.collisionId() != pp.collisionId())
171  continue;
172 
173  if (getDeltaR(p, pp) > cone)
174  continue;
175 
176  double et = p.et();
177  // double pt = p.pt();
178  EtCone += et;
179  }
180  EtPhoton = pp.et();
181 
182  // isolation cuts
183 
184  if (EtCone - EtPhoton > 5)
185  return kFALSE;
186  return kTRUE;
187 }

References EgHLTOffHistBins_cfi::et, fSigParticles, getDeltaR(), mps_fire::i, createfilelist::int, AlCaHLTBitMon_ParallelJobs::p, and createTree::pp.

◆ IsIsolatedPP()

bool HiGammaJetSignalDef::IsIsolatedPP ( const reco::GenParticle pp)

Definition at line 111 of file HiPhotonType.cc.

111  {
112  using namespace std;
113  using namespace edm;
114  using namespace reco;
115 
116  // Check if a given particle is isolated.
117 
118  double EtCone = 0;
119  double EtPhoton = 0;
120  double cone = 0.4;
121 
122  const int maxindex = (int)fSigParticles->size();
123  for (int i = 0; i < maxindex; ++i) {
124  const GenParticle &p = (*fSigParticles)[i];
125 
126  if (p.status() != 1)
127  continue;
128  if (p.collisionId() != pp.collisionId())
129  continue;
130 
131  //int apid= abs(p.pdgId());
132  // if(apid>11 && apid<20)
133  // continue; //get rid of muons and neutrinos
134 
135  if (getDeltaR(p, pp) > cone)
136  continue;
137  double et = p.et();
138  // double pt = p.pt();
139  EtCone += et;
140  }
141 
142  EtPhoton = pp.et();
143 
144  // isolation cuts
145 
146  if (EtCone - EtPhoton > 5)
147  return kFALSE;
148  //if(PtMaxHadron > 4.5+EtPhoton/40)
149  // return kFALSE;
150 
151  return kTRUE;
152 }

References EgHLTOffHistBins_cfi::et, fSigParticles, getDeltaR(), mps_fire::i, createfilelist::int, AlCaHLTBitMon_ParallelJobs::p, and createTree::pp.

Member Data Documentation

◆ fSigParticles

const reco::GenParticleCollection* HiGammaJetSignalDef::fSigParticles
private

Definition at line 38 of file HiPhotonType.h.

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

◆ PI

double HiGammaJetSignalDef::PI

Definition at line 35 of file HiPhotonType.h.

Referenced by getDeltaPhi(), and getDeltaR().

PI
Definition: PayloadInspector.h:20
mps_fire.i
i
Definition: mps_fire.py:428
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
HiGammaJetSignalDef::fSigParticles
const reco::GenParticleCollection * fSigParticles
Definition: HiPhotonType.h:38
reco::Candidate::eta
virtual double eta() const =0
momentum pseudorapidity
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
HiGammaJetSignalDef::getDeltaR
double getDeltaR(const reco::Candidate &track1, const reco::Candidate &track2)
Definition: HiPhotonType.cc:189
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
HLT_FULL_cff.dPhi
dPhi
Definition: HLT_FULL_cff.py:13765
GenParticle
Definition: GenParticle.py:1
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
createfilelist.int
int
Definition: createfilelist.py:10
EgHLTOffHistBins_cfi.et
et
Definition: EgHLTOffHistBins_cfi.py:8
std
Definition: JetResolutionObject.h:76
HLT_FULL_cff.dEta
dEta
Definition: HLT_FULL_cff.py:13764
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
createTree.pp
pp
Definition: createTree.py:17
HiGammaJetSignalDef::PI
double PI
Definition: HiPhotonType.h:35
reco::Candidate::phi
virtual double phi() const =0
momentum azimuthal angle