DQMHcalDiJetsAlCaReco

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259
/*
 * \file DQMHcalPhiSymAlCaReco.cc
 *
 * \author Olga Kodolova
 *
 *
 *
 * Description: Monitoring of Phi Symmetry Calibration Stream
 */

// work on collections

#include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
#include "DataFormats/DetId/interface/DetId.h"
#include "DataFormats/HcalDetId/interface/HcalDetId.h"
#include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
#include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "DataFormats/JetReco/interface/CaloJet.h"
#include "DataFormats/JetReco/interface/CaloJetCollection.h"

#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"

// DQM include files

#include "DQMServices/Core/interface/DQMEDAnalyzer.h"
#include "DQMServices/Core/interface/DQMStore.h"
#include "DQMServices/Core/interface/DQMStore.h"

#include <string>

class DQMHcalDiJetsAlCaReco : public DQMEDAnalyzer {
public:
  DQMHcalDiJetsAlCaReco(const edm::ParameterSet &);
  ~DQMHcalDiJetsAlCaReco() override;

protected:
  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
  void analyze(const edm::Event &e, const edm::EventSetup &c) override;

private:
  int eventCounter_;

  //
  // Monitor elements
  //
  MonitorElement *hiDistrRecHitEnergyEBEE_;
  MonitorElement *hiDistrRecHitEnergyHBHE_;
  MonitorElement *hiDistrRecHitEnergyHF_;
  MonitorElement *hiDistrRecHitEnergyHO_;

  MonitorElement *hiDistrProbeJetEnergy_;
  MonitorElement *hiDistrProbeJetEta_;
  MonitorElement *hiDistrProbeJetPhi_;

  MonitorElement *hiDistrTagJetEnergy_;
  MonitorElement *hiDistrTagJetEta_;
  MonitorElement *hiDistrTagJetPhi_;

  MonitorElement *hiDistrEtThirdJet_;

  /// object to monitor
  edm::EDGetTokenT<reco::CaloJetCollection> jets_;
  edm::EDGetTokenT<EcalRecHitCollection> ec_;
  edm::EDGetTokenT<HBHERecHitCollection> hbhe_;
  edm::EDGetTokenT<HORecHitCollection> ho_;
  edm::EDGetTokenT<HFRecHitCollection> hf_;

  /// DQM folder name
  std::string folderName_;

  /// Write to file
  bool saveToFile_;

  /// Output file name if required
  std::string fileName_;

  bool allowMissingInputs_;
};

// ******************************************
// constructors
// *****************************************

DQMHcalDiJetsAlCaReco::DQMHcalDiJetsAlCaReco(const edm::ParameterSet &iConfig) : eventCounter_(0) {
  //
  // Input from configurator file
  //
  folderName_ = iConfig.getUntrackedParameter<std::string>("FolderName", "ALCAStreamHcalDiJets");

  jets_ = consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("jetsInput"));
  ec_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("ecInput"));
  hbhe_ = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"));
  ho_ = consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoInput"));
  hf_ = consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hfInput"));

  saveToFile_ = iConfig.getUntrackedParameter<bool>("SaveToFile", false);
  fileName_ = iConfig.getUntrackedParameter<std::string>("FileName", "MonitorAlCaHcalDiJets.root");
}

DQMHcalDiJetsAlCaReco::~DQMHcalDiJetsAlCaReco() {}

//--------------------------------------------------------
void DQMHcalDiJetsAlCaReco::bookHistograms(DQMStore::IBooker &ibooker,
                                           edm::Run const & /* iRun*/,
                                           edm::EventSetup const & /* iSetup */) {
  // create and cd into new folder
  ibooker.setCurrentFolder(folderName_);

  // book some histograms 1D
  hiDistrRecHitEnergyEBEE_ = ibooker.book1D("RecHitEnergyEBEE", "the number of hits inside jets", 100, 0, 800);
  hiDistrRecHitEnergyEBEE_->setAxisTitle("E, GeV", 1);
  hiDistrRecHitEnergyEBEE_->setAxisTitle("# rechits", 2);

  hiDistrRecHitEnergyHBHE_ = ibooker.book1D("RecHitEnergyHBHE", "the number of hits inside jets", 100, 0, 800);
  hiDistrRecHitEnergyHBHE_->setAxisTitle("E, GeV", 1);
  hiDistrRecHitEnergyHBHE_->setAxisTitle("# rechits", 2);

  hiDistrRecHitEnergyHF_ = ibooker.book1D("RecHitEnergyHF", "the number of hits inside jets", 150, 0, 1500);
  hiDistrRecHitEnergyHF_->setAxisTitle("E, GeV", 1);
  hiDistrRecHitEnergyHF_->setAxisTitle("# rechits", 2);

  hiDistrRecHitEnergyHO_ = ibooker.book1D("RecHitEnergyHO", "the number of hits inside jets", 100, 0, 100);
  hiDistrRecHitEnergyHO_->setAxisTitle("E, GeV", 1);
  hiDistrRecHitEnergyHO_->setAxisTitle("# rechits", 2);

  hiDistrProbeJetEnergy_ = ibooker.book1D("ProbeJetEnergy", "the energy of probe jets", 250, 0, 2500);
  hiDistrProbeJetEnergy_->setAxisTitle("E, GeV", 1);
  hiDistrProbeJetEnergy_->setAxisTitle("# jets", 2);

  hiDistrProbeJetEta_ = ibooker.book1D("ProbeJetEta", "the number of probe jets", 100, -5., 5.);
  hiDistrProbeJetEta_->setAxisTitle("#eta", 1);
  hiDistrProbeJetEta_->setAxisTitle("# jets", 2);

  hiDistrProbeJetPhi_ = ibooker.book1D("ProbeJetPhi", "the number of probe jets", 50, -3.14, 3.14);
  hiDistrProbeJetPhi_->setAxisTitle("#phi", 1);
  hiDistrProbeJetPhi_->setAxisTitle("# jets", 2);

  hiDistrTagJetEnergy_ = ibooker.book1D("TagJetEnergy", "the energy of tsg jets", 250, 0, 2500);
  hiDistrTagJetEnergy_->setAxisTitle("E, GeV", 1);
  hiDistrTagJetEnergy_->setAxisTitle("# jets", 2);

  hiDistrTagJetEta_ = ibooker.book1D("TagJetEta", "the number of  tag jets", 100, -5., 5.);
  hiDistrTagJetEta_->setAxisTitle("#eta", 1);
  hiDistrTagJetEta_->setAxisTitle("# jets", 2);

  hiDistrTagJetPhi_ = ibooker.book1D("TagJetPhi", "the number of tag jets", 50, -3.14, 3.14);
  hiDistrTagJetPhi_->setAxisTitle("#phi", 1);
  hiDistrTagJetPhi_->setAxisTitle("# jets", 2);

  hiDistrEtThirdJet_ = ibooker.book1D("EtThirdJet", "Et of the third jet", 90, 0, 90);

  //==================================================================================
}

