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