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