#include <RecoEgamma/Examples/interface/SimplePhotonAnalyzer.h>
Description: Get Photon collection from the event and make very basic histos
Definition at line 32 of file SimplePhotonAnalyzer.h.
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_.
00032 { 00033 photonCollectionProducer_ = ps.getParameter<std::string>("phoProducer"); 00034 photonCollection_ = ps.getParameter<std::string>("photonCollection"); 00035 00036 barrelEcalHits_ = ps.getParameter<edm::InputTag>("barrelEcalHits"); 00037 endcapEcalHits_ = ps.getParameter<edm::InputTag>("endcapEcalHits"); 00038 00039 00040 00041 mcProducer_ = ps.getParameter<std::string>("mcProducer"); 00042 //mcCollection_ = ps.getParameter<std::string>("mcCollection"); 00043 vertexProducer_ = ps.getParameter<std::string>("primaryVertexProducer"); 00044 sample_ = ps.getParameter<int>("sample"); 00045 00046 00047 }
SimplePhotonAnalyzer::~SimplePhotonAnalyzer | ( | ) |
SimplePhotonAnalyzer::SimplePhotonAnalyzer | ( | const edm::ParameterSet & | ) | [explicit] |
SimplePhotonAnalyzer::~SimplePhotonAnalyzer | ( | ) |
virtual void SimplePhotonAnalyzer::analyze | ( | const edm::Event & | , | |
const edm::EventSetup & | ||||
) | [virtual] |
Implements edm::EDAnalyzer.
void SimplePhotonAnalyzer::analyze | ( | const edm::Event & | evt, | |
const edm::EventSetup & | es | |||
) | [virtual] |
Get the MC truth
Match reconstructed photon candidates with the nearest generated photonPho;
Plot kinematic disctributions for matched photons
Implements edm::EDAnalyzer.
Definition at line 140 of file SimplePhotonAnalyzer.cc.
References barrelEcalHits_, GenMuonPlsPt100GeV_cfg::cout, deltaPhi(), EcalClusterTools::e3x3(), EcalClusterTools::e5x5(), effEta_, effPhi_, endcapEcalHits_, lat::endl(), relval_parameters_module::energy, reco::Particle::eta(), eta, etaTransformation(), edm::EventSetup::get(), edm::Event::getByLabel(), h1_deltaEta_, h1_deltaEtaSC_, h1_deltaPhi_, h1_deltaPhiSC_, h1_e5x5_unconvBarrel_, h1_e5x5_unconvEndcap_, h1_ePho_convBarrel_, h1_ePho_convEndcap_, h1_pho_E_, h1_pho_Eta_, h1_pho_Phi_, h1_pho_R9Barrel_, h1_pho_R9Endcap_, h1_recEoverTrueEBarrel_, h1_recEoverTrueEEndcap_, h1_recESCoverTrueEBarrel_, h1_recESCoverTrueEEndcap_, h1_scE_, h1_scEt_, h1_scEta_, h1_scPhi_, edm::Event::id(), index, edm::Handle< T >::isValid(), edm::InputTag::label(), mcProducer_, p, phi, reco::Particle::phi(), configurableAnalysis::Photon, photonCollection_, photonCollectionProducer_, pi, funct::pow(), edm::ESHandle< T >::product(), edm::Handle< T >::product(), r9, funct::sqrt(), theCaloTopo_, and vertexProducer_.
00140 { 00141 //======================================================================== 00142 00143 using namespace edm; // needed for all fwk related classes 00144 edm::LogInfo("PhotonAnalyzer") << "Analyzing event number: " << evt.id() << "\n"; 00145 00146 00147 // get the calo topology from the event setup: 00148 edm::ESHandle<CaloTopology> pTopology; 00149 es.get<CaloTopologyRecord>().get(theCaloTopo_); 00150 const CaloTopology *topology = theCaloTopo_.product(); 00151 00152 00153 // Get the corrected photon collection (set in the configuration) which also contains infos about conversions 00154 00155 Handle<reco::PhotonCollection> photonHandle; 00156 evt.getByLabel(photonCollectionProducer_, photonCollection_ , photonHandle); 00157 const reco::PhotonCollection photonCollection = *(photonHandle.product()); 00158 00159 00160 // Get the primary event vertex 00161 Handle<reco::VertexCollection> vertexHandle; 00162 evt.getByLabel(vertexProducer_, vertexHandle); 00163 reco::VertexCollection vertexCollection = *(vertexHandle.product()); 00164 math::XYZPoint vtx(0.,0.,0.); 00165 if (vertexCollection.size()>0) { 00166 std::cout << " I am using the Primary vertex position " << std::endl; 00167 vtx = vertexCollection.begin()->position(); 00168 } 00169 00171 Handle< HepMCProduct > hepProd ; 00172 evt.getByLabel( mcProducer_.c_str(), hepProd ) ; 00173 const HepMC::GenEvent * myGenEvent = hepProd->GetEvent(); 00174 00175 for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) { 00176 if ( !( (*p)->pdg_id() == 22 && (*p)->status()==1 ) ) continue; 00177 00178 // single primary photons or photons from Higgs or RS Graviton 00179 HepMC::GenParticle* mother = 0; 00180 if ( (*p)->production_vertex() ) { 00181 if ( (*p)->production_vertex()->particles_begin(HepMC::parents) != 00182 (*p)->production_vertex()->particles_end(HepMC::parents)) 00183 mother = *((*p)->production_vertex()->particles_begin(HepMC::parents)); 00184 } 00185 if ( ((mother == 0) || ((mother != 0) && (mother->pdg_id() == 25)) 00186 || ((mother != 0) && (mother->pdg_id() == 22)))) { 00187 00188 float minDelta=10000.; 00189 std::vector<reco::Photon> localPhotons; 00190 int index=0; 00191 int iMatch=-1; 00192 00193 float phiPho=(*p)->momentum().phi(); 00194 float etaPho=(*p)->momentum().eta(); 00195 etaPho = etaTransformation(etaPho, (*p)->production_vertex()->position().z()/10. ); 00196 00197 00198 00199 // loop over corrected Photon candidates 00200 for( reco::PhotonCollection::const_iterator iPho = photonCollection.begin(); iPho != photonCollection.end(); iPho++) { 00201 00203 reco::Photon localPho = reco::Photon(*iPho); 00204 // localPho.setVertex(vtx); 00205 localPhotons.push_back(localPho); 00206 00208 float phiClu=localPho.phi(); 00209 float etaClu=localPho.eta(); 00210 float deltaPhi = phiClu-phiPho; 00211 float deltaEta = etaClu-etaPho; 00212 00213 if ( deltaPhi > pi ) deltaPhi -= twopi; 00214 if ( deltaPhi < -pi) deltaPhi += twopi; 00215 deltaPhi=pow(deltaPhi,2); 00216 deltaEta=pow(deltaEta,2); 00217 float delta = sqrt( deltaPhi+deltaEta); 00218 if ( delta<0.1 && delta < minDelta ) { 00219 minDelta=delta; 00220 iMatch=index; 00221 00222 } 00223 index++; 00224 } // End loop over uncorrected photons 00225 00226 double wt=0.; 00227 if ( iMatch>-1 ) wt=1.; 00228 00229 effEta_ ->Fill ( etaPho, wt); 00230 effPhi_ ->Fill ( phiPho, wt); 00231 00233 if (iMatch>-1) { 00234 00235 00236 bool phoIsInBarrel=false; 00237 bool phoIsInEndcap=false; 00238 if ( fabs(localPhotons[iMatch].superCluster()->position().eta() ) < 1.479 ) { 00239 phoIsInBarrel=true; 00240 } else { 00241 phoIsInEndcap=true; 00242 } 00243 edm::Handle<EcalRecHitCollection> ecalRecHitHandle; 00244 00245 00246 if ( phoIsInBarrel ) { 00247 // Get handle to rec hits ecal barrel 00248 evt.getByLabel(barrelEcalHits_, ecalRecHitHandle); 00249 if (!ecalRecHitHandle.isValid()) { 00250 edm::LogError("PhotonProducer") << "Error! Can't get the product "<<barrelEcalHits_.label(); 00251 return; 00252 } 00253 00254 } else if ( phoIsInEndcap ) { 00255 00256 // Get handle to rec hits ecal encap 00257 evt.getByLabel(endcapEcalHits_, ecalRecHitHandle); 00258 if (!ecalRecHitHandle.isValid()) { 00259 edm::LogError("PhotonProducer") << "Error! Can't get the product "<<endcapEcalHits_.label(); 00260 return; 00261 } 00262 00263 00264 } 00265 const EcalRecHitCollection ecalRecHitCollection = *(ecalRecHitHandle.product()); 00266 float e5x5= EcalClusterTools::e5x5( *( localPhotons[iMatch].superCluster()->seed()), &ecalRecHitCollection, &(*topology)); 00267 float e3x3= EcalClusterTools::e3x3( *( localPhotons[iMatch].superCluster()->seed() ), &ecalRecHitCollection, &(*topology)); 00268 float r9 =e3x3/( localPhotons[iMatch].superCluster()->rawEnergy()+ localPhotons[iMatch].superCluster()->preshowerEnergy()); 00269 00270 00271 00272 00273 h1_scE_->Fill( localPhotons[iMatch].superCluster()->energy() ); 00274 h1_scEt_->Fill( localPhotons[iMatch].superCluster()->energy()/cosh(localPhotons[iMatch].superCluster()->position().eta()) ); 00275 h1_scEta_->Fill( localPhotons[iMatch].superCluster()->position().eta() ); 00276 h1_scPhi_->Fill( localPhotons[iMatch].superCluster()->position().phi() ); 00277 00278 float trueEta= (*p)->momentum().eta() ; 00279 trueEta = etaTransformation(trueEta, (*p)->production_vertex()->position().z()/10. ); 00280 h1_deltaEtaSC_ -> Fill( localPhotons[iMatch].superCluster()->eta()- trueEta ); 00281 h1_deltaPhiSC_ -> Fill( localPhotons[iMatch].phi()- (*p)->momentum().phi() ); 00282 00283 h1_pho_E_->Fill( localPhotons[iMatch].energy() ); 00284 h1_pho_Eta_->Fill( localPhotons[iMatch].eta() ); 00285 h1_pho_Phi_->Fill( localPhotons[iMatch].phi() ); 00286 00287 h1_deltaEta_ -> Fill( localPhotons[iMatch].eta()- (*p)->momentum().eta() ); 00288 h1_deltaPhi_ -> Fill( localPhotons[iMatch].phi()- (*p)->momentum().phi() ); 00289 00290 00291 00292 00293 if ( phoIsInBarrel ) { 00294 h1_recEoverTrueEBarrel_ -> Fill (localPhotons[iMatch].energy() / (*p)->momentum().e() ); 00295 h1_recESCoverTrueEBarrel_ -> Fill (localPhotons[iMatch].superCluster()->energy() / (*p)->momentum().e() ); 00296 h1_pho_R9Barrel_->Fill( r9 ); 00297 00298 00299 00300 if ( r9 > 0.93 ) 00301 h1_e5x5_unconvBarrel_ -> Fill ( e5x5/ (*p)->momentum().e() ); 00302 else 00303 h1_ePho_convBarrel_ -> Fill ( localPhotons[iMatch].energy()/ (*p)->momentum().e() ); 00304 00305 } else { 00306 h1_recEoverTrueEEndcap_ -> Fill (localPhotons[iMatch].energy() / (*p)->momentum().e() ); 00307 h1_recESCoverTrueEEndcap_ -> Fill (localPhotons[iMatch].energy() / (*p)->momentum().e() ); 00308 h1_pho_R9Endcap_->Fill( r9 ); 00309 00310 00311 if ( r9 > 0.93 ) 00312 h1_e5x5_unconvEndcap_ -> Fill ( e5x5/ (*p)->momentum().e() ); 00313 else 00314 h1_ePho_convEndcap_ -> Fill ( localPhotons[iMatch].energy()/ (*p)->momentum().e() ); 00315 00316 } 00317 00318 00319 } 00320 00321 00322 00323 00324 } // End loop over MC particles 00325 00326 } 00327 00328 00329 }
virtual void SimplePhotonAnalyzer::beginJob | ( | edm::EventSetup const & | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
void SimplePhotonAnalyzer::beginJob | ( | edm::EventSetup const & | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 60 of file SimplePhotonAnalyzer.cc.
References dPhi(), effEta_, effPhi_, h1_deltaEta_, h1_deltaEtaSC_, h1_deltaPhi_, h1_deltaPhiSC_, h1_e5x5_unconvBarrel_, h1_e5x5_unconvEndcap_, h1_ePho_convBarrel_, h1_ePho_convEndcap_, h1_pho_E_, h1_pho_Eta_, h1_pho_Phi_, h1_pho_R9Barrel_, h1_pho_R9Endcap_, h1_recEoverTrueEBarrel_, h1_recEoverTrueEEndcap_, h1_recESCoverTrueEBarrel_, h1_recESCoverTrueEEndcap_, h1_scE_, h1_scEt_, h1_scEta_, h1_scPhi_, and sample_.
00060 { 00061 //======================================================================== 00062 00063 edm::Service<TFileService> fs; 00064 00065 float hiE=0; 00066 float loE=0; 00067 float hiEt=0; 00068 float loEt=0; 00069 float dPhi=0; 00070 float loRes=0; 00071 float hiRes=0; 00072 if ( sample_ ==1 ) { 00073 loE=0.; 00074 hiE=30.; 00075 loEt=0.; 00076 hiEt=30.; 00077 dPhi=0.2; 00078 loRes=0.; 00079 hiRes=1.2; 00080 } else if ( sample_ ==2 ) { 00081 loE=0.; 00082 hiE=200.; 00083 loEt=0.; 00084 hiEt=50.; 00085 dPhi=0.05; 00086 loRes=0.7; 00087 hiRes=1.2; 00088 } else if (sample_==4) { 00089 loE=0.; 00090 hiE=6000.; 00091 loEt=0.; 00092 hiEt=1200.; 00093 dPhi=0.05; 00094 loRes=0.7; 00095 hiRes=1.2; 00096 } 00097 00098 00099 effEta_ = fs->make<TProfile> ("effEta"," Photon reconstruction efficiency",50,-2.5,2.5); 00100 effPhi_ = fs->make<TProfile> ("effPhi"," Photon reconstruction efficiency",80, -3.14, 3.14); 00101 00102 h1_scEta_ = fs->make<TH1F>("scEta"," SC Eta ",40,-3., 3.); 00103 h1_scPhi_ = fs->make<TH1F>("scPhi"," SC Phi ",40,-3.14, 3.14); 00104 h1_deltaEtaSC_ = fs->make<TH1F>("deltaEtaSC"," SC Eta minus Generated photon Eta ",100,-0.02, 0.02); 00105 h1_deltaPhiSC_ = fs->make<TH1F>("deltaPhiSC"," SC Phi minus Generated photon Phi ",100,-dPhi, dPhi); 00106 h1_deltaEta_ = fs->make<TH1F>("deltaEta"," Reco photon Eta minus Generated photon Eta ",100,-0.2, 0.2); 00107 h1_deltaPhi_ = fs->make<TH1F>("deltaPhi","Reco photon Phi minus Generated photon Phi ",100,-dPhi, dPhi); 00108 h1_pho_Eta_ = fs->make<TH1F>("phoEta","Photon Eta ",40,-3., 3.); 00109 h1_pho_Phi_ = fs->make<TH1F>("phoPhi","Photon Phi ",40,-3.14, 3.14); 00110 00111 00112 h1_scEt_ = fs->make<TH1F>("scEt"," SC Et ",100,loEt,hiEt); 00113 h1_scE_ = fs->make<TH1F>("scE"," SC Energy ",100,loE,hiE); 00114 h1_pho_E_ = fs->make<TH1F>("phoE","Photon Energy ",100,loE,hiE); 00115 00116 h1_e5x5_unconvBarrel_ = fs->make<TH1F>("e5x5_unconvBarrelOverEtrue"," Photon rec/true energy if R9>0.93 Barrel ",100,loRes, hiRes); 00117 h1_e5x5_unconvEndcap_ = fs->make<TH1F>("e5x5_unconvEndcapOverEtrue"," Photon rec/true energy if R9>0.93 Endcap ",100,loRes, hiRes); 00118 h1_ePho_convBarrel_ = fs->make<TH1F>("ePho_convBarrelOverEtrue"," Photon rec/true energy if R9<=0.93 Barrel ",100,loRes, hiRes); 00119 h1_ePho_convEndcap_ = fs->make<TH1F>("ePho_convEndcapOverEtrue"," Photon rec/true energy if R9<=0.93 Endcap ",100,loRes, hiRes); 00120 00121 00122 // 00123 h1_recEoverTrueEBarrel_ = fs->make<TH1F>("recEoverTrueEBarrel"," Reco photon Energy over Generated photon Energy: Barrel ",100,loRes, hiRes); 00124 h1_recEoverTrueEEndcap_ = fs->make<TH1F>("recEoverTrueEEndcap"," Reco photon Energy over Generated photon Energy: Endcap ",100,loRes, hiRes); 00125 h1_recESCoverTrueEBarrel_ = fs->make<TH1F>("recESCoverTrueEBarrel"," Reco SC Energy over Generated photon Energy: Barrel ",100,loRes, hiRes); 00126 h1_recESCoverTrueEEndcap_ = fs->make<TH1F>("recESCoverTrueEEndcap"," Reco SC Energy over Generated photon Energy: Endcap ",100,loRes, hiRes); 00127 00128 // 00129 00130 h1_pho_R9Barrel_ = fs->make<TH1F>("phoR9Barrel","Photon 3x3 energy / SuperCluster energy : Barrel ",100,0.,1.2); 00131 h1_pho_R9Endcap_ = fs->make<TH1F>("phoR9Endcap","Photon 3x3 energy / SuperCluster energy : Endcap ",100,0.,1.2); 00132 00133 00134 00135 }
virtual void SimplePhotonAnalyzer::endJob | ( | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Reimplemented from edm::EDAnalyzer.
Definition at line 372 of file SimplePhotonAnalyzer.cc.
00372 { 00373 //======================================================================== 00374 00375 00376 00377 }
float SimplePhotonAnalyzer::etaTransformation | ( | float | a, | |
float | b | |||
) | [private] |
Definition at line 332 of file SimplePhotonAnalyzer.cc.
References ETA, etaBarrelEndcap, funct::log(), PI, R_ECAL, funct::tan(), TWOPI, and Z_Endcap.
Referenced by analyze().
00332 { 00333 00334 //---Definitions 00335 const float PI = 3.1415927; 00336 const float TWOPI = 2.0*PI; 00337 00338 //---Definitions for ECAL 00339 const float R_ECAL = 136.5; 00340 const float Z_Endcap = 328.0; 00341 const float etaBarrelEndcap = 1.479; 00342 00343 //---ETA correction 00344 00345 float Theta = 0.0 ; 00346 float ZEcal = R_ECAL*sinh(EtaParticle)+Zvertex; 00347 00348 if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal); 00349 if(Theta<0.0) Theta = Theta+PI ; 00350 float ETA = - log(tan(0.5*Theta)); 00351 00352 if( fabs(ETA) > etaBarrelEndcap ) 00353 { 00354 float Zend = Z_Endcap ; 00355 if(EtaParticle<0.0 ) Zend = -Zend ; 00356 float Zlen = Zend - Zvertex ; 00357 float RR = Zlen/sinh(EtaParticle); 00358 Theta = atan(RR/Zend); 00359 if(Theta<0.0) Theta = Theta+PI ; 00360 ETA = - log(tan(0.5*Theta)); 00361 } 00362 //---Return the result 00363 return ETA; 00364 //---end 00365 }
Definition at line 54 of file SimplePhotonAnalyzer.h.
Referenced by analyze(), and SimplePhotonAnalyzer().
std::string SimplePhotonAnalyzer::correctedPhotonCollection_ [private] |
Definition at line 54 of file SimplePhotonAnalyzer.h.
TProfile* SimplePhotonAnalyzer::effEta_ [private] |
TProfile* SimplePhotonAnalyzer::effPhi_ [private] |
Definition at line 55 of file SimplePhotonAnalyzer.h.
Referenced by analyze(), and SimplePhotonAnalyzer().
TH1F* SimplePhotonAnalyzer::h1_corrPho_E_ [private] |
Definition at line 74 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_corrPho_Eta_ [private] |
Definition at line 75 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_corrPho_Phi_ [private] |
Definition at line 76 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_corrPho_R9_ [private] |
Definition at line 77 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_deltaEta_ [private] |
Definition at line 84 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_deltaEta_ [private] |
TH1F* SimplePhotonAnalyzer::h1_deltaEtaSC_ [private] |
TH1F* SimplePhotonAnalyzer::h1_deltaPhi_ [private] |
Definition at line 85 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_deltaPhi_ [private] |
TH1F* SimplePhotonAnalyzer::h1_deltaPhiSC_ [private] |
TH1F* SimplePhotonAnalyzer::h1_e5x5_unconvBarrel_ [private] |
TH1F* SimplePhotonAnalyzer::h1_e5x5_unconvEndcap_ [private] |
TH1F* SimplePhotonAnalyzer::h1_ePho_convBarrel_ [private] |
TH1F* SimplePhotonAnalyzer::h1_ePho_convEndcap_ [private] |
TH1F* SimplePhotonAnalyzer::h1_pho_E_ [private] |
TH1F* SimplePhotonAnalyzer::h1_pho_Eta_ [private] |
TH1F* SimplePhotonAnalyzer::h1_pho_Phi_ [private] |
TH1F* SimplePhotonAnalyzer::h1_pho_R9Barrel_ [private] |
TH1F* SimplePhotonAnalyzer::h1_pho_R9Endcap_ [private] |
TH1F* SimplePhotonAnalyzer::h1_phoE_ [private] |
Definition at line 62 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_phoEta_ [private] |
Definition at line 63 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_phoPhi_ [private] |
Definition at line 64 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_recEoverTrueE_ [private] |
Definition at line 67 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_recEoverTrueEBarrel_ [private] |
TH1F* SimplePhotonAnalyzer::h1_recEoverTrueEEndcap_ [private] |
TH1F* SimplePhotonAnalyzer::h1_recESCoverTrueEBarrel_ [private] |
TH1F* SimplePhotonAnalyzer::h1_recESCoverTrueEEndcap_ [private] |
TH1F* SimplePhotonAnalyzer::h1_scE_ [private] |
Definition at line 66 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_scE_ [private] |
TH1F* SimplePhotonAnalyzer::h1_scEt_ [private] |
TH1F* SimplePhotonAnalyzer::h1_scEta_ [private] |
Definition at line 68 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_scEta_ [private] |
TH1F* SimplePhotonAnalyzer::h1_scPhi_ [private] |
Definition at line 69 of file SimplePhotonAnalyzer.h.
TH1F* SimplePhotonAnalyzer::h1_scPhi_ [private] |
std::string SimplePhotonAnalyzer::mcCollection_ [private] |
Definition at line 47 of file SimplePhotonAnalyzer.h.
std::string SimplePhotonAnalyzer::mcProducer_ [private] |
Definition at line 46 of file SimplePhotonAnalyzer.h.
Referenced by analyze(), and SimplePhotonAnalyzer().
std::string SimplePhotonAnalyzer::outputFile_ [private] |
Definition at line 43 of file SimplePhotonAnalyzer.h.
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 48 of file SimplePhotonAnalyzer.h.
Referenced by analyze(), and SimplePhotonAnalyzer().
std::string SimplePhotonAnalyzer::photonCorrCollectionProducer_ [private] |
Definition at line 49 of file SimplePhotonAnalyzer.h.
TFile* SimplePhotonAnalyzer::rootFile_ [private] |
Definition at line 44 of file SimplePhotonAnalyzer.h.
float SimplePhotonAnalyzer::sample_ [private] |
Definition at line 60 of file SimplePhotonAnalyzer.h.
Referenced by beginJob(), and SimplePhotonAnalyzer().
std::string SimplePhotonAnalyzer::uncorrectedPhotonCollection_ [private] |
Definition at line 50 of file SimplePhotonAnalyzer.h.
std::string SimplePhotonAnalyzer::vertexProducer_ [private] |
Definition at line 51 of file SimplePhotonAnalyzer.h.
Referenced by analyze(), and SimplePhotonAnalyzer().