CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/DQMOffline/CalibCalo/src/DQMHcalDiJetsAlCaReco.cc

Go to the documentation of this file.
00001 /*
00002  * \file DQMHcalPhiSymAlCaReco.cc
00003  *
00004  * \author Olga Kodolova
00005  *        
00006  * $Date: 2010/04/23 16:35:17 $
00007  * $Revision: 1.7 $
00008  *
00009  *
00010  * Description: Monitoring of Phi Symmetry Calibration Stream  
00011 */
00012 
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/ServiceRegistry/interface/Service.h"
00017 
00018 // DQM include files
00019 
00020 #include "DQMServices/Core/interface/MonitorElement.h"
00021 
00022 // work on collections
00023 
00024 
00025 #include "DataFormats/JetReco/interface/CaloJet.h"
00026 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00027 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00028 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00029 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00030 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
00031 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00032 #include "DataFormats/DetId/interface/DetId.h"
00033 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00034 
00035 
00036 // #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00037 // #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00038 // #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00039 
00040 #include "DQMServices/Core/interface/DQMStore.h"
00041 #include "DQMOffline/CalibCalo/src/DQMHcalDiJetsAlCaReco.h"
00042 
00043 using namespace std;
00044 using namespace edm;
00045 using namespace reco;
00046 
00047 // ******************************************
00048 // constructors
00049 // *****************************************
00050 
00051 DQMHcalDiJetsAlCaReco::DQMHcalDiJetsAlCaReco( const edm::ParameterSet& iConfig ) :
00052 eventCounter_(0)
00053 {
00054   dbe_ = Service<DQMStore>().operator->();
00055 //
00056 // Input from configurator file 
00057 //
00058   folderName_ = iConfig.getUntrackedParameter<string>("FolderName","ALCAStreamHcalDiJets");
00059 
00060 
00061   jets_= iConfig.getParameter<edm::InputTag>("jetsInput");
00062   ec_= iConfig.getParameter<edm::InputTag>("ecInput");
00063   hbhe_= iConfig.getParameter<edm::InputTag>("hbheInput");
00064   ho_= iConfig.getParameter<edm::InputTag>("hoInput");
00065   hf_= iConfig.getParameter<edm::InputTag>("hfInput");
00066   
00067   saveToFile_= iConfig.getUntrackedParameter<bool>("SaveToFile",false);
00068   fileName_=  iConfig.getUntrackedParameter<string>("FileName","MonitorAlCaHcalDiJets.root");
00069     
00070 }
00071 
00072 DQMHcalDiJetsAlCaReco::~DQMHcalDiJetsAlCaReco()
00073 {}
00074 
00075 //--------------------------------------------------------
00076 void DQMHcalDiJetsAlCaReco::beginJob(){
00077    
00078 
00079   // create and cd into new folder
00080   dbe_->setCurrentFolder(folderName_);
00081 
00082   // book some histograms 1D
00083   hiDistrRecHitEnergyEBEE_ = dbe_->book1D("RecHitEnergyEBEE", "the number of hits inside jets", 100, 0, 800);
00084   hiDistrRecHitEnergyEBEE_->setAxisTitle("E, GeV", 1);
00085   hiDistrRecHitEnergyEBEE_->setAxisTitle("# rechits", 2);
00086 
00087   hiDistrRecHitEnergyHBHE_ = dbe_->book1D("RecHitEnergyHBHE", "the number of hits inside jets", 100,0, 800);
00088   hiDistrRecHitEnergyHBHE_->setAxisTitle("E, GeV", 1);
00089   hiDistrRecHitEnergyHBHE_->setAxisTitle("# rechits", 2);
00090 
00091   hiDistrRecHitEnergyHF_ = dbe_->book1D("RecHitEnergyHF", "the number of hits inside jets", 150,0, 1500); 
00092   hiDistrRecHitEnergyHF_->setAxisTitle("E, GeV", 1);
00093   hiDistrRecHitEnergyHF_->setAxisTitle("# rechits", 2);
00094 
00095   hiDistrRecHitEnergyHO_ = dbe_->book1D("RecHitEnergyHO", "the number of hits inside jets", 100,0, 100); 
00096   hiDistrRecHitEnergyHO_->setAxisTitle("E, GeV", 1);
00097   hiDistrRecHitEnergyHO_->setAxisTitle("# rechits", 2);
00098 
00099   hiDistrProbeJetEnergy_ = dbe_->book1D("ProbeJetEnergy", "the energy of probe jets", 250,0, 2500); 
00100   hiDistrProbeJetEnergy_->setAxisTitle("E, GeV", 1);
00101   hiDistrProbeJetEnergy_->setAxisTitle("# jets", 2);
00102 
00103   hiDistrProbeJetEta_ = dbe_->book1D("ProbeJetEta", "the number of probe jets", 100, -5., 5.); 
00104   hiDistrProbeJetEta_->setAxisTitle("#eta", 1);
00105   hiDistrProbeJetEta_->setAxisTitle("# jets", 2);
00106 
00107   hiDistrProbeJetPhi_ = dbe_->book1D("ProbeJetPhi", "the number of probe jets", 50, -3.14, 3.14); 
00108   hiDistrProbeJetPhi_->setAxisTitle("#phi", 1);
00109   hiDistrProbeJetPhi_->setAxisTitle("# jets", 2);
00110 
00111   hiDistrTagJetEnergy_ = dbe_->book1D("TagJetEnergy", "the energy of tsg jets", 250,0, 2500); 
00112   hiDistrTagJetEnergy_->setAxisTitle("E, GeV", 1);
00113   hiDistrTagJetEnergy_->setAxisTitle("# jets", 2);
00114 
00115   hiDistrTagJetEta_ = dbe_->book1D("TagJetEta", "the number of  tag jets", 100, -5., 5.); 
00116   hiDistrTagJetEta_->setAxisTitle("#eta", 1);
00117   hiDistrTagJetEta_->setAxisTitle("# jets", 2);
00118 
00119   hiDistrTagJetPhi_ = dbe_->book1D("TagJetPhi", "the number of tag jets", 50, -3.14, 3.14); 
00120   hiDistrTagJetPhi_->setAxisTitle("#phi", 1);
00121   hiDistrTagJetPhi_->setAxisTitle("# jets", 2);
00122 
00123   hiDistrEtThirdJet_ = dbe_->book1D("EtThirdJet", "Et of the third jet", 90, 0, 90); 
00124 
00125 
00126 //==================================================================================
00127 
00128 
00129 }
00130 
00131 //--------------------------------------------------------
00132 void DQMHcalDiJetsAlCaReco::beginRun(const edm::Run& r, const EventSetup& context) {
00133 
00134 }
00135 
00136 //--------------------------------------------------------
00137 void DQMHcalDiJetsAlCaReco::beginLuminosityBlock(const LuminosityBlock& lumiSeg, 
00138      const EventSetup& context) {
00139   
00140 }
00141 
00142 //-------------------------------------------------------------
00143 
00144 void DQMHcalDiJetsAlCaReco::analyze(const Event& iEvent, 
00145                                const EventSetup& iSetup ){  
00146  
00147 
00148    eventCounter_++;
00149 
00150    CaloJet jet1, jet2, jet3;
00151    Float_t etVetoJet; 
00152 
00153    edm::Handle<CaloJetCollection> jets;
00154    iEvent.getByLabel(jets_,jets);
00155    
00156   if(!jets.isValid()){
00157     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't getjet product!" << std::endl;
00158     return ;
00159   }
00160    
00161    
00162    
00163    if(jets->size()>1){
00164     jet1 = (*jets)[0];
00165     jet2 = (*jets)[1];
00166      if(fabs(jet1.eta())>fabs(jet2.eta())){
00167        CaloJet jet = jet1;
00168        jet1 = jet2;
00169        jet2 = jet;
00170      }
00171      //     if(fabs(jet1.eta())>eta_1 || (fabs(jet2.eta())-jet_R) < eta_2){ return;}
00172    } else {return;}
00173 
00174    hiDistrTagJetEnergy_->Fill(jet1.energy());
00175    hiDistrTagJetEta_->Fill(jet1.eta());
00176    hiDistrTagJetPhi_->Fill(jet1.phi());
00177 
00178    hiDistrProbeJetEnergy_->Fill(jet2.energy());
00179    hiDistrProbeJetEta_->Fill(jet2.eta());
00180    hiDistrProbeJetPhi_->Fill(jet2.phi());
00181 
00182    if(jets->size()>2){
00183      jet3 = (*jets)[2];
00184      etVetoJet = jet3.et();
00185      hiDistrEtThirdJet_->Fill(etVetoJet);
00186    } else { etVetoJet = 0.; hiDistrEtThirdJet_->Fill(etVetoJet); }
00187 
00188 
00189 
00190       Handle<EcalRecHitCollection> ec;
00191       iEvent.getByLabel(ec_,ec);
00192       
00193   if(!ec.isValid()){
00194     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ec product!" << std::endl;
00195     return ;
00196   }
00197       
00198       
00199       for(EcalRecHitCollection::const_iterator ecItr = (*ec).begin();
00200                                                 ecItr != (*ec).end(); ++ecItr)
00201       {
00202         hiDistrRecHitEnergyEBEE_->Fill(ecItr->energy()); 
00203       }
00204 
00205 
00206 
00207 
00208       Handle<HBHERecHitCollection> hbhe;
00209       iEvent.getByLabel(hbhe_, hbhe);
00210 
00211   if(!hbhe.isValid()){
00212     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hbhe product!" << std::endl;
00213     return ;
00214   }
00215 
00216 
00217       for(HBHERecHitCollection::const_iterator hbheItr=hbhe->begin();
00218                                                  hbheItr!=hbhe->end(); hbheItr++)
00219       {
00220         hiDistrRecHitEnergyHBHE_->Fill(hbheItr->energy()); 
00221       }
00222 
00223    
00224       Handle<HORecHitCollection> ho;
00225       iEvent.getByLabel(ho_, ho);
00226 
00227   if(!ho.isValid()){
00228     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ho product!" << std::endl;
00229     return ;
00230   }
00231 
00232 
00233       for(HORecHitCollection::const_iterator hoItr=ho->begin();
00234                                                hoItr!=ho->end(); hoItr++)
00235       {
00236          hiDistrRecHitEnergyHO_->Fill(hoItr->energy());
00237 
00238       }
00239 
00240 
00241 
00242 
00243       Handle<HFRecHitCollection> hf;
00244       iEvent.getByLabel(hf_, hf);
00245 
00246   if(!hf.isValid()){
00247     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hf product!" << std::endl;
00248     return ;
00249   }
00250 
00251 
00252       for(HFRecHitCollection::const_iterator hfItr=hf->begin();
00253                                                hfItr!=hf->end(); hfItr++)
00254       {
00255         hiDistrRecHitEnergyHF_->Fill(hfItr->energy()); 
00256       }
00257         
00258 } //analyze
00259 
00260 
00261 
00262 
00263 //--------------------------------------------------------
00264 void DQMHcalDiJetsAlCaReco::endLuminosityBlock(const LuminosityBlock& lumiSeg, 
00265                                           const EventSetup& context) {
00266 }
00267 //--------------------------------------------------------
00268 void DQMHcalDiJetsAlCaReco::endRun(const Run& r, const EventSetup& context){
00269 
00270 }
00271 //--------------------------------------------------------
00272 void DQMHcalDiJetsAlCaReco::endJob(){
00273   
00274   if (saveToFile_) {
00275      dbe_->save(fileName_);
00276   }
00277   
00278 }
00279 
00280