CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
SimplePhotonAnalyzer Class Reference

#include <SimplePhotonAnalyzer.h>

Inheritance diagram for SimplePhotonAnalyzer:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob ()
 
virtual void endJob ()
 
 SimplePhotonAnalyzer (const edm::ParameterSet &)
 
 ~SimplePhotonAnalyzer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

float etaTransformation (float a, float b)
 

Private Attributes

edm::InputTag barrelEcalHits_
 
TProfile * effEta_
 
TProfile * effPhi_
 
edm::InputTag endcapEcalHits_
 
TH1F * h1_deltaEta_
 
TH1F * h1_deltaEtaSC_
 
TH1F * h1_deltaPhi_
 
TH1F * h1_pho_E_
 
TH1F * h1_pho_ecalIsoBarrel_
 
TH1F * h1_pho_ecalIsoEndcap_
 
TH1F * h1_pho_Et_
 
TH1F * h1_pho_Eta_
 
TH1F * h1_pho_hcalIsoBarrel_
 
TH1F * h1_pho_hcalIsoEndcap_
 
TH1F * h1_pho_hOverEBarrel_
 
TH1F * h1_pho_hOverEEndcap_
 
TH1F * h1_pho_Phi_
 
TH1F * h1_pho_R9Barrel_
 
TH1F * h1_pho_R9Endcap_
 
TH1F * h1_pho_sigmaIetaIetaBarrel_
 
TH1F * h1_pho_sigmaIetaIetaEndcap_
 
TH1F * h1_pho_trkIsoBarrel_
 
TH1F * h1_pho_trkIsoEndcap_
 
TH1F * h1_recEoverTrueEBarrel_
 
TH1F * h1_recEoverTrueEEndcap_
 
TH1F * h1_scEta_
 
std::string mcCollection_
 
std::string mcProducer_
 
std::string photonCollection_
 
std::string photonCollectionProducer_
 
float sample_
 
edm::ESHandle< CaloTopologytheCaloTopo_
 
std::string vertexProducer_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Description: Get Photon collection from the event and make very basic histos

Date:
2009/03/06 12:42:16
Revision:
1.10
Author
Nancy Marinelli, U. of Notre Dame, US

Definition at line 36 of file SimplePhotonAnalyzer.h.

Constructor & Destructor Documentation

SimplePhotonAnalyzer::SimplePhotonAnalyzer ( const edm::ParameterSet ps)
explicit

Definition at line 30 of file SimplePhotonAnalyzer.cc.

References barrelEcalHits_, endcapEcalHits_, edm::ParameterSet::getParameter(), mcProducer_, photonCollection_, photonCollectionProducer_, sample_, and vertexProducer_.

32 {
33  photonCollectionProducer_ = ps.getParameter<std::string>("phoProducer");
34  photonCollection_ = ps.getParameter<std::string>("photonCollection");
35 
36  barrelEcalHits_ = ps.getParameter<edm::InputTag>("barrelEcalHits");
37  endcapEcalHits_ = ps.getParameter<edm::InputTag>("endcapEcalHits");
38 
39 
40 
41  mcProducer_ = ps.getParameter<std::string>("mcProducer");
42  //mcCollection_ = ps.getParameter<std::string>("mcCollection");
43  vertexProducer_ = ps.getParameter<std::string>("primaryVertexProducer");
44  sample_ = ps.getParameter<int>("sample");
45 
46 
47 }
T getParameter(std::string const &) const
std::string photonCollectionProducer_
SimplePhotonAnalyzer::~SimplePhotonAnalyzer ( )

Definition at line 51 of file SimplePhotonAnalyzer.cc.

53 {
54 
55 
56 }

Member Function Documentation

void SimplePhotonAnalyzer::analyze ( const edm::Event evt,
const edm::EventSetup es 
)
virtual

Match reconstructed photon candidates with the nearest generated photonPho;

Plot kinematic disctributions for matched photons

Implements edm::EDAnalyzer.

Definition at line 148 of file SimplePhotonAnalyzer.cc.

