Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:15

0001 /*
0002  * \file DQMHcalPhiSymAlCaReco.cc
0003  *
0004  * \author Olga Kodolova
0005  *
0006  *
0007  *
0008  * Description: Monitoring of Phi Symmetry Calibration Stream
0009  */
0010 
0011 // work on collections
0012 
0013 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
0014 #include "DataFormats/DetId/interface/DetId.h"
0015 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0016 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
0017 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
0018 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0019 #include "DataFormats/JetReco/interface/CaloJet.h"
0020 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0021 
0022 #include "FWCore/Framework/interface/Frameworkfwd.h"
0023 #include "FWCore/Framework/interface/Event.h"
0024 #include "FWCore/Framework/interface/EventSetup.h"
0025 #include "FWCore/Framework/interface/MakerMacros.h"
0026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0028 #include "FWCore/ServiceRegistry/interface/Service.h"
0029 
0030 // DQM include files
0031 
0032 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0033 #include "DQMServices/Core/interface/DQMStore.h"
0034 #include "DQMServices/Core/interface/DQMStore.h"
0035 
0036 #include <string>
0037 
0038 class DQMHcalDiJetsAlCaReco : public DQMEDAnalyzer {
0039 public:
0040   DQMHcalDiJetsAlCaReco(const edm::ParameterSet &);
0041   ~DQMHcalDiJetsAlCaReco() override;
0042 
0043 protected:
0044   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0045   void analyze(const edm::Event &e, const edm::EventSetup &c) override;
0046 
0047 private:
0048   int eventCounter_;
0049 
0050   //
0051   // Monitor elements
0052   //
0053   MonitorElement *hiDistrRecHitEnergyEBEE_;
0054   MonitorElement *hiDistrRecHitEnergyHBHE_;
0055   MonitorElement *hiDistrRecHitEnergyHF_;
0056   MonitorElement *hiDistrRecHitEnergyHO_;
0057 
0058   MonitorElement *hiDistrProbeJetEnergy_;
0059   MonitorElement *hiDistrProbeJetEta_;
0060   MonitorElement *hiDistrProbeJetPhi_;
0061 
0062   MonitorElement *hiDistrTagJetEnergy_;
0063   MonitorElement *hiDistrTagJetEta_;
0064   MonitorElement *hiDistrTagJetPhi_;
0065 
0066   MonitorElement *hiDistrEtThirdJet_;
0067 
0068   /// object to monitor
0069   edm::EDGetTokenT<reco::CaloJetCollection> jets_;
0070   edm::EDGetTokenT<EcalRecHitCollection> ec_;
0071   edm::EDGetTokenT<HBHERecHitCollection> hbhe_;
0072   edm::EDGetTokenT<HORecHitCollection> ho_;
0073   edm::EDGetTokenT<HFRecHitCollection> hf_;
0074 
0075   /// DQM folder name
0076   std::string folderName_;
0077 
0078   /// Write to file
0079   bool saveToFile_;
0080 
0081   /// Output file name if required
0082   std::string fileName_;
0083 
0084   bool allowMissingInputs_;
0085 };
0086 
0087 // ******************************************
0088 // constructors
0089 // *****************************************
0090 
0091 DQMHcalDiJetsAlCaReco::DQMHcalDiJetsAlCaReco(const edm::ParameterSet &iConfig) : eventCounter_(0) {
0092   //
0093   // Input from configurator file
0094   //
0095   folderName_ = iConfig.getUntrackedParameter<std::string>("FolderName", "ALCAStreamHcalDiJets");
0096 
0097   jets_ = consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("jetsInput"));
0098   ec_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("ecInput"));
0099   hbhe_ = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"));
0100   ho_ = consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoInput"));
0101   hf_ = consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hfInput"));
0102 
0103   saveToFile_ = iConfig.getUntrackedParameter<bool>("SaveToFile", false);
0104   fileName_ = iConfig.getUntrackedParameter<std::string>("FileName", "MonitorAlCaHcalDiJets.root");
0105 }
0106 
0107 DQMHcalDiJetsAlCaReco::~DQMHcalDiJetsAlCaReco() {}
0108 
0109 //--------------------------------------------------------
0110 void DQMHcalDiJetsAlCaReco::bookHistograms(DQMStore::IBooker &ibooker,
0111                                            edm::Run const & /* iRun*/,
0112                                            edm::EventSetup const & /* iSetup */) {
0113   // create and cd into new folder
0114   ibooker.setCurrentFolder(folderName_);
0115 
0116   // book some histograms 1D
0117   hiDistrRecHitEnergyEBEE_ = ibooker.book1D("RecHitEnergyEBEE", "the number of hits inside jets", 100, 0, 800);
0118   hiDistrRecHitEnergyEBEE_->setAxisTitle("E, GeV", 1);
0119   hiDistrRecHitEnergyEBEE_->setAxisTitle("# rechits", 2);
0120 
0121   hiDistrRecHitEnergyHBHE_ = ibooker.book1D("RecHitEnergyHBHE", "the number of hits inside jets", 100, 0, 800);
0122   hiDistrRecHitEnergyHBHE_->setAxisTitle("E, GeV", 1);
0123   hiDistrRecHitEnergyHBHE_->setAxisTitle("# rechits", 2);
0124 
0125   hiDistrRecHitEnergyHF_ = ibooker.book1D("RecHitEnergyHF", "the number of hits inside jets", 150, 0, 1500);
0126   hiDistrRecHitEnergyHF_->setAxisTitle("E, GeV", 1);
0127   hiDistrRecHitEnergyHF_->setAxisTitle("# rechits", 2);
0128 
0129   hiDistrRecHitEnergyHO_ = ibooker.book1D("RecHitEnergyHO", "the number of hits inside jets", 100, 0, 100);
0130   hiDistrRecHitEnergyHO_->setAxisTitle("E, GeV", 1);
0131   hiDistrRecHitEnergyHO_->setAxisTitle("# rechits", 2);
0132 
0133   hiDistrProbeJetEnergy_ = ibooker.book1D("ProbeJetEnergy", "the energy of probe jets", 250, 0, 2500);
0134   hiDistrProbeJetEnergy_->setAxisTitle("E, GeV", 1);
0135   hiDistrProbeJetEnergy_->setAxisTitle("# jets", 2);
0136 
0137   hiDistrProbeJetEta_ = ibooker.book1D("ProbeJetEta", "the number of probe jets", 100, -5., 5.);
0138   hiDistrProbeJetEta_->setAxisTitle("#eta", 1);
0139   hiDistrProbeJetEta_->setAxisTitle("# jets", 2);
0140 
0141   hiDistrProbeJetPhi_ = ibooker.book1D("ProbeJetPhi", "the number of probe jets", 50, -3.14, 3.14);
0142   hiDistrProbeJetPhi_->setAxisTitle("#phi", 1);
0143   hiDistrProbeJetPhi_->setAxisTitle("# jets", 2);
0144 
0145   hiDistrTagJetEnergy_ = ibooker.book1D("TagJetEnergy", "the energy of tsg jets", 250, 0, 2500);
0146   hiDistrTagJetEnergy_->setAxisTitle("E, GeV", 1);
0147   hiDistrTagJetEnergy_->setAxisTitle("# jets", 2);
0148 
0149   hiDistrTagJetEta_ = ibooker.book1D("TagJetEta", "the number of  tag jets", 100, -5., 5.);
0150   hiDistrTagJetEta_->setAxisTitle("#eta", 1);
0151   hiDistrTagJetEta_->setAxisTitle("# jets", 2);
0152 
0153   hiDistrTagJetPhi_ = ibooker.book1D("TagJetPhi", "the number of tag jets", 50, -3.14, 3.14);
0154   hiDistrTagJetPhi_->setAxisTitle("#phi", 1);
0155   hiDistrTagJetPhi_->setAxisTitle("# jets", 2);
0156 
0157   hiDistrEtThirdJet_ = ibooker.book1D("EtThirdJet", "Et of the third jet", 90, 0, 90);
0158 
0159   //==================================================================================
0160 }
0161 
0162 //-------------------------------------------------------------
0163 
0164 void DQMHcalDiJetsAlCaReco::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0165   eventCounter_++;
0166 
0167   reco::CaloJet jet1, jet2, jet3;
0168   Float_t etVetoJet;
0169 
0170   edm::Handle<reco::CaloJetCollection> jets;
0171   iEvent.getByToken(jets_, jets);
0172 
0173   if (!jets.isValid()) {
0174     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't getjet product!" << std::endl;
0175     return;
0176   }
0177 
0178   if (jets->size() > 1) {
0179     jet1 = (*jets)[0];
0180     jet2 = (*jets)[1];
0181     if (fabs(jet1.eta()) > fabs(jet2.eta())) {
0182       reco::CaloJet jet = jet1;
0183       jet1 = jet2;
0184       jet2 = jet;
0185     }
0186     //     if(fabs(jet1.eta())>eta_1 || (fabs(jet2.eta())-jet_R) < eta_2){
0187     //     return;}
0188   } else {
0189     return;
0190   }
0191 
0192   hiDistrTagJetEnergy_->Fill(jet1.energy());
0193   hiDistrTagJetEta_->Fill(jet1.eta());
0194   hiDistrTagJetPhi_->Fill(jet1.phi());
0195 
0196   hiDistrProbeJetEnergy_->Fill(jet2.energy());
0197   hiDistrProbeJetEta_->Fill(jet2.eta());
0198   hiDistrProbeJetPhi_->Fill(jet2.phi());
0199 
0200   if (jets->size() > 2) {
0201     jet3 = (*jets)[2];
0202     etVetoJet = jet3.et();
0203     hiDistrEtThirdJet_->Fill(etVetoJet);
0204   } else {
0205     etVetoJet = 0.;
0206     hiDistrEtThirdJet_->Fill(etVetoJet);
0207   }
0208 
0209   edm::Handle<EcalRecHitCollection> ec;
0210   iEvent.getByToken(ec_, ec);
0211 
0212   if (!ec.isValid()) {
0213     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ec product!" << std::endl;
0214     return;
0215   }
0216 
0217   for (EcalRecHitCollection::const_iterator ecItr = (*ec).begin(); ecItr != (*ec).end(); ++ecItr) {
0218     hiDistrRecHitEnergyEBEE_->Fill(ecItr->energy());
0219   }
0220 
0221   edm::Handle<HBHERecHitCollection> hbhe;
0222   iEvent.getByToken(hbhe_, hbhe);
0223 
0224   if (!hbhe.isValid()) {
0225     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hbhe product!" << std::endl;
0226     return;
0227   }
0228 
0229   for (HBHERecHitCollection::const_iterator hbheItr = hbhe->begin(); hbheItr != hbhe->end(); hbheItr++) {
0230     hiDistrRecHitEnergyHBHE_->Fill(hbheItr->energy());
0231   }
0232 
0233   edm::Handle<HORecHitCollection> ho;
0234   iEvent.getByToken(ho_, ho);
0235 
0236   if (!ho.isValid()) {
0237     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ho product!" << std::endl;
0238     return;
0239   }
0240 
0241   for (HORecHitCollection::const_iterator hoItr = ho->begin(); hoItr != ho->end(); hoItr++) {
0242     hiDistrRecHitEnergyHO_->Fill(hoItr->energy());
0243   }
0244 
0245   edm::Handle<HFRecHitCollection> hf;
0246   iEvent.getByToken(hf_, hf);
0247 
0248   if (!hf.isValid()) {
0249     LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hf product!" << std::endl;
0250     return;
0251   }
0252 
0253   for (HFRecHitCollection::const_iterator hfItr = hf->begin(); hfItr != hf->end(); hfItr++) {
0254     hiDistrRecHitEnergyHF_->Fill(hfItr->energy());
0255   }
0256 
0257 }  // analyze
0258 
0259 DEFINE_FWK_MODULE(DQMHcalDiJetsAlCaReco);