CMS 3D CMS Logo

HGCalTriggerValidator.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <iostream>
4 #include <string>
5 #include <vector>
6 
7 // user include files
13 
18 
21 
27 
32 
33 //
34 // class declaration
35 //
36 
37 struct Histograms {
38  //histogram tc related
49 
50  //histogram cl related
59 
60  //histogram multicl related
67  // cluster shower shapes
83 
84  //histogram tower related
94 };
95 
96 class HGCalTriggerValidator : public DQMGlobalEDAnalyzer<Histograms> {
97 public:
98  explicit HGCalTriggerValidator(const edm::ParameterSet &);
99  ~HGCalTriggerValidator() override = default;
100 
101 private:
102  void dqmBeginRun(edm::Run const &, edm::EventSetup const &, Histograms &) const override;
103  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &, Histograms &) const override;
104  void dqmAnalyze(edm::Event const &, edm::EventSetup const &, Histograms const &) const override;
105 
106 private:
107  // ----------member data ---------------------------
113 
114  std::unique_ptr<HGCalTriggerClusterIdentificationBase> id_;
115 
116  std::shared_ptr<HGCalTriggerTools> triggerTools_;
117 };
118 
120  : trigger_cells_token_{consumes<l1t::HGCalTriggerCellBxCollection>(
121  iConfig.getParameter<edm::InputTag>("TriggerCells"))},
122  clusters_token_{consumes<l1t::HGCalClusterBxCollection>(iConfig.getParameter<edm::InputTag>("Clusters"))},
123  multiclusters_token_{
124  consumes<l1t::HGCalMulticlusterBxCollection>(iConfig.getParameter<edm::InputTag>("Multiclusters"))},
125  towers_token_{consumes<l1t::HGCalTowerBxCollection>(iConfig.getParameter<edm::InputTag>("Towers"))},
126  triggerGeomToken_(esConsumes<HGCalTriggerGeometryBase, CaloGeometryRecord, edm::Transition::BeginRun>()),
127  id_{HGCalTriggerClusterIdentificationFactory::get()->create("HGCalTriggerClusterIdentificationBDT")} {
128  id_->initialize(iConfig.getParameter<edm::ParameterSet>("EGIdentification"));
129 
130  triggerTools_ = std::make_shared<HGCalTriggerTools>();
131 }
132 
134  edm::EventSetup const &iSetup,
135  Histograms &histograms) const {
136  auto const triggerGeometry = iSetup.getHandle(triggerGeomToken_);
137  triggerTools_->setGeometry(triggerGeometry.product());
138 }
139 
141  edm::Run const &,
142  edm::EventSetup const &iSetup,
143  Histograms &histograms) const {
144  iBooker.cd();
145  iBooker.setCurrentFolder("HGCALTPG");
146 
147  //initiating histograms
148  // trigger cells
149  histograms.h_tc_n_ = iBooker.book1D("tc_n", "trigger cell number; number", 400, 0, 400);
150  histograms.h_tc_mipPt_ = iBooker.book1D("tc_mipPt", "trigger cell mipPt; mipPt", 400, 0, 400);
151  histograms.h_tc_pt_ = iBooker.book1D("tc_pt", "trigger cell pt; pt [GeV]", 15, 0, 15);
152  histograms.h_tc_energy_ = iBooker.book1D("tc_energy", "trigger cell energy; energy [GeV]", 70, 0, 70);
153  histograms.h_tc_eta_ = iBooker.book1D("tc_eta", "trigger cell eta; eta", 60, -3.14, 3.14);
154  histograms.h_tc_phi_ = iBooker.book1D("tc_phi", "trigger cell phi; phi", 60, -3.14, 3.14);
155  histograms.h_tc_x_ = iBooker.book1D("tc_x", "trigger cell x; x [cm]", 500, -250, 250);
156  histograms.h_tc_y_ = iBooker.book1D("tc_y", "trigger cell y; y [cm]", 500, -250, 250);
157  histograms.h_tc_z_ = iBooker.book1D("tc_z", "trigger cell z; z [cm]", 1100, -550, 550);
158  histograms.h_tc_layer_ = iBooker.book1D("tc_layer", "trigger cell layer; layer", 50, 0, 50);
159 
160  // cluster 2D histograms
161  histograms.h_cl_n_ = iBooker.book1D("cl_n", "cluster2D number; number", 80, 0, 80);
162  histograms.h_cl_mipPt_ = iBooker.book1D("cl_mipPt", "cluster2D mipPt; mipPt", 600, 0, 600);
163  histograms.h_cl_pt_ = iBooker.book1D("cl_pt", "cluster2D pt; pt [GeV]", 20, 0, 20);
164  histograms.h_cl_energy_ = iBooker.book1D("cl_energy", "cluster2D energy; energy [GeV]", 80, 0, 80);
165  histograms.h_cl_eta_ = iBooker.book1D("cl_eta", "cluster2D eta; eta", 60, -3.14, 3.14);
166  histograms.h_cl_phi_ = iBooker.book1D("cl_phi", "cluster2D phi; phi", 60, -3.14, 3.14);
167  histograms.h_cl_cells_n_ = iBooker.book1D("cl_cells_n", "cluster2D cells_n; cells_n", 16, 0, 16);
168  histograms.h_cl_layer_ = iBooker.book1D("cl_layer", "cluster2D layer; layer", 50, 0, 50);
169 
170  // multiclusters
171  histograms.h_cl3d_n_ = iBooker.book1D("cl3d_n", "cl3duster3D number; number", 12, 0, 12);
172  histograms.h_cl3d_pt_ = iBooker.book1D("cl3d_pt", "cl3duster3D pt; pt [GeV]", 50, 0, 50);
173  histograms.h_cl3d_energy_ = iBooker.book1D("cl3d_energy", "cl3duster3D energy; energy [GeV]", 80, 0, 80);
174  histograms.h_cl3d_eta_ = iBooker.book1D("cl3d_eta", "cl3duster3D eta; eta", 60, -3.14, 3.14);
175  histograms.h_cl3d_phi_ = iBooker.book1D("cl3d_phi", "cl3duster3D phi; phi", 60, -3.14, 3.14);
176  histograms.h_cl3d_clusters_n_ = iBooker.book1D("cl3d_clusters_n", "cl3duster3D clusters_n; clusters_n", 30, 0, 30);
177  // cluster shower shapes
178  histograms.h_cl3d_showerlength_ =
179  iBooker.book1D("cl3d_showerlength", "cl3duster3D showerlength; showerlength", 50, 0, 50);
180  histograms.h_cl3d_coreshowerlength_ =
181  iBooker.book1D("cl3d_coreshowerlength", "cl3duster3D coreshowerlength; coreshowerlength", 16, 0, 16);
182  histograms.h_cl3d_firstlayer_ = iBooker.book1D("cl3d_firstlayer", "cl3duster3D firstlayer; firstlayer", 50, 0, 50);
183  histograms.h_cl3d_maxlayer_ = iBooker.book1D("cl3d_maxlayer", "cl3duster3D maxlayer; maxlayer", 50, 0, 50);
184  histograms.h_cl3d_seetot_ = iBooker.book1D("cl3d_seetot", "cl3duster3D seetot; seetot", 50, 0, 0.05);
185  histograms.h_cl3d_seemax_ = iBooker.book1D("cl3d_seemax", "cl3duster3D seemax; seemax", 40, 0, 0.04);
186  histograms.h_cl3d_spptot_ = iBooker.book1D("cl3d_spptot", "cl3duster3D spptot; spptot", 800, 0, 0.08);
187  histograms.h_cl3d_sppmax_ = iBooker.book1D("cl3d_sppmax", "cl3duster3D sppmax; sppmax", 800, 0, 0.08);
188  histograms.h_cl3d_szz_ = iBooker.book1D("cl3d_szz", "cl3duster3D szz; szz", 50, 0, 50);
189  histograms.h_cl3d_srrtot_ = iBooker.book1D("cl3d_srrtot", "cl3duster3D srrtot; srrtot", 800, 0, 0.008);
190  histograms.h_cl3d_srrmax_ = iBooker.book1D("cl3d_srrmax", "cl3duster3D srrmax; srrmax", 900, 0, 0.009);
191  histograms.h_cl3d_srrmean_ = iBooker.book1D("cl3d_srrmean", "cl3duster3D srrmean; srrmean", 800, 0, 0.008);
192  histograms.h_cl3d_emaxe_ = iBooker.book1D("cl3d_emaxe", "cl3duster3D emaxe; emaxe", 15, 0, 1.5);
193  histograms.h_cl3d_bdteg_ = iBooker.book1D("cl3d_bdteg", "cl3duster3D bdteg; bdteg", 30, -0.7, 0.4);
194  histograms.h_cl3d_quality_ = iBooker.book1D("cl3d_quality", "cl3duster3D quality; quality", 20, 0, 2);
195 
196  // towers
197  histograms.h_tower_n_ = iBooker.book1D("tower_n", "tower n; number", 400, 1200, 1600);
198  histograms.h_tower_pt_ = iBooker.book1D("tower_pt", "tower pt; pt [GeV]", 50, 0, 50);
199  histograms.h_tower_energy_ = iBooker.book1D("tower_energy", "tower energy; energy [GeV]", 200, 0, 200);
200  histograms.h_tower_eta_ = iBooker.book1D("tower_eta", "tower eta; eta", 60, -3.14, 3.14);
201  histograms.h_tower_phi_ = iBooker.book1D("tower_phi", "tower phi; phi", 60, -3.14, 3.14);
202  histograms.h_tower_etEm_ = iBooker.book1D("tower_etEm", "tower etEm; etEm", 50, 0, 50);
203  histograms.h_tower_etHad_ = iBooker.book1D("tower_etHad", "tower etHad; etHad", 30, 0, 0.3);
204  histograms.h_tower_iEta_ = iBooker.book1D("tower_iEta", "tower iEta; iEta", 20, 0, 20);
205  histograms.h_tower_iPhi_ = iBooker.book1D("tower_iPhi", "tower iPhi; iPhi", 80, 0, 80);
206 }
207 
209  edm::EventSetup const &iSetup,
210  Histograms const &histograms) const {
211  int tc_n = 0;
212  int cl_n = 0;
213  int cl3d_n = 0;
214  int tower_n = 0;
215 
216  // retrieve trigger cells
218  iEvent.getByToken(trigger_cells_token_, trigger_cells_h);
219  const l1t::HGCalTriggerCellBxCollection &trigger_cells = *trigger_cells_h;
220 
221  if (trigger_cells_h.isValid()) {
222  for (auto tc_itr = trigger_cells.begin(0); tc_itr != trigger_cells.end(0); tc_itr++) {
223  tc_n++;
224  DetId id(tc_itr->detId());
225  histograms.h_tc_pt_->Fill(tc_itr->pt());
226  histograms.h_tc_mipPt_->Fill(tc_itr->mipPt());
227  histograms.h_tc_energy_->Fill(tc_itr->energy());
228  histograms.h_tc_eta_->Fill(tc_itr->eta());
229  histograms.h_tc_phi_->Fill(tc_itr->phi());
230  histograms.h_tc_x_->Fill(tc_itr->position().x());
231  histograms.h_tc_y_->Fill(tc_itr->position().y());
232  histograms.h_tc_z_->Fill(tc_itr->position().z());
233  histograms.h_tc_layer_->Fill(triggerTools_->layerWithOffset(id));
234  }
235  }
236  histograms.h_tc_n_->Fill(tc_n);
237 
238  // retrieve clusters
240  iEvent.getByToken(clusters_token_, clusters_h);
241  const l1t::HGCalClusterBxCollection &clusters = *clusters_h;
242 
243  if (clusters_h.isValid()) {
244  for (auto cl_itr = clusters.begin(0); cl_itr != clusters.end(0); cl_itr++) {
245  cl_n++;
246  histograms.h_cl_mipPt_->Fill(cl_itr->mipPt());
247  histograms.h_cl_pt_->Fill(cl_itr->pt());
248  histograms.h_cl_energy_->Fill(cl_itr->energy());
249  histograms.h_cl_eta_->Fill(cl_itr->eta());
250  histograms.h_cl_phi_->Fill(cl_itr->phi());
251  histograms.h_cl_layer_->Fill(triggerTools_->layerWithOffset(cl_itr->detId()));
252  histograms.h_cl_cells_n_->Fill(cl_itr->constituents().size());
253  }
254  }
255  histograms.h_cl_n_->Fill(cl_n);
256 
257  // retrieve clusters 3D
259  iEvent.getByToken(multiclusters_token_, multiclusters_h);
260  const l1t::HGCalMulticlusterBxCollection &multiclusters = *multiclusters_h;
261 
262  if (multiclusters_h.isValid()) {
263  for (auto cl3d_itr = multiclusters.begin(0); cl3d_itr != multiclusters.end(0); cl3d_itr++) {
264  cl3d_n++;
265  histograms.h_cl3d_pt_->Fill(cl3d_itr->pt());
266  histograms.h_cl3d_energy_->Fill(cl3d_itr->energy());
267  histograms.h_cl3d_eta_->Fill(cl3d_itr->eta());
268  histograms.h_cl3d_phi_->Fill(cl3d_itr->phi());
269  histograms.h_cl3d_clusters_n_->Fill(cl3d_itr->constituents().size());
270  // cluster shower shapes
271  histograms.h_cl3d_showerlength_->Fill(cl3d_itr->showerLength());
272  histograms.h_cl3d_coreshowerlength_->Fill(cl3d_itr->coreShowerLength());
273  histograms.h_cl3d_firstlayer_->Fill(cl3d_itr->firstLayer());
274  histograms.h_cl3d_maxlayer_->Fill(cl3d_itr->maxLayer());
275  histograms.h_cl3d_seetot_->Fill(cl3d_itr->sigmaEtaEtaTot());
276  histograms.h_cl3d_seemax_->Fill(cl3d_itr->sigmaEtaEtaMax());
277  histograms.h_cl3d_spptot_->Fill(cl3d_itr->sigmaPhiPhiTot());
278  histograms.h_cl3d_sppmax_->Fill(cl3d_itr->sigmaPhiPhiMax());
279  histograms.h_cl3d_szz_->Fill(cl3d_itr->sigmaZZ());
280  histograms.h_cl3d_srrtot_->Fill(cl3d_itr->sigmaRRTot());
281  histograms.h_cl3d_srrmax_->Fill(cl3d_itr->sigmaRRMax());
282  histograms.h_cl3d_srrmean_->Fill(cl3d_itr->sigmaRRMean());
283  histograms.h_cl3d_emaxe_->Fill(cl3d_itr->eMax() / cl3d_itr->energy());
284  histograms.h_cl3d_bdteg_->Fill(id_->value(*cl3d_itr));
285  histograms.h_cl3d_quality_->Fill(cl3d_itr->hwQual());
286  }
287  }
288  histograms.h_cl3d_n_->Fill(cl3d_n);
289 
290  // retrieve towers
292  iEvent.getByToken(towers_token_, towers_h);
293  const l1t::HGCalTowerBxCollection &towers = *towers_h;
294 
295  if (towers_h.isValid()) {
296  for (auto tower_itr = towers.begin(0); tower_itr != towers.end(0); tower_itr++) {
297  tower_n++;
298  histograms.h_tower_pt_->Fill(tower_itr->pt());
299  histograms.h_tower_energy_->Fill(tower_itr->energy());
300  histograms.h_tower_eta_->Fill(tower_itr->eta());
301  histograms.h_tower_phi_->Fill(tower_itr->phi());
302  histograms.h_tower_etEm_->Fill(tower_itr->etEm());
303  histograms.h_tower_etHad_->Fill(tower_itr->etHad());
304  histograms.h_tower_iEta_->Fill(tower_itr->id().iEta());
305  histograms.h_tower_iPhi_->Fill(tower_itr->id().iPhi());
306  }
307  }
308  histograms.h_tower_n_->Fill(tower_n);
309 }
310 
312 
313 //define this as a plug-in
dqm::reco::MonitorElement * h_tower_etHad_
const edm::EDGetToken clusters_token_
dqm::reco::MonitorElement * h_cl3d_pt_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
dqm::reco::MonitorElement * h_cl_cells_n_
dqm::reco::MonitorElement * h_tower_energy_
dqm::reco::MonitorElement * h_tc_z_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
dqm::reco::MonitorElement * h_cl3d_srrtot_
HGCalTriggerValidator(const edm::ParameterSet &)
dqm::reco::MonitorElement * h_cl_phi_
dqm::reco::MonitorElement * h_cl3d_clusters_n_
dqm::reco::MonitorElement * h_cl3d_n_
dqm::reco::MonitorElement * h_cl3d_eta_
~HGCalTriggerValidator() override=default
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &, Histograms &) const override
dqm::reco::MonitorElement * h_tower_n_
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, Histograms const &) const override
void dqmBeginRun(edm::Run const &, edm::EventSetup const &, Histograms &) const override
dqm::reco::MonitorElement * h_tc_x_
dqm::reco::MonitorElement * h_cl_eta_
dqm::reco::MonitorElement * h_cl_energy_
dqm::reco::MonitorElement * h_cl3d_srrmax_
dqm::reco::MonitorElement * h_tc_layer_
dqm::reco::MonitorElement * h_cl3d_coreshowerlength_
const_iterator begin(int bx) const
std::unique_ptr< HGCalTriggerClusterIdentificationBase > id_
dqm::reco::MonitorElement * h_cl3d_szz_
dqm::reco::MonitorElement * h_tc_mipPt_
int iEvent
Definition: GenABIO.cc:224
dqm::reco::MonitorElement * h_tc_phi_
dqm::reco::MonitorElement * h_cl3d_seemax_
dqm::reco::MonitorElement * h_tower_eta_
dqm::reco::MonitorElement * h_tc_pt_
dqm::reco::MonitorElement * h_cl_layer_
dqm::reco::MonitorElement * h_cl_mipPt_
dqm::reco::MonitorElement * h_cl3d_spptot_
const edm::EDGetToken towers_token_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
dqm::reco::MonitorElement * h_cl_n_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
dqm::reco::MonitorElement * h_cl3d_sppmax_
dqm::reco::MonitorElement * h_cl_pt_
dqm::reco::MonitorElement * h_cl3d_bdteg_
Definition: DetId.h:17
dqm::reco::MonitorElement * h_tower_iEta_
dqm::reco::MonitorElement * h_cl3d_showerlength_
const edm::EDGetToken multiclusters_token_
dqm::reco::MonitorElement * h_cl3d_quality_
dqm::reco::MonitorElement * h_tc_energy_
dqm::reco::MonitorElement * h_cl3d_firstlayer_
bool isValid() const
Definition: HandleBase.h:70
dqm::reco::MonitorElement * h_cl3d_energy_
const edm::ESGetToken< HGCalTriggerGeometryBase, CaloGeometryRecord > triggerGeomToken_
const_iterator end(int bx) const
dqm::reco::MonitorElement * h_cl3d_maxlayer_
std::shared_ptr< HGCalTriggerTools > triggerTools_
dqm::reco::MonitorElement * h_tower_iPhi_
const edm::EDGetToken trigger_cells_token_
dqm::reco::MonitorElement * h_tc_y_
dqm::reco::MonitorElement * h_tower_etEm_
dqm::reco::MonitorElement * h_tc_n_
dqm::reco::MonitorElement * h_tc_eta_
dqm::reco::MonitorElement * h_cl3d_emaxe_
#define get
dqm::reco::MonitorElement * h_cl3d_srrmean_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
dqm::reco::MonitorElement * h_tower_pt_
dqm::reco::MonitorElement * h_tower_phi_
dqm::reco::MonitorElement * h_cl3d_seetot_
Definition: Run.h:45
dqm::reco::MonitorElement * h_cl3d_phi_