CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DQM/Physics/src/QcdPhotonsDQM.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2011/01/21 10:21:39 $
00005  *  $Revision: 1.27 $
00006  *  \author Michael B. Anderson, University of Wisconsin Madison
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 
00019 #include "DataFormats/Common/interface/Handle.h"
00020 
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 
00023 // Physics Objects
00024 #include "DataFormats/EgammaCandidates/interface/Photon.h"
00025 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
00026 #include "DataFormats/JetReco/interface/CaloJet.h"
00027 
00028 // Vertex
00029 #include "DataFormats/VertexReco/interface/Vertex.h"
00030 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00031 
00032 // For removing ECAL Spikes
00033 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h"
00034 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00035 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
00036 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgoRcd.h"
00037 #include "FWCore/ServiceRegistry/interface/Service.h"
00038 
00039 #include "FWCore/Framework/interface/ESHandle.h"
00040 #include "CondFormats/EcalObjects/interface/EcalCondObjectContainer.h"
00041 
00042 //geometry
00043 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00044 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00045 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00046 #include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
00047 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00048 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00049 #include "Geometry/EcalAlgo/interface/EcalPreshowerGeometry.h"
00050 
00051 // Math stuff
00052 #include "DataFormats/Math/interface/deltaR.h"
00053 #include "DataFormats/Math/interface/deltaPhi.h"
00054 
00055 #include <vector>
00056 
00057 #include <string>
00058 #include <cmath>
00059 using namespace std;
00060 using namespace edm;
00061 using namespace reco;
00062 
00063 
00064 
00065 QcdPhotonsDQM::QcdPhotonsDQM(const ParameterSet& parameters) {
00066   // Get parameters from configuration file
00067   theTriggerPathToPass        = parameters.getParameter<string>("triggerPathToPass");
00068   thePlotTheseTriggersToo     = parameters.getParameter<vector<string> >("plotTheseTriggersToo");
00069   theHltMenu                  = parameters.getParameter<string>("hltMenu");
00070   theTriggerResultsCollection = parameters.getParameter<string>("triggerResultsCollection");
00071   thePhotonCollectionLabel    = parameters.getParameter<InputTag>("photonCollection");
00072   theCaloJetCollectionLabel   = parameters.getParameter<InputTag>("caloJetCollection");
00073   theVertexCollectionLabel    = parameters.getParameter<InputTag>("vertexCollection");
00074   theMinCaloJetPt             = parameters.getParameter<double>("minCaloJetPt");
00075   theMinPhotonEt              = parameters.getParameter<double>("minPhotonEt");
00076   theRequirePhotonFound       = parameters.getParameter<bool>("requirePhotonFound");
00077   thePlotPhotonMaxEt          = parameters.getParameter<double>("plotPhotonMaxEt");
00078   thePlotPhotonMaxEta         = parameters.getParameter<double>("plotPhotonMaxEta");
00079   thePlotJetMaxEta            = parameters.getParameter<double>("plotJetMaxEta");
00080   // just to initialize
00081   isValidHltConfig_ = false;
00082 
00083 }
00084 
00085 QcdPhotonsDQM::~QcdPhotonsDQM() { 
00086 }
00087 
00088 
00089 void QcdPhotonsDQM::beginJob() {
00090  
00091   logTraceName = "QcdPhotonAnalyzer";
00092 
00093   LogTrace(logTraceName)<<"Parameters initialization";
00094   theDbe = Service<DQMStore>().operator->();
00095  
00096   theDbe->setCurrentFolder("Physics/QcdPhotons");  // Use folder with name of PAG
00097 
00098   std::stringstream aStringStream;
00099   std::string aString;
00100   aStringStream << theMinCaloJetPt;
00101   aString = aStringStream.str();
00102 
00103   // Monitor of triggers passed
00104   int numOfTriggersToMonitor = thePlotTheseTriggersToo.size();
00105   h_triggers_passed = theDbe->book1D("h_triggers_passed", "Events passing these trigger paths", numOfTriggersToMonitor, 0, numOfTriggersToMonitor);
00106   for (int i=0; i<numOfTriggersToMonitor; i++) {
00107     h_triggers_passed->setBinLabel(i+1,thePlotTheseTriggersToo[i]);
00108   }
00109 
00110   // Keep the number of plots and number of bins to a minimum!
00111   h_photon_et           = theDbe->book1D("h_photon_et",       "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt);
00112   h_photon_eta          = theDbe->book1D("h_photon_eta",      "#gamma with highest E_{T};#eta(#gamma)", 40, -thePlotPhotonMaxEta, thePlotPhotonMaxEta);
00113   h_photon_count_bar    = theDbe->book1D("h_photon_count_bar","Number of #gamma's passing selection (Barrel);Number of #gamma's", 8, -0.5, 7.5);
00114   h_photon_count_end    = theDbe->book1D("h_photon_count_end","Number of #gamma's passing selection (Endcap);Number of #gamma's", 8, -0.5, 7.5);
00115 
00116   h_jet_pt             = theDbe->book1D("h_jet_pt",   "Jet with highest p_{T} (from "+theCaloJetCollectionLabel.label()+");p_{T}(1^{st} jet) (GeV)",    20, 0., thePlotPhotonMaxEt);
00117   h_jet_eta             = theDbe->book1D("h_jet_eta", "Jet with highest p_{T} (from "+theCaloJetCollectionLabel.label()+");#eta(1^{st} jet)", 20, -thePlotJetMaxEta, thePlotJetMaxEta);
00118   h_deltaPhi_photon_jet = theDbe->book1D("h_deltaPhi_photon_jet", "#Delta#phi between Highest E_{T} #gamma and jet;#Delta#phi(#gamma,1^{st} jet)", 20, 0, 3.1415926);
00119   h_deltaPhi_jet_jet2   = theDbe->book1D("h_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);
00120   h_deltaEt_photon_jet  = theDbe->book1D("h_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);
00121   h_jet_count           = theDbe->book1D("h_jet_count", "Number of "+theCaloJetCollectionLabel.label()+" (p_{T} > "+aString+" GeV);Number of Jets", 8, -0.5, 7.5);
00122   h_jet2_pt             = theDbe->book1D("h_jet2_pt",   "Jet with 2^{nd} highest p_{T} (from "+theCaloJetCollectionLabel.label()+");p_{T}(2^{nd} jet) (GeV)",    20, 0., thePlotPhotonMaxEt);
00123   h_jet2_eta            = theDbe->book1D("h_jet2_eta", "Jet with 2^{nd} highest p_{T} (from "+theCaloJetCollectionLabel.label()+");#eta(2^{nd} jet)", 20, -thePlotJetMaxEta, thePlotJetMaxEta);
00124   h_jet2_ptOverPhotonEt = theDbe->book1D("h_jet2_ptOverPhotonEt", "p_{T}(2^{nd} highest jet) / E_{T}(#gamma);p_{T}(2^{nd} Jet)/E_{T}(#gamma)", 20, 0.0, 4.0);
00125   h_deltaPhi_photon_jet2 = theDbe->book1D("h_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);
00126   h_deltaR_jet_jet2      = theDbe->book1D("h_deltaR_jet_jet2", "#DeltaR between Highest Jet and 2^{nd} Highest;#DeltaR(1^{st} jet,2^{nd} jet)", 30, 0, 6.0);
00127   h_deltaR_photon_jet2   = theDbe->book1D("h_deltaR_photon_jet2", "#DeltaR between Highest E_{T} #gamma and 2^{nd} jet;#DeltaR(#gamma, 2^{nd} jet)", 30, 0, 6.0);
00128 
00129   // Photon Et for different jet configurations
00130   Float_t bins_et[] = {15,20,30,50,80};
00131   int num_bins_et = 4;
00132   h_photon_et_jetcs = theDbe->book1D("h_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);
00133   h_photon_et_jetco = theDbe->book1D("h_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);
00134   h_photon_et_jetfs = theDbe->book1D("h_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);
00135   h_photon_et_jetfo = theDbe->book1D("h_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);
00136   h_photon_et_jetcs->getTH1F()->Sumw2();
00137   h_photon_et_jetco->getTH1F()->Sumw2();
00138   h_photon_et_jetfs->getTH1F()->Sumw2();
00139   h_photon_et_jetfo->getTH1F()->Sumw2();
00140   // Ratio of the above Photon Et distributions
00141   h_photon_et_ratio_co_cs = theDbe->book1D("h_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);
00142   h_photon_et_ratio_fo_fs = theDbe->book1D("h_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);
00143   h_photon_et_ratio_cs_fs = theDbe->book1D("h_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);
00144   h_photon_et_ratio_co_fs = theDbe->book1D("h_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);
00145   h_photon_et_ratio_cs_fo = theDbe->book1D("h_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);
00146   h_photon_et_ratio_co_fo = theDbe->book1D("h_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);
00147   h_photon_et_ratio_co_cs->getTH1F()->Sumw2();
00148   h_photon_et_ratio_fo_fs->getTH1F()->Sumw2();
00149   h_photon_et_ratio_cs_fs->getTH1F()->Sumw2();
00150   h_photon_et_ratio_co_fs->getTH1F()->Sumw2();
00151   h_photon_et_ratio_cs_fo->getTH1F()->Sumw2();
00152   h_photon_et_ratio_co_fo->getTH1F()->Sumw2();
00153 }
00154 
00155 
00159 void QcdPhotonsDQM::beginRun( const edm::Run &r, const edm::EventSetup &iSetup ) {
00160 
00161   // passed as parameter to HLTConfigProvider::init(), not yet used
00162   bool isConfigChanged = false;
00163 
00164   // isValidHltConfig_ used to short-circuit analyze() in case of problems
00165   isValidHltConfig_ = hltConfigProvider_.init( r, iSetup, theHltMenu, isConfigChanged );
00166 
00167   num_events_in_run = 0;
00168 }
00169 
00170 
00171 void QcdPhotonsDQM::analyze(const Event& iEvent, const EventSetup& iSetup) {
00172   num_events_in_run++;
00173 
00174   // short-circuit if hlt problems
00175   if( ! isValidHltConfig_ ) return;
00176   
00177   LogTrace(logTraceName)<<"Analysis of event # ";
00178 
00180   // Did event pass HLT paths?
00181   Handle<TriggerResults> HLTresults;
00182   iEvent.getByLabel(InputTag(theTriggerResultsCollection, "", theHltMenu), HLTresults); 
00183   if (!HLTresults.isValid()) return;
00184 
00185   unsigned int triggerIndex; // index of trigger path
00186   bool passed_HLT;
00187 
00188   // See if event passed trigger paths
00189   //  increment that bin in the trigger plot
00190   for (unsigned int i=0; i<thePlotTheseTriggersToo.size(); i++) {
00191     passed_HLT = false;
00192     triggerIndex = hltConfigProvider_.triggerIndex(thePlotTheseTriggersToo[i]);
00193     if (triggerIndex < HLTresults->size()) passed_HLT = HLTresults->accept(triggerIndex);
00194     if (passed_HLT) h_triggers_passed->Fill(i);
00195   }
00196 
00197   // Quit if the event did not pass the HLT path we care about
00198   passed_HLT = false;
00199   triggerIndex = hltConfigProvider_.triggerIndex(theTriggerPathToPass); // index of trigger path
00200   if (triggerIndex < HLTresults->size()) passed_HLT = HLTresults->accept(triggerIndex);
00201   if (!passed_HLT) return;
00203 
00204 
00206   // Does event have valid vertex?
00207   // Get the primary event vertex
00208   Handle<VertexCollection> vertexHandle;
00209   iEvent.getByLabel(theVertexCollectionLabel, vertexHandle);
00210   VertexCollection vertexCollection = *(vertexHandle.product());
00211   double vtx_ndof = -1.0;
00212   double vtx_z    = 0.0;
00213   bool   vtx_isFake = true;
00214   if (vertexCollection.size()>0) {
00215     vtx_ndof = vertexCollection.begin()->ndof();
00216     vtx_z    = vertexCollection.begin()->z();
00217     vtx_isFake = false;
00218   }
00219   if (vtx_isFake || fabs(vtx_z)>15 || vtx_ndof<4) return;
00221 
00222 
00224   // Did the event pass certain L1 Technical Trigger bits?
00225   // It's probably beam halo
00226   //  TODO: ADD code
00228 
00229   // grab photons
00230   Handle<PhotonCollection> photonCollection;
00231   iEvent.getByLabel(thePhotonCollectionLabel, photonCollection);
00232 
00233   // If photon collection is empty, exit
00234   if (!photonCollection.isValid()) return;
00235 
00236   // For finding spikes
00237   Handle<EcalRecHitCollection> EBReducedRecHits;
00238   iEvent.getByLabel("reducedEcalRecHitsEB", EBReducedRecHits);
00239   Handle<EcalRecHitCollection> EEReducedRecHits;
00240   iEvent.getByLabel("reducedEcalRecHitsEE", EEReducedRecHits); 
00241   EcalClusterLazyTools lazyTool(iEvent, iSetup, InputTag("reducedEcalRecHitsEB"), InputTag("reducedEcalRecHitsEE") );
00242   
00243 
00244   // Find the highest et "decent" photon
00245   float photon_et  = -9.0;
00246   float photon_eta = -9.0;
00247   float photon_phi = -9.0;
00248   bool  photon_passPhotonID = false;
00249   bool  found_lead_pho = false;
00250   int   photon_count_bar = 0;
00251   int   photon_count_end = 0;
00252   // Assumption: reco photons are ordered by Et
00253   for (PhotonCollection::const_iterator recoPhoton = photonCollection->begin(); recoPhoton!=photonCollection->end(); recoPhoton++){
00254 
00255     // stop looping over photons once we get to too low Et
00256     if ( recoPhoton->et() < theMinPhotonEt ) break;
00257 
00258     //  Ignore ECAL Spikes
00259     const reco::CaloClusterPtr  seed = recoPhoton->superCluster()->seed();
00260     DetId id = lazyTool.getMaximum(*seed).first; // Cluster shape variables
00261     float time  = -999., outOfTimeChi2 = -999., chi2 = -999.;
00262     int   flags=-1, severity = -1; 
00263     const EcalRecHitCollection & rechits = ( recoPhoton->isEB() ? *EBReducedRecHits : *EEReducedRecHits); 
00264     EcalRecHitCollection::const_iterator it = rechits.find( id );
00265     if( it != rechits.end() ) {
00266       time = it->time(); 
00267       outOfTimeChi2 = it->outOfTimeChi2();
00268       chi2 = it->chi2();
00269       flags = it->recoFlag();
00270 
00271       edm::ESHandle<EcalSeverityLevelAlgo> sevlv;
00272       iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
00273       severity = sevlv->severityLevel( id, rechits);
00274     }
00275     bool isNotSpike = ((recoPhoton->isEB() && (severity!=3 && severity!=4 ) && (flags != 2) ) || recoPhoton->isEE());
00276     if (!isNotSpike) continue;  // move on to next photon
00277     // END of determining ECAL Spikes
00278 
00279     bool pho_current_passPhotonID = false;
00280     bool pho_current_isEB = recoPhoton->isEB();
00281     bool pho_current_isEE = recoPhoton->isEE();
00282 
00283     if ( pho_current_isEB && (recoPhoton->sigmaIetaIeta() < 0.01 || recoPhoton->hadronicOverEm() < 0.05) ) {
00284       // Photon object in barrel passes photon ID
00285       pho_current_passPhotonID = true;
00286       photon_count_bar++;
00287     } else if ( pho_current_isEE && (recoPhoton->hadronicOverEm() < 0.05) ) {
00288       // Photon object in endcap passes photon ID
00289       pho_current_passPhotonID = true;
00290       photon_count_end++;
00291     }
00292 
00293     if (!found_lead_pho) {
00294       found_lead_pho = true;
00295       photon_passPhotonID = pho_current_passPhotonID;
00296       photon_et  = recoPhoton->et();
00297       photon_eta = recoPhoton->eta();
00298       photon_phi = recoPhoton->phi();
00299     }
00300   }
00301   
00302   // If user requires a photon to be found, but none is, return.
00303   //   theRequirePhotonFound should pretty much always be set to 'True'
00304   //    except when running on qcd monte carlo just to see the jets.
00305   if ( theRequirePhotonFound && (!photon_passPhotonID || photon_et<theMinPhotonEt) ) return;
00306 
00307 
00309   // Find the highest et jet
00310   Handle<CaloJetCollection> caloJetCollection;
00311   iEvent.getByLabel (theCaloJetCollectionLabel,caloJetCollection);
00312   if (!caloJetCollection.isValid()) return;
00313 
00314   float jet_pt    = -8.0;
00315   float jet_eta   = -8.0;
00316   float jet_phi   = -8.0;
00317   int   jet_count = 0;
00318   float jet2_pt   = -9.0;
00319   float jet2_eta  = -9.0;
00320   float jet2_phi  = -9.0;
00321   // Assumption: jets are ordered by Et
00322   for (CaloJetCollection::const_iterator i_calojet = caloJetCollection->begin(); i_calojet != caloJetCollection->end(); i_calojet++) {
00323 
00324     float jet_current_pt = i_calojet->pt();
00325 
00326     // don't care about jets that overlap with the lead photon
00327     if ( deltaR(i_calojet->eta(), i_calojet->phi(), photon_eta, photon_phi) < 0.5 ) continue;
00328     // stop looping over jets once we get to too low Et
00329     if (jet_current_pt < theMinCaloJetPt) break;
00330 
00331     jet_count++;
00332     if (jet_current_pt > jet_pt) {
00333       jet2_pt  = jet_pt;  // 2nd highest jet get's et from current highest
00334       jet2_eta = jet_eta;
00335       jet2_phi = jet_phi;
00336       jet_pt   = jet_current_pt; // current highest jet gets et from the new highest
00337       jet_eta  = i_calojet->eta();
00338       jet_phi  = i_calojet->phi();
00339     } else if (jet_current_pt > jet2_pt) {
00340       jet2_pt  = jet_current_pt;
00341       jet2_eta = i_calojet->eta();
00342       jet2_phi = i_calojet->phi();
00343     }
00344   }
00346 
00347 
00349   // Fill histograms if a jet found
00350   // NOTE: if a photon was required to be found, but wasn't
00351   //        we wouldn't have made it to this point in the code
00352   if ( jet_pt > 0.0 ) {
00353 
00354     // Photon Plots
00355     h_photon_et       ->Fill( photon_et  );
00356     h_photon_eta      ->Fill( photon_eta );
00357     h_photon_count_bar->Fill( photon_count_bar );
00358     h_photon_count_end->Fill( photon_count_end );
00359 
00360     // Photon Et hists for different orientations to the jet
00361     if ( fabs(photon_eta)<1.45 && photon_passPhotonID ) {  // Lead photon is in barrel
00362       if (fabs(jet_eta)<1.45){                          //   jet is in barrel
00363         if (photon_eta*jet_eta>0) {
00364           h_photon_et_jetcs->Fill(photon_et);
00365         } else {
00366           h_photon_et_jetco->Fill(photon_et);
00367         }
00368       } else if (jet_eta>1.55 && jet_eta<2.5) {         // jet is in endcap
00369         if (photon_eta*jet_eta>0) {
00370           h_photon_et_jetfs->Fill(photon_et);
00371         } else {
00372           h_photon_et_jetfo->Fill(photon_et);
00373         }
00374       }
00375     } // END of Lead Photon is in Barrel
00376 
00377     // Jet Plots
00378     h_jet_pt       ->Fill( jet_pt     );
00379     h_jet_eta      ->Fill( jet_eta    );
00380     h_jet_count    ->Fill( jet_count  );
00381     h_deltaPhi_photon_jet   ->Fill( abs(deltaPhi(photon_phi, jet_phi)) );
00382     if ( abs(deltaPhi(photon_phi,jet_phi))>2.8 ) h_deltaEt_photon_jet->Fill( (photon_et-jet_pt)/photon_et );
00383 
00384     // 2nd Highest Jet Plots
00385     if ( jet2_pt  > 0.0 ) {
00386       h_jet2_pt             ->Fill( jet2_pt  );
00387       h_jet2_eta            ->Fill( jet2_eta );
00388       h_jet2_ptOverPhotonEt ->Fill( jet2_pt/photon_et );
00389       h_deltaPhi_photon_jet2->Fill( abs(deltaPhi(photon_phi, jet2_phi)) );
00390       h_deltaPhi_jet_jet2   ->Fill( abs(deltaPhi(   jet_phi, jet2_phi)) );
00391       h_deltaR_jet_jet2     ->Fill( deltaR(   jet_eta,    jet_phi, jet2_eta, jet2_phi) );
00392       h_deltaR_photon_jet2  ->Fill( deltaR(photon_eta, photon_phi, jet2_eta, jet2_phi) );
00393     }
00394   } 
00395   // End of Filling histograms
00397 }
00398 
00399 
00400 void QcdPhotonsDQM::endJob(void) {}
00401 
00402 void QcdPhotonsDQM::endRun(const edm::Run& run, const edm::EventSetup& es) {
00403   if (num_events_in_run>0) { 
00404     h_triggers_passed->getTH1F()->Scale(1.0/num_events_in_run);
00405   }
00406   h_photon_et_ratio_co_cs->getTH1F()->Divide( h_photon_et_jetco->getTH1F(), h_photon_et_jetcs->getTH1F() );
00407   h_photon_et_ratio_fo_fs->getTH1F()->Divide( h_photon_et_jetfo->getTH1F(), h_photon_et_jetfs->getTH1F() );
00408   h_photon_et_ratio_cs_fs->getTH1F()->Divide( h_photon_et_jetcs->getTH1F(), h_photon_et_jetfs->getTH1F() );
00409   h_photon_et_ratio_co_fs->getTH1F()->Divide( h_photon_et_jetco->getTH1F(), h_photon_et_jetfs->getTH1F() );
00410   h_photon_et_ratio_cs_fo->getTH1F()->Divide( h_photon_et_jetcs->getTH1F(), h_photon_et_jetfo->getTH1F() );
00411   h_photon_et_ratio_co_fo->getTH1F()->Divide( h_photon_et_jetco->getTH1F(), h_photon_et_jetfo->getTH1F() );
00412 }
00413