//-------------------------------------------------------------

void DQMHcalDiJetsAlCaReco::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
  eventCounter_++;

  reco::CaloJet jet1, jet2, jet3;
  Float_t etVetoJet;

  edm::Handle<reco::CaloJetCollection> jets;
  iEvent.getByToken(jets_, jets);

  if (!jets.isValid()) {
    LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't getjet product!" << std::endl;
    return;
  }

  if (jets->size() > 1) {
    jet1 = (*jets)[0];
    jet2 = (*jets)[1];
    if (fabs(jet1.eta()) > fabs(jet2.eta())) {
      reco::CaloJet jet = jet1;
      jet1 = jet2;
      jet2 = jet;
    }
    //     if(fabs(jet1.eta())>eta_1 || (fabs(jet2.eta())-jet_R) < eta_2){
    //     return;}
  } else {
    return;
  }

  hiDistrTagJetEnergy_->Fill(jet1.energy());
  hiDistrTagJetEta_->Fill(jet1.eta());
  hiDistrTagJetPhi_->Fill(jet1.phi());

  hiDistrProbeJetEnergy_->Fill(jet2.energy());
  hiDistrProbeJetEta_->Fill(jet2.eta());
  hiDistrProbeJetPhi_->Fill(jet2.phi());

  if (jets->size() > 2) {
    jet3 = (*jets)[2];
    etVetoJet = jet3.et();
    hiDistrEtThirdJet_->Fill(etVetoJet);
  } else {
    etVetoJet = 0.;
    hiDistrEtThirdJet_->Fill(etVetoJet);
  }

  edm::Handle<EcalRecHitCollection> ec;
  iEvent.getByToken(ec_, ec);

  if (!ec.isValid()) {
    LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ec product!" << std::endl;
    return;
  }

  for (EcalRecHitCollection::const_iterator ecItr = (*ec).begin(); ecItr != (*ec).end(); ++ecItr) {
    hiDistrRecHitEnergyEBEE_->Fill(ecItr->energy());
  }

  edm::Handle<HBHERecHitCollection> hbhe;
  iEvent.getByToken(hbhe_, hbhe);

  if (!hbhe.isValid()) {
    LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hbhe product!" << std::endl;
    return;
  }

  for (HBHERecHitCollection::const_iterator hbheItr = hbhe->begin(); hbheItr != hbhe->end(); hbheItr++) {
    hiDistrRecHitEnergyHBHE_->Fill(hbheItr->energy());
  }

  edm::Handle<HORecHitCollection> ho;
  iEvent.getByToken(ho_, ho);

  if (!ho.isValid()) {
    LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ho product!" << std::endl;
    return;
  }

  for (HORecHitCollection::const_iterator hoItr = ho->begin(); hoItr != ho->end(); hoItr++) {
    hiDistrRecHitEnergyHO_->Fill(hoItr->energy());
  }

  edm::Handle<HFRecHitCollection> hf;
  iEvent.getByToken(hf_, hf);

  if (!hf.isValid()) {
    LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hf product!" << std::endl;
    return;
  }

  for (HFRecHitCollection::const_iterator hfItr = hf->begin(); hfItr != hf->end(); hfItr++) {
    hiDistrRecHitEnergyHF_->Fill(hfItr->energy());
  }

}  // analyze

DEFINE_FWK_MODULE(DQMHcalDiJetsAlCaReco);