References delta, Geom::deltaPhi(), reco::Photon::ecalRecHitSumEtConeDR04(), effEta_, effPhi_, reco::LeafCandidate::energy(), reco::LeafCandidate::et(), eta(), reco::LeafCandidate::eta(), etaTransformation(), configurableAnalysis::GenParticle, edm::EventSetup::get(), edm::Event::getByLabel(), h1_deltaEta_, h1_deltaEtaSC_, h1_deltaPhi_, h1_pho_E_, h1_pho_ecalIsoBarrel_, h1_pho_ecalIsoEndcap_, h1_pho_Et_, h1_pho_Eta_, h1_pho_hcalIsoBarrel_, h1_pho_hcalIsoEndcap_, h1_pho_hOverEBarrel_, h1_pho_hOverEEndcap_, h1_pho_Phi_, h1_pho_R9Barrel_, h1_pho_R9Endcap_, h1_pho_sigmaIetaIetaBarrel_, h1_pho_sigmaIetaIetaEndcap_, h1_pho_trkIsoBarrel_, h1_pho_trkIsoEndcap_, h1_recEoverTrueEBarrel_, h1_recEoverTrueEEndcap_, h1_scEta_, reco::Photon::hadronicOverEm(), reco::Photon::hcalTowerSumEtConeDR04(), edm::EventBase::id(), getHLTprescales::index, mcProducer_, L1TEmulatorMonitor_cff::p, parents, reco::LeafCandidate::phi(), photonCollection_, photonCollectionProducer_, pi, funct::pow(), reco::Photon::r9(), reco::Photon::sigmaIetaIeta(), mathSSE::sqrt(), reco::Photon::superCluster(), theCaloTopo_, and reco::Photon::trkSumPtSolidConeDR04().

