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
|
/*
* \file DQMHcalIterativePhiSymAlCaReco.cc
*
* \author Sunanda Banerjee
*
*
*
* Description: Monitoring of Iterative Phi Symmetry Calibration Stream
*/
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/Frameworkfwd.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"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
// DQM include files
#include "DQMServices/Core/interface/DQMStore.h"
#include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
// work on collections
#include "DataFormats/HcalDetId/interface/HcalDetId.h"
#include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
#include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
#include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
#include <string>
class DQMHcalIterativePhiSymAlCaReco : public DQMOneEDAnalyzer<> {
public:
DQMHcalIterativePhiSymAlCaReco(const edm::ParameterSet &);
~DQMHcalIterativePhiSymAlCaReco() override;
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
protected:
// void beginRun(const edm::Run& r, const edm::EventSetup& c);
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
void analyze(const edm::Event &e, const edm::EventSetup &c) override;
void dqmEndRun(const edm::Run &r, const edm::EventSetup &c) override {}
private:
static constexpr int maxDepth_ = 7;
//
// Monitor elements
//
MonitorElement *hiDistr2D_[maxDepth_];
MonitorElement *hiDistrHBHEsize1D_;
MonitorElement *hiDistrHFsize1D_;
MonitorElement *hiDistrHOsize1D_;
std::string folderName_;
int hiDistr_y_nbin_;
int hiDistr_x_nbin_;
double hiDistr_y_min_;
double hiDistr_y_max_;
double hiDistr_x_min_;
double hiDistr_x_max_;
/// object to monitor
edm::EDGetTokenT<HBHERecHitCollection> tok_hbhe_;
edm::EDGetTokenT<HORecHitCollection> tok_ho_;
edm::EDGetTokenT<HFRecHitCollection> tok_hf_;
};
// ******************************************
// constructors
// *****************************************
DQMHcalIterativePhiSymAlCaReco::DQMHcalIterativePhiSymAlCaReco(const edm::ParameterSet &ps) {
//
// Input from configurator file
//
folderName_ = ps.getParameter<std::string>("folderName");
// histogram parameters
tok_ho_ = consumes<HORecHitCollection>(ps.getParameter<edm::InputTag>("hoInput"));
tok_hf_ = consumes<HFRecHitCollection>(ps.getParameter<edm::InputTag>("hfInput"));
tok_hbhe_ = consumes<HBHERecHitCollection>(ps.getParameter<edm::InputTag>("hbheInput"));
// Distribution of rechits in iPhi, iEta
hiDistr_y_nbin_ = ps.getUntrackedParameter<int>("hiDistr_y_nbin", 72);
hiDistr_y_min_ = ps.getUntrackedParameter<double>("hiDistr_y_min", 0.5);
hiDistr_y_max_ = ps.getUntrackedParameter<double>("hiDistr_y_max", 72.5);
hiDistr_x_nbin_ = ps.getUntrackedParameter<int>("hiDistr_x_nbin", 83);
hiDistr_x_min_ = ps.getUntrackedParameter<double>("hiDistr_x_min", -41.5);
hiDistr_x_max_ = ps.getUntrackedParameter<double>("hiDistr_x_max", 41.5);
}
DQMHcalIterativePhiSymAlCaReco::~DQMHcalIterativePhiSymAlCaReco() {}
void DQMHcalIterativePhiSymAlCaReco::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("folderName", "ALCAStreamHcalIterativePhiSym");
desc.add<edm::InputTag>("hbheInput", edm::InputTag("hbhereco"));
desc.add<edm::InputTag>("hfInput", edm::InputTag("hfreco"));
desc.add<edm::InputTag>("hoInput", edm::InputTag("horeco"));
desc.addUntracked<int>("hiDistr_y_nbin", 72);
desc.addUntracked<double>("hiDistr_y_min", 0.5);
desc.addUntracked<double>("hiDistr_y_max", 72.5);
desc.addUntracked<int>("hiDistr_x_nbin", 83);
desc.addUntracked<double>("hiDistr_x_min", -41.5);
desc.addUntracked<double>("hiDistr_x_max", 41.5);
descriptions.add("dqmHcalIterativePhiSymAlCaReco", desc);
}
//--------------------------------------------------------
void DQMHcalIterativePhiSymAlCaReco::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
hiDistrHBHEsize1D_ = ibooker.book1D("DistrHBHEsize", "Size of HBHE Collection", 100, 0.0, 10000.0);
hiDistrHFsize1D_ = ibooker.book1D("DistrHFsize", "Size of HF Collection", 100, 0.0, 10000.0);
hiDistrHOsize1D_ = ibooker.book1D("DistrHOsize", "Size of HO Collection", 100, 0.0, 3000.0);
// Eta-phi occupancy
for (int k = 0; k < maxDepth_; ++k) {
char name[20], title[20];
sprintf(name, "MBdepth%d", (k + 1));
sprintf(title, "Depth %d", (k + 1));
hiDistr2D_[k] = ibooker.book2D(
name, title, hiDistr_x_nbin_, hiDistr_x_min_, hiDistr_x_max_, hiDistr_y_nbin_, hiDistr_y_min_, hiDistr_y_max_);
hiDistr2D_[k]->setAxisTitle("i#phi ", 2);
hiDistr2D_[k]->setAxisTitle("i#eta ", 1);
}
}
//--------------------------------------------------------
//-------------------------------------------------------------
void DQMHcalIterativePhiSymAlCaReco::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
// First HBHE RecHits
edm::Handle<HBHERecHitCollection> hbheMB;
iEvent.getByToken(tok_hbhe_, hbheMB);
if (!hbheMB.isValid()) {
edm::LogWarning("DQMHcal") << "DQMHcalIterativePhiSymAlCaReco: Error! can't get HBHE RecHit Collection!";
} else {
const HBHERecHitCollection Hithbhe = *(hbheMB.product());
hiDistrHBHEsize1D_->Fill(Hithbhe.size());
for (HBHERecHitCollection::const_iterator hbheItr = Hithbhe.begin(); hbheItr != Hithbhe.end(); hbheItr++) {
HcalDetId hid = HcalDetId(hbheItr->detid());
int id = hid.depth() - 1;
if ((id >= 0) && (id < maxDepth_))
hiDistr2D_[id]->Fill(hid.ieta(), hid.iphi(), hbheItr->energy());
}
}
// Then HF RecHits
edm::Handle<HFRecHitCollection> hfMB;
iEvent.getByToken(tok_hf_, hfMB);
if (!hfMB.isValid()) {
edm::LogWarning("DQMHcal") << "DQMHcalIterativePhiSymAlCaReco: Error! can't get HF RecHit Collection!";
} else {
const HFRecHitCollection Hithf = *(hfMB.product());
hiDistrHFsize1D_->Fill(Hithf.size());
for (HFRecHitCollection::const_iterator hfItr = Hithf.begin(); hfItr != Hithf.end(); hfItr++) {
HcalDetId hid = HcalDetId(hfItr->detid());
int id = hid.depth() - 1;
if ((id >= 0) && (id < maxDepth_))
hiDistr2D_[id]->Fill(hid.ieta(), hid.iphi(), hfItr->energy());
}
}
// And finally HO RecHits
edm::Handle<HORecHitCollection> hoMB;
iEvent.getByToken(tok_ho_, hoMB);
if (!hoMB.isValid()) {
edm::LogWarning("DQMHcal") << "DQMHcalIterativePhiSymAlCaReco: Error! can't get HO RecHit Collection!";
} else {
const HORecHitCollection Hitho = *(hoMB.product());
hiDistrHOsize1D_->Fill(Hitho.size());
for (HORecHitCollection::const_iterator hoItr = Hitho.begin(); hoItr != Hitho.end(); hoItr++) {
HcalDetId hid = HcalDetId(hoItr->detid());
int id = hid.depth() - 1;
if ((id >= 0) && (id < maxDepth_))
hiDistr2D_[id]->Fill(hid.ieta(), hid.iphi(), hoItr->energy());
}
}
} // analyze
DEFINE_FWK_MODULE(DQMHcalIterativePhiSymAlCaReco);
|