CMS 3D CMS Logo

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