148  {
149 //========================================================================
150 
151  using namespace edm; // needed for all fwk related classes
152  edm::LogInfo("PhotonAnalyzer") << "Analyzing event number: " << evt.id() << "\n";
153 
154 
155  // get the calo topology from the event setup:
156  edm::ESHandle<CaloTopology> pTopology;
158 
159 
160 
161  // Get the corrected photon collection (set in the configuration) which also contains infos about conversions
162 
163  Handle<reco::PhotonCollection> photonHandle;
165  const reco::PhotonCollection photonCollection = *(photonHandle.product());
166 
167  Handle< HepMCProduct > hepProd ;
168  evt.getByLabel( mcProducer_.c_str(), hepProd ) ;
169  const HepMC::GenEvent * myGenEvent = hepProd->GetEvent();
170 
171 
172  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
173  if ( !( (*p)->pdg_id() == 22 && (*p)->status()==1 ) ) continue;
174 
175  // single primary photons or photons from Higgs or RS Graviton
176  HepMC::GenParticle* mother = 0;
177  if ( (*p)->production_vertex() ) {
178  if ( (*p)->production_vertex()->particles_begin(HepMC::parents) !=
179  (*p)->production_vertex()->particles_end(HepMC::parents))
180  mother = *((*p)->production_vertex()->particles_begin(HepMC::parents));
181  }
182  if ( ((mother == 0) || ((mother != 0) && (mother->pdg_id() == 25))
183  || ((mother != 0) && (mother->pdg_id() == 22)))) {
184 
185  float minDelta=10000.;
186  std::vector<reco::Photon> localPhotons;
187  int index=0;
188  int iMatch=-1;
189 
190  float phiPho=(*p)->momentum().phi();
191  float etaPho=(*p)->momentum().eta();
192  etaPho = etaTransformation(etaPho, (*p)->production_vertex()->position().z()/10. );
193 
194  bool matched=false;
195  // loop Photon candidates
196  for( reco::PhotonCollection::const_iterator iPho = photonCollection.begin(); iPho != photonCollection.end(); iPho++) {
197  reco::Photon localPho = reco::Photon(*iPho);
198  localPhotons.push_back(localPho);
199 
201  float phiClu=localPho.phi();
202  float etaClu=localPho.eta();
203  float deltaPhi = phiClu-phiPho;
204  float deltaEta = etaClu-etaPho;
205 
206  if ( deltaPhi > pi ) deltaPhi -= twopi;
207  if ( deltaPhi < -pi) deltaPhi += twopi;
208  deltaPhi=std::pow(deltaPhi,2);
209  deltaEta=std::pow(deltaEta,2);
210  float delta = sqrt( deltaPhi+deltaEta);
211  if ( delta<0.1 && delta < minDelta ) {
212  minDelta=delta;
213  iMatch=index;
214 
215  }
216  index++;
217  } // End loop over photons
218 
219  double wt=0.;
220  if ( iMatch>-1 ) matched = true;
221 
223  if (matched ) {
224  wt=1.;
225 
226  effEta_ ->Fill ( etaPho, wt);
227  effPhi_ ->Fill ( phiPho, wt);
228  reco::Photon matchingPho = localPhotons[iMatch];
229 
230  bool phoIsInBarrel=false;
231  bool phoIsInEndcap=false;
232  if ( fabs(matchingPho.superCluster()->position().eta() ) < 1.479 ) {
233  phoIsInBarrel=true;
234  } else {
235  phoIsInEndcap=true;
236  }
237  edm::Handle<EcalRecHitCollection> ecalRecHitHandle;
238 
239 
240  h1_scEta_->Fill( matchingPho.superCluster()->position().eta() );
241  float trueEta= (*p)->momentum().eta() ;
242  trueEta = etaTransformation(trueEta, (*p)->production_vertex()->position().z()/10. );
243  h1_deltaEtaSC_ -> Fill( localPhotons[iMatch].superCluster()->eta()- trueEta );
244 
245  float photonE = matchingPho.energy();
246  float photonEt= matchingPho.et() ;
247  float photonEta= matchingPho.eta() ;
248  float photonPhi= matchingPho.phi() ;
249 
250  float r9 = matchingPho.r9();
251  float sigmaIetaIeta = matchingPho.sigmaIetaIeta();
252  float hOverE = matchingPho.hadronicOverEm();
253  float ecalIso = matchingPho.ecalRecHitSumEtConeDR04();
254  float hcalIso = matchingPho.hcalTowerSumEtConeDR04();
255  float trkIso = matchingPho.trkSumPtSolidConeDR04();
256 
257 
258  h1_pho_E_->Fill( photonE );
259  h1_pho_Et_->Fill( photonEt );
260  h1_pho_Eta_->Fill( photonEta );
261  h1_pho_Phi_->Fill( photonPhi );
262 
263  h1_deltaEta_ -> Fill( photonEta - (*p)->momentum().eta() );
264  h1_deltaPhi_ -> Fill( photonPhi - (*p)->momentum().phi() );
265 
266  if ( phoIsInBarrel ) {
267  h1_recEoverTrueEBarrel_ -> Fill ( photonE / (*p)->momentum().e() );
268  h1_pho_R9Barrel_->Fill( r9 );
269  h1_pho_sigmaIetaIetaBarrel_->Fill ( sigmaIetaIeta );
270  h1_pho_hOverEBarrel_-> Fill ( hOverE );
271  h1_pho_ecalIsoBarrel_-> Fill ( ecalIso );
272  h1_pho_hcalIsoBarrel_-> Fill ( hcalIso );
273  h1_pho_trkIsoBarrel_-> Fill ( trkIso );
274 
275  } else {
276  h1_recEoverTrueEEndcap_ -> Fill ( photonE / (*p)->momentum().e() );
277  h1_pho_R9Endcap_->Fill( r9 );
278  h1_pho_sigmaIetaIetaEndcap_->Fill ( sigmaIetaIeta );
279  h1_pho_hOverEEndcap_-> Fill ( hOverE );
280  h1_pho_ecalIsoEndcap_-> Fill ( ecalIso );
281  h1_pho_hcalIsoEndcap_-> Fill ( hcalIso );
282  h1_pho_trkIsoEndcap_-> Fill ( trkIso );
283 
284 
285  }
286 
287 
288  } // reco photon matching MC truth
289 
290 
291 
292 
293  } // End loop over MC particles
294 
295  }
296 
297 
298 }
dbl * delta
Definition: mlp_gen.cc:36
float hcalTowerSumEtConeDR04() const
Hcal isolation sum.
Definition: Photon.h:245
TPRegexp parents
Definition: eve_filter.cc:24
float etaTransformation(float a, float b)
reco::SuperClusterRef superCluster() const
Retrieve photonCore attributes.
Definition: Photon.h:54
virtual double et() const
transverse energy
double deltaPhi(float phi1, float phi2)
Definition: VectorUtil.h:30
float trkSumPtSolidConeDR04() const
Definition: Photon.h:251
float ecalRecHitSumEtConeDR04() const
Definition: Photon.h:243
T eta() const
virtual double eta() const
momentum pseudorapidity
virtual double energy() const
energy
edm::ESHandle< CaloTopology > theCaloTopo_
T sqrt(T t)
Definition: SSEVec.h:28
float sigmaIetaIeta() const
Definition: Photon.h:167
float hadronicOverEm() const
the total hadronic over electromagnetic fraction
Definition: Photon.h:154
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
const T & get() const
Definition: EventSetup.h:55
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
edm::EventID id() const
Definition: EventBase.h:56
float r9() const
Definition: Photon.h:170
double pi
virtual double phi() const
momentum azimuthal angle
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::string photonCollectionProducer_
void SimplePhotonAnalyzer::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 60 of file SimplePhotonAnalyzer.cc.

