CMS 3D CMS Logo

HGCalTriggerNtupleHGCTriggerCells.cc
Go to the documentation of this file.
1 
19 
21 public:
24  void initialize(TTree&, const edm::ParameterSet&, edm::ConsumesCollector&&) final;
25  void fill(const edm::Event& e, const edm::EventSetup& es) final;
26 
27 private:
28  double calibrate(double, int, unsigned);
29  void simhits(const edm::Event& e,
30  std::unordered_map<uint32_t, double>& simhits_ee,
31  std::unordered_map<uint32_t, double>& simhits_fh,
32  std::unordered_map<uint32_t, double>& simhits_bh);
33  void clear() final;
34 
36 
43  double keV2fC_;
44  std::vector<double> fcPerMip_;
45  std::vector<double> layerWeights_;
46  std::vector<double> thicknessCorrections_;
48 
49  static const unsigned kPanelOffset_ = 0;
50  static const unsigned kPanelMask_ = 0x7F;
51  static const unsigned kSectorOffset_ = 7;
52  static const unsigned kSectorMask_ = 0x7;
53 
54  int tc_n_;
55  std::vector<uint32_t> tc_id_;
56  std::vector<int> tc_subdet_;
57  std::vector<int> tc_side_;
58  std::vector<int> tc_layer_;
59  std::vector<int> tc_panel_number_;
60  std::vector<int> tc_panel_sector_;
61  std::vector<int> tc_wafer_;
62  std::vector<int> tc_waferu_;
63  std::vector<int> tc_waferv_;
64  std::vector<int> tc_wafertype_;
65  std::vector<int> tc_cell_;
66  std::vector<int> tc_cellu_;
67  std::vector<int> tc_cellv_;
68  std::vector<uint32_t> tc_data_;
69  std::vector<uint32_t> tc_uncompressedCharge_;
70  std::vector<uint32_t> tc_compressedCharge_;
71  std::vector<float> tc_mipPt_;
72  std::vector<float> tc_pt_;
73  std::vector<float> tc_energy_;
74  std::vector<float> tc_simenergy_;
75  std::vector<float> tc_eta_;
76  std::vector<float> tc_phi_;
77  std::vector<float> tc_x_;
78  std::vector<float> tc_y_;
79  std::vector<float> tc_z_;
80  std::vector<uint32_t> tc_cluster_id_;
81  std::vector<uint32_t> tc_multicluster_id_;
82  std::vector<float> tc_multicluster_pt_;
83  std::vector<int> tc_genparticle_index_;
84 
86 };
87 
89 
91  : HGCalTriggerNtupleBase(conf),
92  fill_simenergy_(conf.getParameter<bool>("FillSimEnergy")),
93  fill_truthmap_(conf.getParameter<bool>("FillTruthMap")),
94  filter_cells_in_multiclusters_(conf.getParameter<bool>("FilterCellsInMulticlusters")),
95  keV2fC_(conf.getParameter<double>("keV2fC")),
96  fcPerMip_(conf.getParameter<std::vector<double>>("fcPerMip")),
97  layerWeights_(conf.getParameter<std::vector<double>>("layerWeights")),
98  thicknessCorrections_(conf.getParameter<std::vector<double>>("thicknessCorrections")) {}
99 
101  const edm::ParameterSet& conf,
102  edm::ConsumesCollector&& collector) {
104  collector.consumes<l1t::HGCalTriggerCellBxCollection>(conf.getParameter<edm::InputTag>("TriggerCells"));
106  collector.consumes<l1t::HGCalMulticlusterBxCollection>(conf.getParameter<edm::InputTag>("Multiclusters"));
107 
108  if (fill_simenergy_) {
109  simhits_ee_token_ = collector.consumes<edm::PCaloHitContainer>(conf.getParameter<edm::InputTag>("eeSimHits"));
110  simhits_fh_token_ = collector.consumes<edm::PCaloHitContainer>(conf.getParameter<edm::InputTag>("fhSimHits"));
111  simhits_bh_token_ = collector.consumes<edm::PCaloHitContainer>(conf.getParameter<edm::InputTag>("bhSimHits"));
112  }
113 
114  if (fill_truthmap_)
116  collector.consumes<CaloToCellsMap>(conf.getParameter<edm::InputTag>("caloParticlesToCells"));
117 
118  std::string prefix(conf.getUntrackedParameter<std::string>("Prefix", "tc"));
119 
120  std::string bname;
121  auto withPrefix([&prefix, &bname](char const* vname) -> char const* {
122  bname = prefix + "_" + vname;
123  return bname.c_str();
124  });
125 
126  tree.Branch(withPrefix("n"), &tc_n_, (prefix + "_n/I").c_str());
127  tree.Branch(withPrefix("id"), &tc_id_);
128  tree.Branch(withPrefix("subdet"), &tc_subdet_);
129  tree.Branch(withPrefix("zside"), &tc_side_);
130  tree.Branch(withPrefix("layer"), &tc_layer_);
131  tree.Branch(withPrefix("wafer"), &tc_wafer_);
132  tree.Branch(withPrefix("waferu"), &tc_waferu_);
133  tree.Branch(withPrefix("waferv"), &tc_waferv_);
134  tree.Branch(withPrefix("wafertype"), &tc_wafertype_);
135  tree.Branch(withPrefix("panel_number"), &tc_panel_number_);
136  tree.Branch(withPrefix("panel_sector"), &tc_panel_sector_);
137  tree.Branch(withPrefix("cell"), &tc_cell_);
138  tree.Branch(withPrefix("cellu"), &tc_cellu_);
139  tree.Branch(withPrefix("cellv"), &tc_cellv_);
140  tree.Branch(withPrefix("data"), &tc_data_);
141  tree.Branch(withPrefix("uncompressedCharge"), &tc_uncompressedCharge_);
142  tree.Branch(withPrefix("compressedCharge"), &tc_compressedCharge_);
143  tree.Branch(withPrefix("pt"), &tc_pt_);
144  tree.Branch(withPrefix("mipPt"), &tc_mipPt_);
145  tree.Branch(withPrefix("energy"), &tc_energy_);
146  if (fill_simenergy_)
147  tree.Branch(withPrefix("simenergy"), &tc_simenergy_);
148  tree.Branch(withPrefix("eta"), &tc_eta_);
149  tree.Branch(withPrefix("phi"), &tc_phi_);
150  tree.Branch(withPrefix("x"), &tc_x_);
151  tree.Branch(withPrefix("y"), &tc_y_);
152  tree.Branch(withPrefix("z"), &tc_z_);
153  tree.Branch(withPrefix("cluster_id"), &tc_cluster_id_);
154  tree.Branch(withPrefix("multicluster_id"), &tc_multicluster_id_);
155  tree.Branch(withPrefix("multicluster_pt"), &tc_multicluster_pt_);
156  if (fill_truthmap_)
157  tree.Branch(withPrefix("genparticle_index"), &tc_genparticle_index_);
158 }
159 
161  // retrieve trigger cells
163  e.getByToken(trigger_cells_token_, trigger_cells_h);
164  const l1t::HGCalTriggerCellBxCollection& trigger_cells = *trigger_cells_h;
165 
166  // retrieve clusters
168  e.getByToken(multiclusters_token_, multiclusters_h);
169  const l1t::HGCalMulticlusterBxCollection& multiclusters = *multiclusters_h;
170 
171  // retrieve geometry
172  es.get<CaloGeometryRecord>().get(geometry_);
173 
174  // sim hit association
175  std::unordered_map<uint32_t, double> simhits_ee;
176  std::unordered_map<uint32_t, double> simhits_fh;
177  std::unordered_map<uint32_t, double> simhits_bh;
178  if (fill_simenergy_)
179  simhits(e, simhits_ee, simhits_fh, simhits_bh);
180 
181  edm::Handle<CaloToCellsMap> caloparticles_map_h;
182  std::unordered_map<uint32_t, unsigned> cell_to_genparticle;
183  if (fill_truthmap_) {
184  e.getByToken(caloparticles_map_token_, caloparticles_map_h);
185  for (auto& keyval : *caloparticles_map_h) {
186  for (auto& tcref : keyval.val)
187  cell_to_genparticle.emplace(tcref->detId(), keyval.key->g4Tracks().at(0).genpartIndex() - 1);
188  }
189  }
190 
191  // Associate cells to clusters
192  std::unordered_map<uint32_t, uint32_t> cell2cluster;
193  std::unordered_map<uint32_t, l1t::HGCalMulticlusterBxCollection::const_iterator> cell2multicluster;
194  for (auto mcl_itr = multiclusters.begin(0); mcl_itr != multiclusters.end(0); mcl_itr++) {
195  // loop on 2D clusters inside 3D clusters
196  for (const auto& cl_ptr : mcl_itr->constituents()) {
197  // loop on TC inside 2D clusters
198  for (const auto& tc_ptr : cl_ptr.second->constituents()) {
199  cell2cluster.emplace(tc_ptr.second->detId(), cl_ptr.second->detId());
200  cell2multicluster.emplace(tc_ptr.second->detId(), mcl_itr);
201  }
202  }
203  }
204 
206 
207  clear();
208  for (auto tc_itr = trigger_cells.begin(0); tc_itr != trigger_cells.end(0); tc_itr++) {
209  if (tc_itr->hwPt() > 0) {
210  auto cl_itr = cell2cluster.find(tc_itr->detId());
211  auto mcl_itr = cell2multicluster.find(tc_itr->detId());
212  uint32_t cl_id = (cl_itr != cell2cluster.end() ? cl_itr->second : 0);
213  uint32_t mcl_id = (mcl_itr != cell2multicluster.end() ? mcl_itr->second->detId() : 0);
214  float mcl_pt = (mcl_itr != cell2multicluster.end() ? mcl_itr->second->pt() : 0.);
215  // Filter cells not included in a multicluster, if requested
216  if (filter_cells_in_multiclusters_ && mcl_id == 0)
217  continue;
218  tc_n_++;
219  // hardware data
220  DetId id(tc_itr->detId());
222  int panel_sector = -999;
223  int panel_number = -999;
224  if (panelId.det() == DetId::Forward) {
225  HGCalDetId panelIdHGCal(panelId);
226  if (panelId.subdetId() == ForwardSubdetector::HGCHEB) {
227  panel_number = panelIdHGCal.wafer();
228  } else {
229  panel_sector = (panelIdHGCal.wafer() >> kSectorOffset_) & kSectorMask_;
230  panel_number = (panelIdHGCal.wafer() >> kPanelOffset_) & kPanelMask_;
231  }
232  } else if (id.det() == DetId::HGCalHSc) {
233  HGCScintillatorDetId panelIdSci(panelId);
234  panel_sector = panelIdSci.iphi();
235  panel_number = panelIdSci.ietaAbs();
236  }
237  tc_panel_number_.emplace_back(panel_number);
238  tc_panel_sector_.emplace_back(panel_sector);
239  tc_id_.emplace_back(tc_itr->detId());
240  tc_side_.emplace_back(triggerTools_.zside(id));
241  tc_layer_.emplace_back(triggerTools_.layerWithOffset(id));
242  // V9 detids
243  if (id.det() == DetId::HGCalTrigger) {
244  HGCalTriggerDetId idv9(id);
245  tc_subdet_.emplace_back(idv9.subdet());
246  tc_waferu_.emplace_back(idv9.waferU());
247  tc_waferv_.emplace_back(idv9.waferV());
248  tc_wafertype_.emplace_back(idv9.type());
249  tc_cellu_.emplace_back(idv9.triggerCellU());
250  tc_cellv_.emplace_back(idv9.triggerCellV());
251  } else if (id.det() == DetId::HGCalHSc) {
252  HGCScintillatorDetId idv9(id);
253  tc_subdet_.emplace_back(idv9.subdet());
254  tc_waferu_.emplace_back(-999);
255  tc_waferv_.emplace_back(-999);
256  tc_wafertype_.emplace_back(idv9.type());
257  tc_cellu_.emplace_back(idv9.ietaAbs());
258  tc_cellv_.emplace_back(idv9.iphi());
259  }
260  // V8 detids
261  else {
262  HGCalDetId idv8(id);
263  tc_subdet_.emplace_back(id.subdetId());
264  tc_wafer_.emplace_back(idv8.wafer());
265  tc_wafertype_.emplace_back(idv8.waferType());
266  tc_cell_.emplace_back(idv8.cell());
267  }
268  tc_data_.emplace_back(tc_itr->hwPt());
269  tc_uncompressedCharge_.emplace_back(tc_itr->uncompressedCharge());
270  tc_compressedCharge_.emplace_back(tc_itr->compressedCharge());
271  tc_mipPt_.emplace_back(tc_itr->mipPt());
272  // physical values
273  tc_pt_.emplace_back(tc_itr->pt());
274  tc_energy_.emplace_back(tc_itr->energy());
275  tc_eta_.emplace_back(tc_itr->eta());
276  tc_phi_.emplace_back(tc_itr->phi());
277  tc_x_.emplace_back(tc_itr->position().x());
278  tc_y_.emplace_back(tc_itr->position().y());
279  tc_z_.emplace_back(tc_itr->position().z());
280  // Links between TC and clusters
281  tc_cluster_id_.emplace_back(cl_id);
282  tc_multicluster_id_.emplace_back(mcl_id);
283  tc_multicluster_pt_.emplace_back(mcl_pt);
284 
285  if (fill_simenergy_) {
286  double energy = 0;
287  unsigned layer = triggerTools_.layerWithOffset(id);
288  // search for simhit for all the cells inside the trigger cell
289  for (uint32_t c_id : geometry_->getCellsFromTriggerCell(id)) {
291  if (triggerTools_.isEm(id)) {
292  auto itr = simhits_ee.find(c_id);
293  if (itr != simhits_ee.end())
294  energy += calibrate(itr->second, thickness, layer);
295  } else if (triggerTools_.isSilicon(id)) {
296  auto itr = simhits_fh.find(c_id);
297  if (itr != simhits_fh.end())
298  energy += calibrate(itr->second, thickness, layer);
299  } else {
300  auto itr = simhits_bh.find(c_id);
301  if (itr != simhits_bh.end())
302  energy += itr->second;
303  }
304  }
305  tc_simenergy_.emplace_back(energy);
306  }
307  }
308 
309  if (fill_truthmap_) {
310  auto itr(cell_to_genparticle.find(tc_itr->detId()));
311  if (itr == cell_to_genparticle.end())
312  tc_genparticle_index_.push_back(-1);
313  else
314  tc_genparticle_index_.push_back(itr->second);
315  }
316  }
317 }
318 
319 double HGCalTriggerNtupleHGCTriggerCells::calibrate(double energy, int thickness, unsigned layer) {
320  double fcPerMip = fcPerMip_[thickness];
322  double layerWeight = layerWeights_[layer];
323  double TeV2GeV = 1.e3;
324  return energy * keV2fC_ / fcPerMip * layerWeight * TeV2GeV / thicknessCorrection;
325 }
326 
328  std::unordered_map<uint32_t, double>& simhits_ee,
329  std::unordered_map<uint32_t, double>& simhits_fh,
330  std::unordered_map<uint32_t, double>& simhits_bh) {
332  e.getByToken(simhits_ee_token_, ee_simhits_h);
333  const edm::PCaloHitContainer& ee_simhits = *ee_simhits_h;
335  e.getByToken(simhits_fh_token_, fh_simhits_h);
336  const edm::PCaloHitContainer& fh_simhits = *fh_simhits_h;
338  e.getByToken(simhits_bh_token_, bh_simhits_h);
339  const edm::PCaloHitContainer& bh_simhits = *bh_simhits_h;
340 
341  //EE
342  for (const auto& simhit : ee_simhits) {
343  DetId id = triggerTools_.simToReco(simhit.id(), geometry_->eeTopology());
344  if (id.rawId() == 0)
345  continue;
346  auto itr_insert = simhits_ee.emplace(id, 0.);
347  itr_insert.first->second += simhit.energy();
348  }
349  // FH
350  for (const auto& simhit : fh_simhits) {
351  DetId id = triggerTools_.simToReco(simhit.id(), geometry_->fhTopology());
352  if (id.rawId() == 0)
353  continue;
354  auto itr_insert = simhits_fh.emplace(id, 0.);
355  itr_insert.first->second += simhit.energy();
356  }
357  // BH
358  for (const auto& simhit : bh_simhits) {
360  : triggerTools_.simToReco(simhit.id(), geometry_->bhTopology()));
361  if (id.rawId() == 0)
362  continue;
363  auto itr_insert = simhits_bh.emplace(id, 0.);
364  itr_insert.first->second += simhit.energy();
365  }
366 }
367 
369  tc_n_ = 0;
370  tc_id_.clear();
371  tc_subdet_.clear();
372  tc_side_.clear();
373  tc_layer_.clear();
374  tc_wafer_.clear();
375  tc_waferu_.clear();
376  tc_waferv_.clear();
377  tc_wafertype_.clear();
378  tc_panel_number_.clear();
379  tc_panel_sector_.clear();
380  tc_cell_.clear();
381  tc_cellu_.clear();
382  tc_cellv_.clear();
383  tc_data_.clear();
384  tc_uncompressedCharge_.clear();
385  tc_compressedCharge_.clear();
386  tc_mipPt_.clear();
387  tc_pt_.clear();
388  tc_energy_.clear();
389  tc_simenergy_.clear();
390  tc_eta_.clear();
391  tc_phi_.clear();
392  tc_x_.clear();
393  tc_y_.clear();
394  tc_z_.clear();
395  tc_cluster_id_.clear();
396  tc_multicluster_id_.clear();
397  tc_multicluster_pt_.clear();
398  tc_genparticle_index_.clear();
399 }
void initialize(TTree &, const edm::ParameterSet &, edm::ConsumesCollector &&) final
const_iterator end(int bx) const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::ESHandle< HGCalTriggerGeometryBase > geometry_
std::vector< PCaloHit > PCaloHitContainer
virtual geom_set getCellsFromTriggerCell(const unsigned cell_det_id) const =0
void eventSetup(const edm::EventSetup &)
HGCalTriggerSubdetector subdet() const
get the subdetector
void fill(const edm::Event &e, const edm::EventSetup &es) final
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const HGCalTopology & eeTopology() const
int triggerCellU() const
get the cell #&#39;s in u,v or in x,y
HGCalTriggerNtupleHGCTriggerCells(const edm::ParameterSet &conf)
int type() const
get the type
unsigned layerWithOffset(const DetId &) const
DetId::Detector subdet() const
get the subdetector
int type() const
get the type
int zside(const DetId &) const
const HcalTopology & bhTopology() const
int iphi() const
get the phi index
int wafer() const
get the wafer #
Definition: HGCalDetId.h:40
int triggerCellV() const
Definition: DetId.h:17
int thicknessIndex(const DetId &, bool tc=false) const
bool isSilicon(const DetId &) const
int cell() const
get the absolute value of the cell #&#39;s in x and y
Definition: HGCalDetId.h:37
virtual unsigned getModuleFromTriggerCell(const unsigned trigger_cell_det_id) const =0
const HGCalTopology & hscTopology() const
bool isEm(const DetId &) const
T get() const
Definition: EventSetup.h:73
int waferType() const
get the wafer type
Definition: HGCalDetId.h:43
#define DEFINE_EDM_PLUGIN(factory, type, name)
DetId simToReco(const DetId &, const HGCalTopology &) const
Definition: tree.py:1
void simhits(const edm::Event &e, std::unordered_map< uint32_t, double > &simhits_ee, std::unordered_map< uint32_t, double > &simhits_fh, std::unordered_map< uint32_t, double > &simhits_bh)
const_iterator begin(int bx) const
edm::AssociationMap< edm::OneToMany< CaloParticleCollection, l1t::HGCalTriggerCellBxCollection > > CaloToCellsMap
const HGCalTopology & fhTopology() const