CMS 3D CMS Logo

DQMHcalDiJetsAlCaReco.cc
Go to the documentation of this file.
1 /*
2  * \file DQMHcalPhiSymAlCaReco.cc
3  *
4  * \author Olga Kodolova
5  *
6  *
7  *
8  * Description: Monitoring of Phi Symmetry Calibration Stream
9  */
10 
11 // work on collections
12 
21 
29 
30 // DQM include files
31 
35 
36 #include <string>
37 
39 public:
41  ~DQMHcalDiJetsAlCaReco() override;
42 
43 protected:
44  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
45  void analyze(const edm::Event &e, const edm::EventSetup &c) override;
46 
47 private:
49 
50  //
51  // Monitor elements
52  //
57 
61 
65 
67 
74 
77 
80 
83 
85 };
86 
87 // ******************************************
88 // constructors
89 // *****************************************
90 
92  //
93  // Input from configurator file
94  //
95  folderName_ = iConfig.getUntrackedParameter<std::string>("FolderName", "ALCAStreamHcalDiJets");
96 
97  jets_ = consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("jetsInput"));
98  ec_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("ecInput"));
99  hbhe_ = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"));
100  ho_ = consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoInput"));
101  hf_ = consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hfInput"));
102 
103  saveToFile_ = iConfig.getUntrackedParameter<bool>("SaveToFile", false);
104  fileName_ = iConfig.getUntrackedParameter<std::string>("FileName", "MonitorAlCaHcalDiJets.root");
105 }
106 
108 
109 //--------------------------------------------------------
111  edm::Run const & /* iRun*/,
112  edm::EventSetup const & /* iSetup */) {
113  // create and cd into new folder
114  ibooker.setCurrentFolder(folderName_);
115 
116  // book some histograms 1D
117  hiDistrRecHitEnergyEBEE_ = ibooker.book1D("RecHitEnergyEBEE", "the number of hits inside jets", 100, 0, 800);
119  hiDistrRecHitEnergyEBEE_->setAxisTitle("# rechits", 2);
120 
121  hiDistrRecHitEnergyHBHE_ = ibooker.book1D("RecHitEnergyHBHE", "the number of hits inside jets", 100, 0, 800);
123  hiDistrRecHitEnergyHBHE_->setAxisTitle("# rechits", 2);
124 
125  hiDistrRecHitEnergyHF_ = ibooker.book1D("RecHitEnergyHF", "the number of hits inside jets", 150, 0, 1500);
126  hiDistrRecHitEnergyHF_->setAxisTitle("E, GeV", 1);
127  hiDistrRecHitEnergyHF_->setAxisTitle("# rechits", 2);
128 
129  hiDistrRecHitEnergyHO_ = ibooker.book1D("RecHitEnergyHO", "the number of hits inside jets", 100, 0, 100);
130  hiDistrRecHitEnergyHO_->setAxisTitle("E, GeV", 1);
131  hiDistrRecHitEnergyHO_->setAxisTitle("# rechits", 2);
132 
133  hiDistrProbeJetEnergy_ = ibooker.book1D("ProbeJetEnergy", "the energy of probe jets", 250, 0, 2500);
134  hiDistrProbeJetEnergy_->setAxisTitle("E, GeV", 1);
135  hiDistrProbeJetEnergy_->setAxisTitle("# jets", 2);
136 
137  hiDistrProbeJetEta_ = ibooker.book1D("ProbeJetEta", "the number of probe jets", 100, -5., 5.);
138  hiDistrProbeJetEta_->setAxisTitle("#eta", 1);
139  hiDistrProbeJetEta_->setAxisTitle("# jets", 2);
140 
141  hiDistrProbeJetPhi_ = ibooker.book1D("ProbeJetPhi", "the number of probe jets", 50, -3.14, 3.14);
142  hiDistrProbeJetPhi_->setAxisTitle("#phi", 1);
143  hiDistrProbeJetPhi_->setAxisTitle("# jets", 2);
144 
145  hiDistrTagJetEnergy_ = ibooker.book1D("TagJetEnergy", "the energy of tsg jets", 250, 0, 2500);
146  hiDistrTagJetEnergy_->setAxisTitle("E, GeV", 1);
147  hiDistrTagJetEnergy_->setAxisTitle("# jets", 2);
148 
149  hiDistrTagJetEta_ = ibooker.book1D("TagJetEta", "the number of tag jets", 100, -5., 5.);
150  hiDistrTagJetEta_->setAxisTitle("#eta", 1);
151  hiDistrTagJetEta_->setAxisTitle("# jets", 2);
152 
153  hiDistrTagJetPhi_ = ibooker.book1D("TagJetPhi", "the number of tag jets", 50, -3.14, 3.14);
154  hiDistrTagJetPhi_->setAxisTitle("#phi", 1);
155  hiDistrTagJetPhi_->setAxisTitle("# jets", 2);
156 
157  hiDistrEtThirdJet_ = ibooker.book1D("EtThirdJet", "Et of the third jet", 90, 0, 90);
158 
159  //==================================================================================
160 }
161 
162 //-------------------------------------------------------------
163 
165  eventCounter_++;
166 
167  reco::CaloJet jet1, jet2, jet3;
168  Float_t etVetoJet;
169 
171  iEvent.getByToken(jets_, jets);
172 
173  if (!jets.isValid()) {
174  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't getjet product!" << std::endl;
175  return;
176  }
177 
178  if (jets->size() > 1) {
179  jet1 = (*jets)[0];
180  jet2 = (*jets)[1];
181  if (fabs(jet1.eta()) > fabs(jet2.eta())) {
182  reco::CaloJet jet = jet1;
183  jet1 = jet2;
184  jet2 = jet;
185  }
186  // if(fabs(jet1.eta())>eta_1 || (fabs(jet2.eta())-jet_R) < eta_2){
187  // return;}
188  } else {
189  return;
190  }
191 
193  hiDistrTagJetEta_->Fill(jet1.eta());
194  hiDistrTagJetPhi_->Fill(jet1.phi());
195 
197  hiDistrProbeJetEta_->Fill(jet2.eta());
198  hiDistrProbeJetPhi_->Fill(jet2.phi());
199 
200  if (jets->size() > 2) {
201  jet3 = (*jets)[2];
202  etVetoJet = jet3.et();
203  hiDistrEtThirdJet_->Fill(etVetoJet);
204  } else {
205  etVetoJet = 0.;
206  hiDistrEtThirdJet_->Fill(etVetoJet);
207  }
208 
210  iEvent.getByToken(ec_, ec);
211 
212  if (!ec.isValid()) {
213  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ec product!" << std::endl;
214  return;
215  }
216 
217  for (EcalRecHitCollection::const_iterator ecItr = (*ec).begin(); ecItr != (*ec).end(); ++ecItr) {
218  hiDistrRecHitEnergyEBEE_->Fill(ecItr->energy());
219  }
220 
222  iEvent.getByToken(hbhe_, hbhe);
223 
224  if (!hbhe.isValid()) {
225  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hbhe product!" << std::endl;
226  return;
227  }
228 
229  for (HBHERecHitCollection::const_iterator hbheItr = hbhe->begin(); hbheItr != hbhe->end(); hbheItr++) {
230  hiDistrRecHitEnergyHBHE_->Fill(hbheItr->energy());
231  }
232 
234  iEvent.getByToken(ho_, ho);
235 
236  if (!ho.isValid()) {
237  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get ho product!" << std::endl;
238  return;
239  }
240 
241  for (HORecHitCollection::const_iterator hoItr = ho->begin(); hoItr != ho->end(); hoItr++) {
242  hiDistrRecHitEnergyHO_->Fill(hoItr->energy());
243  }
244 
246  iEvent.getByToken(hf_, hf);
247 
248  if (!hf.isValid()) {
249  LogDebug("") << "DQMHcalDiJetsAlCaReco: Error! can't get hf product!" << std::endl;
250  return;
251  }
252 
253  for (HFRecHitCollection::const_iterator hfItr = hf->begin(); hfItr != hf->end(); hfItr++) {
254  hiDistrRecHitEnergyHF_->Fill(hfItr->energy());
255  }
256 
257 } // analyze
258 
edm::EDGetTokenT< reco::CaloJetCollection > jets_
object to monitor
MonitorElement * hiDistrTagJetEnergy_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void analyze(const edm::Event &e, const edm::EventSetup &c) override
MonitorElement * hiDistrRecHitEnergyHF_
Jets made from CaloTowers.
Definition: CaloJet.h:27
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::string fileName_
Output file name if required.
MonitorElement * hiDistrTagJetPhi_
MonitorElement * hiDistrTagJetEta_
std::vector< EcalRecHit >::const_iterator const_iterator
MonitorElement * hiDistrRecHitEnergyEBEE_
T getUntrackedParameter(std::string const &, T const &) const
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< HORecHitCollection > ho_
edm::EDGetTokenT< HFRecHitCollection > hf_
MonitorElement * hiDistrProbeJetEta_
DQMHcalDiJetsAlCaReco(const edm::ParameterSet &)
MonitorElement * hiDistrProbeJetPhi_
edm::EDGetTokenT< HBHERecHitCollection > hbhe_
bool isValid() const
Definition: HandleBase.h:70
bool saveToFile_
Write to file.
double et() const final
transverse energy
edm::EDGetTokenT< EcalRecHitCollection > ec_
MonitorElement * hiDistrRecHitEnergyHBHE_
MonitorElement * hiDistrProbeJetEnergy_
MonitorElement * hiDistrRecHitEnergyHO_
std::string folderName_
DQM folder name.
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
double phi() const final
momentum azimuthal angle
MonitorElement * hiDistrEtThirdJet_
Definition: Run.h:45
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
#define LogDebug(id)
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
double energy() const final
energy
double eta() const final
momentum pseudorapidity