References dPhi(), effEta_, effPhi_, h1_deltaEta_, h1_deltaEtaSC_, h1_deltaPhi_, h1_pho_E_, h1_pho_ecalIsoBarrel_, h1_pho_ecalIsoEndcap_, h1_pho_Et_, h1_pho_Eta_, h1_pho_hcalIsoBarrel_, h1_pho_hcalIsoEndcap_, h1_pho_hOverEBarrel_, h1_pho_hOverEEndcap_, h1_pho_Phi_, h1_pho_R9Barrel_, h1_pho_R9Endcap_, h1_pho_sigmaIetaIetaBarrel_, h1_pho_sigmaIetaIetaEndcap_, h1_pho_trkIsoBarrel_, h1_pho_trkIsoEndcap_, h1_recEoverTrueEBarrel_, h1_recEoverTrueEEndcap_, h1_scEta_, ExpressReco_HICollisions_FallBack::hiRes, TFileDirectory::make(), and sample_.

60  {
61 //========================================================================
62 
64 
65  float hiE=0;
66  float loE=0;
67  float hiEt=0;
68  float loEt=0;
69  float dPhi=0;
70  float loRes=0;
71  float hiRes=0;
72  if ( sample_ ==1 ) {
73  loE=0.;
74  hiE=30.;
75  loEt=0.;
76  hiEt=30.;
77  dPhi=0.2;
78  loRes=0.;
79  hiRes=1.2;
80  } else if ( sample_ ==2 ) {
81  loE=0.;
82  hiE=200.;
83  loEt=0.;
84  hiEt=50.;
85  dPhi=0.05;
86  loRes=0.7;
87  hiRes=1.2;
88  } else if ( sample_ ==3 ) {
89  loE=0.;
90  hiE=500.;
91  loEt=0.;
92  hiEt=500.;
93  dPhi=0.05;
94  loRes=0.7;
95  hiRes=1.2;
96  } else if (sample_==4) {
97  loE=0.;
98  hiE=6000.;
99  loEt=0.;
100  hiEt=1200.;
101  dPhi=0.05;
102  loRes=0.7;
103  hiRes=1.2;
104  }
105 
106 
107  effEta_ = fs->make<TProfile> ("effEta"," Photon reconstruction efficiency",50,-2.5,2.5);
108  effPhi_ = fs->make<TProfile> ("effPhi"," Photon reconstruction efficiency",80, -3.14, 3.14);
109 
110  h1_deltaEta_ = fs->make<TH1F>("deltaEta"," Reco photon Eta minus Generated photon Eta ",100,-0.2, 0.2);
111  h1_deltaPhi_ = fs->make<TH1F>("deltaPhi","Reco photon Phi minus Generated photon Phi ",100,-dPhi, dPhi);
112  h1_pho_Eta_ = fs->make<TH1F>("phoEta","Photon Eta ",40,-3., 3.);
113  h1_pho_Phi_ = fs->make<TH1F>("phoPhi","Photon Phi ",40,-3.14, 3.14);
114  h1_pho_E_ = fs->make<TH1F>("phoE","Photon Energy ",100,loE,hiE);
115  h1_pho_Et_ = fs->make<TH1F>("phoEt","Photon Et ",100,loEt,hiEt);
116 
117  h1_scEta_ = fs->make<TH1F>("scEta"," SC Eta ",40,-3., 3.);
118  h1_deltaEtaSC_ = fs->make<TH1F>("deltaEtaSC"," SC Eta minus Generated photon Eta ",100,-0.02, 0.02);
119 
120  //
121  h1_recEoverTrueEBarrel_ = fs->make<TH1F>("recEoverTrueEBarrel"," Reco photon Energy over Generated photon Energy: Barrel ",100,loRes, hiRes);
122  h1_recEoverTrueEEndcap_ = fs->make<TH1F>("recEoverTrueEEndcap"," Reco photon Energy over Generated photon Energy: Endcap ",100,loRes, hiRes);
123 
124  //
125 
126  h1_pho_R9Barrel_ = fs->make<TH1F>("phoR9Barrel","Photon 3x3 energy / SuperCluster energy : Barrel ",100,0.,1.2);
127  h1_pho_R9Endcap_ = fs->make<TH1F>("phoR9Endcap","Photon 3x3 energy / SuperCluster energy : Endcap ",100,0.,1.2);
128  h1_pho_sigmaIetaIetaBarrel_ = fs->make<TH1F>("sigmaIetaIetaBarrel", "sigmaIetaIeta: Barrel",100,0., 0.05) ;
129  h1_pho_sigmaIetaIetaEndcap_ = fs->make<TH1F>("sigmaIetaIetaEndcap" , "sigmaIetaIeta: Endcap",100,0., 0.1) ;
130  h1_pho_hOverEBarrel_ = fs->make<TH1F>("hOverEBarrel", "H/E: Barrel",100,0., 0.1) ;
131  h1_pho_hOverEEndcap_ = fs->make<TH1F>("hOverEEndcap", "H/E: Endcap",100,0., 0.1) ;
132  h1_pho_ecalIsoBarrel_ = fs->make<TH1F>("ecalIsolBarrel", "isolation et sum in Ecal: Barrel",100,0., 100.) ;
133  h1_pho_ecalIsoEndcap_ = fs->make<TH1F>("ecalIsolEndcap", "isolation et sum in Ecal: Endcap",100,0., 100.) ;
134  h1_pho_hcalIsoBarrel_ = fs->make<TH1F>("hcalIsolBarrel", "isolation et sum in Hcal: Barrel",100,0., 100.) ;
135  h1_pho_hcalIsoEndcap_ = fs->make<TH1F>("hcalIsolEndcap", "isolation et sum in Hcal: Endcap",100,0., 100.) ;
136  h1_pho_trkIsoBarrel_ = fs->make<TH1F>("trkIsolBarrel", "isolation pt sum in the tracker: Barrel",100,0., 100.) ;
137  h1_pho_trkIsoEndcap_ = fs->make<TH1F>("trkIsolEndcap", "isolation pt sum in the tracker: Endcap",100,0., 100.) ;
138 
139 
140 
141 
142 
143 }
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
T * make() const
make new ROOT object
void SimplePhotonAnalyzer::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 341 of file SimplePhotonAnalyzer.cc.

