00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "DQM/Physics/src/QcdPhotonsDQM.h"
00010
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/Framework/interface/EventSetup.h"
00013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00014
00015 #include "DQMServices/Core/interface/DQMStore.h"
00016 #include "DQMServices/Core/interface/MonitorElement.h"
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 #include "FWCore/Common/interface/TriggerNames.h"
00019
00020 #include "DataFormats/Common/interface/Handle.h"
00021
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023
00024
00025 #include "DataFormats/EgammaCandidates/interface/Photon.h"
00026 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
00027 #include "DataFormats/JetReco/interface/Jet.h"
00028
00029
00030 #include "DataFormats/VertexReco/interface/Vertex.h"
00031 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00032
00033
00034 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h"
00035 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00036 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
00037 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgoRcd.h"
00038 #include "FWCore/ServiceRegistry/interface/Service.h"
00039
00040 #include "FWCore/Framework/interface/ESHandle.h"
00041 #include "CondFormats/EcalObjects/interface/EcalCondObjectContainer.h"
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053 #include "DataFormats/Math/interface/deltaR.h"
00054 #include "DataFormats/Math/interface/deltaPhi.h"
00055
00056 #include <vector>
00057
00058 #include <string>
00059 #include <cmath>
00060 using namespace std;
00061 using namespace edm;
00062 using namespace reco;
00063
00064
00065 QcdPhotonsDQM::QcdPhotonsDQM(const ParameterSet& parameters) {
00066
00067 theTriggerPathToPass_ = parameters.getParameter<string>("triggerPathToPass");
00068 thePlotTheseTriggersToo_ = parameters.getParameter<vector<string> >("plotTheseTriggersToo");
00069 trigTag_ = parameters.getUntrackedParameter<edm::InputTag>("trigTag");
00070 thePhotonCollectionLabel_ = parameters.getParameter<InputTag>("photonCollection");
00071 theJetCollectionLabel_ = parameters.getParameter<InputTag>("jetCollection");
00072 theVertexCollectionLabel_ = parameters.getParameter<InputTag>("vertexCollection");
00073 theMinJetPt_ = parameters.getParameter<double>("minJetPt");
00074 theMinPhotonEt_ = parameters.getParameter<double>("minPhotonEt");
00075 theRequirePhotonFound_ = parameters.getParameter<bool>("requirePhotonFound");
00076 thePlotPhotonMaxEt_ = parameters.getParameter<double>("plotPhotonMaxEt");
00077 thePlotPhotonMaxEta_ = parameters.getParameter<double>("plotPhotonMaxEta");
00078 thePlotJetMaxEta_ = parameters.getParameter<double>("plotJetMaxEta");
00079
00080 isValidHltConfig_ = false;
00081
00082
00083 h_deltaEt_photon_jet = 0;
00084 h_deltaPhi_jet_jet2 = 0;
00085 h_deltaPhi_photon_jet = 0;
00086 h_deltaPhi_photon_jet2 = 0;
00087 h_deltaR_jet_jet2 = 0;
00088 h_deltaR_photon_jet2 = 0;
00089 h_jet2_eta = 0;
00090 h_jet2_pt = 0;
00091 h_jet2_ptOverPhotonEt = 0;
00092 h_jet_count = 0;
00093 h_jet_eta = 0;
00094 h_jet_pt = 0;
00095 h_photon_count_bar = 0;
00096 h_photon_count_end = 0;
00097 h_photon_et = 0;
00098 h_photon_et_beforeCuts = 0;
00099 h_photon_et_jetco = 0;
00100 h_photon_et_jetcs = 0;
00101 h_photon_et_jetfo = 0;
00102 h_photon_et_jetfs = 0;
00103 h_photon_et_ratio_co_cs = 0;
00104 h_photon_et_ratio_co_fo = 0;
00105 h_photon_et_ratio_co_fs = 0;
00106 h_photon_et_ratio_cs_fo = 0;
00107 h_photon_et_ratio_cs_fs = 0;
00108 h_photon_et_ratio_fo_fs = 0;
00109 h_photon_eta = 0;
00110 h_triggers_passed = 0;
00111
00112 theDbe = Service<DQMStore>().operator->();
00113
00114 }
00115
00116 QcdPhotonsDQM::~QcdPhotonsDQM() {
00117 }
00118
00119
00120 void QcdPhotonsDQM::beginJob() {
00121
00122 logTraceName = "QcdPhotonAnalyzer";
00123
00124 LogTrace(logTraceName)<<"Parameters initialization";
00125
00126 theDbe->setCurrentFolder("Physics/QcdPhotons");
00127
00128 std::stringstream aStringStream;
00129 std::string aString;
00130 aStringStream << theMinJetPt_;
00131 aString = aStringStream.str();
00132
00133
00134 int numOfTriggersToMonitor = thePlotTheseTriggersToo_.size();
00135 h_triggers_passed = theDbe->book1D("triggers_passed", "Events passing these trigger paths", numOfTriggersToMonitor, 0, numOfTriggersToMonitor);
00136 for (int i=0; i<numOfTriggersToMonitor; i++) {
00137 h_triggers_passed->setBinLabel(i+1,thePlotTheseTriggersToo_[i]);
00138 }
00139
00140
00141 h_photon_et_beforeCuts = theDbe->book1D("photon_et_beforeCuts", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
00142 h_photon_et = theDbe->book1D("photon_et", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
00143 h_photon_eta = theDbe->book1D("photon_eta", "#gamma with highest E_{T};#eta(#gamma)", 40, -thePlotPhotonMaxEta_, thePlotPhotonMaxEta_);
00144 h_photon_count_bar = theDbe->book1D("photon_count_bar","Number of #gamma's passing selection (Barrel);Number of #gamma's", 8, -0.5, 7.5);
00145 h_photon_count_end = theDbe->book1D("photon_count_end","Number of #gamma's passing selection (Endcap);Number of #gamma's", 8, -0.5, 7.5);
00146
00147 h_jet_pt = theDbe->book1D("jet_pt", "Jet with highest p_{T} (from "+theJetCollectionLabel_.label()+");p_{T}(1^{st} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
00148 h_jet_eta = theDbe->book1D("jet_eta", "Jet with highest p_{T} (from "+theJetCollectionLabel_.label()+");#eta(1^{st} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
00149 h_deltaPhi_photon_jet = theDbe->book1D("deltaPhi_photon_jet", "#Delta#phi between Highest E_{T} #gamma and jet;#Delta#phi(#gamma,1^{st} jet)", 20, 0, 3.1415926);
00150 h_deltaPhi_jet_jet2 = theDbe->book1D("deltaPhi_jet_jet2", "#Delta#phi between Highest E_{T} jet and 2^{nd} jet;#Delta#phi(1^{st} jet,2^{nd} jet)", 20, 0, 3.1415926);
00151 h_deltaEt_photon_jet = theDbe->book1D("deltaEt_photon_jet", "(E_{T}(#gamma)-p_{T}(jet))/E_{T}(#gamma) when #Delta#phi(#gamma,1^{st} jet) > 2.8;#DeltaE_{T}(#gamma,1^{st} jet)/E_{T}(#gamma)", 20, -1.0, 1.0);
00152 h_jet_count = theDbe->book1D("jet_count", "Number of "+theJetCollectionLabel_.label()+" (p_{T} > "+aString+" GeV);Number of Jets", 8, -0.5, 7.5);
00153 h_jet2_pt = theDbe->book1D("jet2_pt", "Jet with 2^{nd} highest p_{T} (from "+theJetCollectionLabel_.label()+");p_{T}(2^{nd} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
00154 h_jet2_eta = theDbe->book1D("jet2_eta", "Jet with 2^{nd} highest p_{T} (from "+theJetCollectionLabel_.label()+");#eta(2^{nd} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
00155 h_jet2_ptOverPhotonEt = theDbe->book1D("jet2_ptOverPhotonEt", "p_{T}(2^{nd} highest jet) / E_{T}(#gamma);p_{T}(2^{nd} Jet)/E_{T}(#gamma)", 20, 0.0, 4.0);
00156 h_deltaPhi_photon_jet2 = theDbe->book1D("deltaPhi_photon_jet2","#Delta#phi between Highest E_{T} #gamma and 2^{nd} highest jet;#Delta#phi(#gamma,2^{nd} jet)", 20, 0, 3.1415926);
00157 h_deltaR_jet_jet2 = theDbe->book1D("deltaR_jet_jet2", "#DeltaR between Highest Jet and 2^{nd} Highest;#DeltaR(1^{st} jet,2^{nd} jet)", 30, 0, 6.0);
00158 h_deltaR_photon_jet2 = theDbe->book1D("deltaR_photon_jet2", "#DeltaR between Highest E_{T} #gamma and 2^{nd} jet;#DeltaR(#gamma, 2^{nd} jet)", 30, 0, 6.0);
00159
00160
00161 Float_t bins_et[] = {15,20,30,50,80};
00162 int num_bins_et = 4;
00163 h_photon_et_jetcs = theDbe->book1D("photon_et_jetcs", "#gamma with highest E_{T} (#eta(jet)<1.45, #eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
00164 h_photon_et_jetco = theDbe->book1D("photon_et_jetco", "#gamma with highest E_{T} (#eta(jet)<1.45, #eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
00165 h_photon_et_jetfs = theDbe->book1D("photon_et_jetfs", "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, #eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
00166 h_photon_et_jetfo = theDbe->book1D("photon_et_jetfo", "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, #eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
00167 h_photon_et_jetcs->getTH1F()->Sumw2();
00168 h_photon_et_jetco->getTH1F()->Sumw2();
00169 h_photon_et_jetfs->getTH1F()->Sumw2();
00170 h_photon_et_jetfo->getTH1F()->Sumw2();
00171
00172 h_photon_et_ratio_co_cs = theDbe->book1D("photon_et_ratio_00_co_cs", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
00173 h_photon_et_ratio_fo_fs = theDbe->book1D("photon_et_ratio_01_fo_fs", "D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)<0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
00174 h_photon_et_ratio_cs_fs = theDbe->book1D("photon_et_ratio_02_cs_fs", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)>0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
00175 h_photon_et_ratio_co_fs = theDbe->book1D("photon_et_ratio_03_co_fs", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
00176 h_photon_et_ratio_cs_fo = theDbe->book1D("photon_et_ratio_04_cs_fo", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)>0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)<0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
00177 h_photon_et_ratio_co_fo = theDbe->book1D("photon_et_ratio_05_co_fo", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)<0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
00178 h_photon_et_ratio_co_cs->getTH1F()->Sumw2();
00179 h_photon_et_ratio_fo_fs->getTH1F()->Sumw2();
00180 h_photon_et_ratio_cs_fs->getTH1F()->Sumw2();
00181 h_photon_et_ratio_co_fs->getTH1F()->Sumw2();
00182 h_photon_et_ratio_cs_fo->getTH1F()->Sumw2();
00183 h_photon_et_ratio_co_fo->getTH1F()->Sumw2();
00184 }
00185
00186
00190 void QcdPhotonsDQM::beginRun( const edm::Run &iRun, const edm::EventSetup &iSet ) {
00191
00192
00193 bool isConfigChanged = false;
00194
00195
00196 isValidHltConfig_ = hltConfigProvider_.init( iRun, iSet, "HLT", isConfigChanged );
00197
00198 num_events_in_run = 0;
00199 }
00200
00201
00202 void QcdPhotonsDQM::analyze(const Event& iEvent, const EventSetup& iSetup) {
00203 num_events_in_run++;
00204
00206
00207
00208 LogTrace(logTraceName)<<"Analysis of event # ";
00209
00211
00212 Handle<TriggerResults> HLTresults;
00213 iEvent.getByLabel(trigTag_, HLTresults);
00214 if (!HLTresults.isValid()) {
00215
00216 return;
00217 }
00218 const edm::TriggerNames & trigNames = iEvent.triggerNames(*HLTresults);
00219
00220 bool passed_HLT=false;
00221
00222
00223
00224 for (unsigned int i=0; i<thePlotTheseTriggersToo_.size(); i++) {
00225 passed_HLT = false;
00226 for (unsigned int ti=0; (ti<trigNames.size()) && !passed_HLT; ++ti) {
00227 size_t pos = trigNames.triggerName(ti).find(thePlotTheseTriggersToo_[i]);
00228 if (pos==0) passed_HLT = HLTresults->accept(ti);
00229 }
00230 if (passed_HLT) h_triggers_passed->Fill(i);
00231 }
00232
00233
00234
00235 Handle<PhotonCollection> photonCollection;
00236 iEvent.getByLabel(thePhotonCollectionLabel_, photonCollection);
00237
00238
00239 if (!photonCollection.isValid()) return;
00240
00241
00242
00243 passed_HLT = false;
00244 {
00245
00246 for (unsigned int ti=0; ti<trigNames.size(); ++ti) {
00247 size_t pos = trigNames.triggerName(ti).find(theTriggerPathToPass_);
00248 if (pos==0) {
00249 passed_HLT = HLTresults->accept(ti);
00250
00251 break;
00252 }
00253 }
00254
00255
00256 for (PhotonCollection::const_iterator recoPhoton = photonCollection->begin(); recoPhoton!=photonCollection->end(); recoPhoton++){
00257
00258 if ( recoPhoton->et() < theMinPhotonEt_ ) break;
00259
00260 h_photon_et_beforeCuts->Fill(recoPhoton->et());
00261 break;
00262 }
00263
00264 if (!passed_HLT) {
00265 return;
00266 }
00267 }
00268
00270
00271
00272
00274
00275
00276 Handle<VertexCollection> vertexHandle;
00277 iEvent.getByLabel(theVertexCollectionLabel_, vertexHandle);
00278 VertexCollection vertexCollection = *(vertexHandle.product());
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289 int nvvertex = 0;
00290 for (unsigned int i=0; i<vertexCollection.size(); ++i) {
00291 if (vertexCollection[i].isValid()) nvvertex++;
00292 }
00293 if (nvvertex==0) return;
00294
00296
00297
00298
00299
00301
00302
00303
00305
00306
00307
00308 Handle<EcalRecHitCollection> EBReducedRecHits;
00309 iEvent.getByLabel("reducedEcalRecHitsEB", EBReducedRecHits);
00310 Handle<EcalRecHitCollection> EEReducedRecHits;
00311 iEvent.getByLabel("reducedEcalRecHitsEE", EEReducedRecHits);
00312 EcalClusterLazyTools lazyTool(iEvent, iSetup, InputTag("reducedEcalRecHitsEB"), InputTag("reducedEcalRecHitsEE") );
00313
00314
00315
00316 float photon_et = -9.0;
00317 float photon_eta = -9.0;
00318 float photon_phi = -9.0;
00319 bool photon_passPhotonID = false;
00320 bool found_lead_pho = false;
00321 int photon_count_bar = 0;
00322 int photon_count_end = 0;
00323
00324 for (PhotonCollection::const_iterator recoPhoton = photonCollection->begin(); recoPhoton!=photonCollection->end(); recoPhoton++){
00325
00326
00327 if ( recoPhoton->et() < theMinPhotonEt_ ) break;
00328
00329
00330 const reco::CaloClusterPtr seed = recoPhoton->superCluster()->seed();
00331 DetId id = lazyTool.getMaximum(*seed).first;
00332
00333 int flags=-1, severity = -1;
00334 const EcalRecHitCollection & rechits = ( recoPhoton->isEB() ? *EBReducedRecHits : *EEReducedRecHits);
00335 EcalRecHitCollection::const_iterator it = rechits.find( id );
00336 if( it != rechits.end() ) {
00337
00338
00339
00340 flags = it->recoFlag();
00341
00342 edm::ESHandle<EcalSeverityLevelAlgo> sevlv;
00343 iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
00344 severity = sevlv->severityLevel( id, rechits);
00345 }
00346 bool isNotSpike = ((recoPhoton->isEB() && (severity!=3 && severity!=4 ) && (flags != 2) ) || recoPhoton->isEE());
00347 if (!isNotSpike) continue;
00348
00349
00350 bool pho_current_passPhotonID = false;
00351 bool pho_current_isEB = recoPhoton->isEB();
00352 bool pho_current_isEE = recoPhoton->isEE();
00353
00354 if ( pho_current_isEB && (recoPhoton->sigmaIetaIeta() < 0.01 || recoPhoton->hadronicOverEm() < 0.05) ) {
00355
00356 pho_current_passPhotonID = true;
00357 photon_count_bar++;
00358 } else if ( pho_current_isEE && (recoPhoton->hadronicOverEm() < 0.05) ) {
00359
00360 pho_current_passPhotonID = true;
00361 photon_count_end++;
00362 }
00363
00364 if (!found_lead_pho) {
00365 found_lead_pho = true;
00366 photon_passPhotonID = pho_current_passPhotonID;
00367 photon_et = recoPhoton->et();
00368 photon_eta = recoPhoton->eta();
00369 photon_phi = recoPhoton->phi();
00370 }
00371 }
00372
00373
00374
00375
00376 if ( theRequirePhotonFound_ && (!photon_passPhotonID || photon_et<theMinPhotonEt_) ) return;
00377
00378
00380
00381 Handle<View<Jet> > jetCollection;
00382 iEvent.getByLabel (theJetCollectionLabel_,jetCollection);
00383 if (!jetCollection.isValid()) return;
00384
00385 float jet_pt = -8.0;
00386 float jet_eta = -8.0;
00387 float jet_phi = -8.0;
00388 int jet_count = 0;
00389 float jet2_pt = -9.0;
00390 float jet2_eta = -9.0;
00391 float jet2_phi = -9.0;
00392
00393 for (unsigned int i_jet = 0; i_jet < jetCollection->size(); i_jet++) {
00394 const Jet* jet = & jetCollection->at(i_jet);
00395
00396 float jet_current_pt = jet->pt();
00397
00398
00399 if ( deltaR(jet->eta(), jet->phi(), photon_eta, photon_phi) < 0.5 ) continue;
00400
00401 if (jet_current_pt < theMinJetPt_) break;
00402
00403 jet_count++;
00404 if (jet_current_pt > jet_pt) {
00405 jet2_pt = jet_pt;
00406 jet2_eta = jet_eta;
00407 jet2_phi = jet_phi;
00408 jet_pt = jet_current_pt;
00409 jet_eta = jet->eta();
00410 jet_phi = jet->phi();
00411 } else if (jet_current_pt > jet2_pt) {
00412 jet2_pt = jet_current_pt;
00413 jet2_eta = jet->eta();
00414 jet2_phi = jet->phi();
00415 }
00416 }
00418
00419
00421
00422
00423
00424 if ( jet_pt > 0.0 ) {
00425
00426
00427 h_photon_et ->Fill( photon_et );
00428 h_photon_eta ->Fill( photon_eta );
00429 h_photon_count_bar->Fill( photon_count_bar );
00430 h_photon_count_end->Fill( photon_count_end );
00431
00432
00433 if ( fabs(photon_eta)<1.45 && photon_passPhotonID ) {
00434 if (fabs(jet_eta)<1.45){
00435 if (photon_eta*jet_eta>0) {
00436 h_photon_et_jetcs->Fill(photon_et);
00437 } else {
00438 h_photon_et_jetco->Fill(photon_et);
00439 }
00440 } else if (jet_eta>1.55 && jet_eta<2.5) {
00441 if (photon_eta*jet_eta>0) {
00442 h_photon_et_jetfs->Fill(photon_et);
00443 } else {
00444 h_photon_et_jetfo->Fill(photon_et);
00445 }
00446 }
00447 }
00448
00449
00450 h_jet_pt ->Fill( jet_pt );
00451 h_jet_eta ->Fill( jet_eta );
00452 h_jet_count ->Fill( jet_count );
00453 h_deltaPhi_photon_jet ->Fill( abs(deltaPhi(photon_phi, jet_phi)) );
00454 if ( abs(deltaPhi(photon_phi,jet_phi))>2.8 ) h_deltaEt_photon_jet->Fill( (photon_et-jet_pt)/photon_et );
00455
00456
00457 if ( jet2_pt > 0.0 ) {
00458 h_jet2_pt ->Fill( jet2_pt );
00459 h_jet2_eta ->Fill( jet2_eta );
00460 h_jet2_ptOverPhotonEt ->Fill( jet2_pt/photon_et );
00461 h_deltaPhi_photon_jet2->Fill( abs(deltaPhi(photon_phi, jet2_phi)) );
00462 h_deltaPhi_jet_jet2 ->Fill( abs(deltaPhi( jet_phi, jet2_phi)) );
00463 h_deltaR_jet_jet2 ->Fill( deltaR( jet_eta, jet_phi, jet2_eta, jet2_phi) );
00464 h_deltaR_photon_jet2 ->Fill( deltaR(photon_eta, photon_phi, jet2_eta, jet2_phi) );
00465 }
00466 }
00467
00469 }
00470
00471
00472 void QcdPhotonsDQM::endJob(void) {}
00473
00474 void QcdPhotonsDQM::endRun(const edm::Run& run, const edm::EventSetup& es) {
00475 if (num_events_in_run>0) {
00476 h_triggers_passed->getTH1F()->Scale(1.0/num_events_in_run);
00477 }
00478 h_photon_et_ratio_co_cs->getTH1F()->Divide( h_photon_et_jetco->getTH1F(), h_photon_et_jetcs->getTH1F() );
00479 h_photon_et_ratio_fo_fs->getTH1F()->Divide( h_photon_et_jetfo->getTH1F(), h_photon_et_jetfs->getTH1F() );
00480 h_photon_et_ratio_cs_fs->getTH1F()->Divide( h_photon_et_jetcs->getTH1F(), h_photon_et_jetfs->getTH1F() );
00481 h_photon_et_ratio_co_fs->getTH1F()->Divide( h_photon_et_jetco->getTH1F(), h_photon_et_jetfs->getTH1F() );
00482 h_photon_et_ratio_cs_fo->getTH1F()->Divide( h_photon_et_jetcs->getTH1F(), h_photon_et_jetfo->getTH1F() );
00483 h_photon_et_ratio_co_fo->getTH1F()->Divide( h_photon_et_jetco->getTH1F(), h_photon_et_jetfo->getTH1F() );
00484 }
00485