341  {
342 //========================================================================
343 
344 
345 
346 }
float SimplePhotonAnalyzer::etaTransformation ( float  a,
float  b 
)
private

Definition at line 301 of file SimplePhotonAnalyzer.cc.

References ETA, etaBarrelEndcap, funct::log(), PI, R_ECAL, funct::tan(), and Z_Endcap.

Referenced by analyze().

301  {
302 
303  //---Definitions
304  const float PI = 3.1415927;
305  //UNUSED const float TWOPI = 2.0*PI;
306 
307  //---Definitions for ECAL
308  const float R_ECAL = 136.5;
309  const float Z_Endcap = 328.0;
310  const float etaBarrelEndcap = 1.479;
311 
312  //---ETA correction
313 
314  float Theta = 0.0 ;
315  float ZEcal = R_ECAL*sinh(EtaParticle)+Zvertex;
316 
317  if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
318  if(Theta<0.0) Theta = Theta+PI ;
319  float ETA = - log(tan(0.5*Theta));
320 
321  if( fabs(ETA) > etaBarrelEndcap )
322  {
323  float Zend = Z_Endcap ;
324  if(EtaParticle<0.0 ) Zend = -Zend ;
325  float Zlen = Zend - Zvertex ;
326  float RR = Zlen/sinh(EtaParticle);
327  Theta = atan(RR/Zend);
328  if(Theta<0.0) Theta = Theta+PI ;
329  ETA = - log(tan(0.5*Theta));
330  }
331  //---Return the result
332  return ETA;
333  //---end
334 }
#define PI
#define ETA
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Log< T >::type log(const T &t)
Definition: Log.h:22
static const float etaBarrelEndcap
static const float Z_Endcap
static const float R_ECAL

Member Data Documentation

edm::InputTag SimplePhotonAnalyzer::barrelEcalHits_
private

Definition at line 54 of file SimplePhotonAnalyzer.h.

Referenced by SimplePhotonAnalyzer().

TProfile* SimplePhotonAnalyzer::effEta_
private

Definition at line 63 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TProfile* SimplePhotonAnalyzer::effPhi_
private

Definition at line 64 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

edm::InputTag SimplePhotonAnalyzer::endcapEcalHits_
private

Definition at line 55 of file SimplePhotonAnalyzer.h.

Referenced by SimplePhotonAnalyzer().

TH1F* SimplePhotonAnalyzer::h1_deltaEta_
private

Definition at line 89 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_deltaEtaSC_
private

Definition at line 67 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_deltaPhi_
private

Definition at line 90 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_E_
private

Definition at line 68 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_ecalIsoBarrel_
private

Definition at line 78 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_ecalIsoEndcap_
private

Definition at line 79 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_Et_
private

Definition at line 69 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_Eta_
private

Definition at line 70 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_hcalIsoBarrel_
private

Definition at line 80 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_hcalIsoEndcap_
private

Definition at line 81 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_hOverEBarrel_
private

Definition at line 76 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_hOverEEndcap_
private

Definition at line 77 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_Phi_
private

Definition at line 71 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_R9Barrel_
private

Definition at line 72 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_R9Endcap_
private

Definition at line 73 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_sigmaIetaIetaBarrel_
private

Definition at line 74 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_sigmaIetaIetaEndcap_
private

Definition at line 75 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_trkIsoBarrel_
private

Definition at line 82 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_pho_trkIsoEndcap_
private

Definition at line 83 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_recEoverTrueEBarrel_
private

Definition at line 87 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_recEoverTrueEEndcap_
private

Definition at line 88 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* SimplePhotonAnalyzer::h1_scEta_
private

Definition at line 66 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and beginJob().

std::string SimplePhotonAnalyzer::mcCollection_
private

Definition at line 50 of file SimplePhotonAnalyzer.h.

std::string SimplePhotonAnalyzer::mcProducer_
private

Definition at line 49 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and SimplePhotonAnalyzer().

std::string SimplePhotonAnalyzer::photonCollection_
private

Definition at line 52 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and SimplePhotonAnalyzer().

std::string SimplePhotonAnalyzer::photonCollectionProducer_
private

Definition at line 51 of file SimplePhotonAnalyzer.h.

Referenced by analyze(), and SimplePhotonAnalyzer().

float SimplePhotonAnalyzer::sample_
private

Definition at line 60 of file SimplePhotonAnalyzer.h.

Referenced by beginJob(), and SimplePhotonAnalyzer().

edm::ESHandle<CaloTopology> SimplePhotonAnalyzer::theCaloTopo_
private

Definition at line 57 of file SimplePhotonAnalyzer.h.

Referenced by analyze().

std::string SimplePhotonAnalyzer::vertexProducer_
private

Definition at line 59 of file SimplePhotonAnalyzer.h.

Referenced by SimplePhotonAnalyzer().