CMS 3D CMS Logo

HGCalTriggerNtupleHGCDigis.cc
Go to the documentation of this file.
8 
11 
17 
18 
20 {
21 
22  public:
25  void initialize(TTree&, const edm::ParameterSet&, edm::ConsumesCollector&&) final;
26  void fill(const edm::Event& e, const edm::EventSetup& es) final;
27 
28  private:
29  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);
30  void clear() final;
31 
35 
37 
38  int hgcdigi_n_ ;
39  std::vector<int> hgcdigi_id_;
40  std::vector<int> hgcdigi_subdet_;
41  std::vector<int> hgcdigi_side_;
42  std::vector<int> hgcdigi_layer_;
43  std::vector<int> hgcdigi_wafer_;
44  std::vector<int> hgcdigi_wafertype_ ;
45  std::vector<int> hgcdigi_cell_;
46  std::vector<float> hgcdigi_eta_;
47  std::vector<float> hgcdigi_phi_;
48  std::vector<float> hgcdigi_z_;
49  std::vector<uint32_t> hgcdigi_data_;
50  std::vector<int> hgcdigi_isadc_;
51  std::vector<float> hgcdigi_simenergy_;
52 
53  int bhdigi_n_ ;
54  std::vector<int> bhdigi_id_;
55  std::vector<int> bhdigi_subdet_;
56  std::vector<int> bhdigi_side_;
57  std::vector<int> bhdigi_layer_;
58  std::vector<int> bhdigi_ieta_;
59  std::vector<int> bhdigi_iphi_;
60  std::vector<float> bhdigi_eta_;
61  std::vector<float> bhdigi_phi_;
62  std::vector<float> bhdigi_z_;
63  std::vector<uint32_t> bhdigi_data_;
64  std::vector<float> bhdigi_simenergy_;
65 
67 
68 };
69 
72  "HGCalTriggerNtupleHGCDigis" );
73 
74 
77 {
78  is_Simhit_comp_ = conf.getParameter<bool>("isSimhitComp");
79 
80 }
81 
82 void
84 initialize(TTree& tree, const edm::ParameterSet& conf, edm::ConsumesCollector&& collector)
85 {
86 
87  ee_token_ = collector.consumes<HGCEEDigiCollection>(conf.getParameter<edm::InputTag>("HGCDigisEE"));
88  fh_token_ = collector.consumes<HGCHEDigiCollection>(conf.getParameter<edm::InputTag>("HGCDigisFH"));
89  bh_token_ = collector.consumes<HGCBHDigiCollection>(conf.getParameter<edm::InputTag>("HGCDigisBH"));
90  if (is_Simhit_comp_) {
91  SimHits_inputee_ = collector.consumes<edm::PCaloHitContainer>(conf.getParameter<edm::InputTag>("eeSimHits"));
92  SimHits_inputfh_ = collector.consumes<edm::PCaloHitContainer>(conf.getParameter<edm::InputTag>("fhSimHits"));
93  SimHits_inputbh_ = collector.consumes<edm::PCaloHitContainer>(conf.getParameter<edm::InputTag>("bhSimHits"));
94  }
95  tree.Branch("hgcdigi_n", &hgcdigi_n_, "hgcdigi_n/I");
96  tree.Branch("hgcdigi_id", &hgcdigi_id_);
97  tree.Branch("hgcdigi_subdet", &hgcdigi_subdet_);
98  tree.Branch("hgcdigi_zside", &hgcdigi_side_);
99  tree.Branch("hgcdigi_layer", &hgcdigi_layer_);
100  tree.Branch("hgcdigi_wafer", &hgcdigi_wafer_);
101  tree.Branch("hgcdigi_wafertype", &hgcdigi_wafertype_);
102  tree.Branch("hgcdigi_cell", &hgcdigi_cell_);
103  tree.Branch("hgcdigi_eta", &hgcdigi_eta_);
104  tree.Branch("hgcdigi_phi", &hgcdigi_phi_);
105  tree.Branch("hgcdigi_z", &hgcdigi_z_);
106  tree.Branch("hgcdigi_data", &hgcdigi_data_);
107  tree.Branch("hgcdigi_isadc", &hgcdigi_isadc_);
108  if (is_Simhit_comp_) tree.Branch("hgcdigi_simenergy", &hgcdigi_simenergy_);
109 
110  tree.Branch("bhdigi_n", &bhdigi_n_, "bhdigi_n/I");
111  tree.Branch("bhdigi_id", &bhdigi_id_);
112  tree.Branch("bhdigi_subdet", &bhdigi_subdet_);
113  tree.Branch("bhdigi_zside", &bhdigi_side_);
114  tree.Branch("bhdigi_layer", &bhdigi_layer_);
115  tree.Branch("bhdigi_ieta", &bhdigi_ieta_);
116  tree.Branch("bhdigi_iphi", &bhdigi_iphi_);
117  tree.Branch("bhdigi_eta", &bhdigi_eta_);
118  tree.Branch("bhdigi_phi", &bhdigi_phi_);
119  tree.Branch("bhdigi_z", &bhdigi_z_);
120  tree.Branch("bhdigi_data", &bhdigi_data_);
121  if (is_Simhit_comp_) tree.Branch("bhdigi_simenergy", &bhdigi_simenergy_);
122 }
123 
124 void
126 fill(const edm::Event& e, const edm::EventSetup& es)
127 {
129 
131  e.getByToken(ee_token_, ee_digis_h);
132  const HGCEEDigiCollection& ee_digis = *ee_digis_h;
134  e.getByToken(fh_token_, fh_digis_h);
135  const HGCHEDigiCollection& fh_digis = *fh_digis_h;
137  e.getByToken(bh_token_, bh_digis_h);
138  const HGCBHDigiCollection& bh_digis = *bh_digis_h;
139 
141 
142  // sim hit association
143  std::unordered_map<uint32_t, double> simhits_ee;
144  std::unordered_map<uint32_t, double> simhits_fh;
145  std::unordered_map<uint32_t, double> simhits_bh;
146  if (is_Simhit_comp_) simhits(e, simhits_ee, simhits_fh, simhits_bh);
147 
148  clear();
149  hgcdigi_n_ = ee_digis.size() + fh_digis.size();
150  hgcdigi_id_.reserve(hgcdigi_n_);
151  hgcdigi_subdet_.reserve(hgcdigi_n_);
152  hgcdigi_side_.reserve(hgcdigi_n_);
153  hgcdigi_layer_.reserve(hgcdigi_n_);
154  hgcdigi_wafer_.reserve(hgcdigi_n_);
156  hgcdigi_cell_.reserve(hgcdigi_n_);
157  hgcdigi_eta_.reserve(hgcdigi_n_);
158  hgcdigi_phi_.reserve(hgcdigi_n_);
159  hgcdigi_z_.reserve(hgcdigi_n_);
160  hgcdigi_data_.reserve(hgcdigi_n_);
161  hgcdigi_isadc_.reserve(hgcdigi_n_);
163 
164  bhdigi_n_ = bh_digis.size();
165  bhdigi_id_.reserve(bhdigi_n_);
166  bhdigi_subdet_.reserve(bhdigi_n_);
167  bhdigi_side_.reserve(bhdigi_n_);
168  bhdigi_layer_.reserve(bhdigi_n_);
169  bhdigi_ieta_.reserve(bhdigi_n_);
170  bhdigi_iphi_.reserve(bhdigi_n_);
171  bhdigi_eta_.reserve(bhdigi_n_);
172  bhdigi_phi_.reserve(bhdigi_n_);
173  bhdigi_z_.reserve(bhdigi_n_);
175 
176  const int kIntimeSample = 2;
177  for(const auto& digi : ee_digis)
178  {
179  const HGCalDetId id(digi.id());
180  hgcdigi_id_.emplace_back(id.rawId());
182  hgcdigi_side_.emplace_back(id.zside());
183  hgcdigi_layer_.emplace_back(triggerTools_.layerWithOffset(id));
184  hgcdigi_wafer_.emplace_back(id.wafer());
185  hgcdigi_wafertype_.emplace_back(id.waferType());
186  hgcdigi_cell_.emplace_back(id.cell());
187  GlobalPoint cellpos = triggerGeometry_->eeGeometry()->getPosition(id.rawId());
188  hgcdigi_eta_.emplace_back(cellpos.eta());
189  hgcdigi_phi_.emplace_back(cellpos.phi());
190  hgcdigi_z_.emplace_back(cellpos.z());
191  hgcdigi_data_.emplace_back(digi[kIntimeSample].data());
192  int is_adc=0;
193  if (!(digi[kIntimeSample].mode())) is_adc =1;
194  hgcdigi_isadc_.emplace_back(is_adc);
195  if (is_Simhit_comp_) {
196  double hit_energy=0;
197  auto itr = simhits_ee.find(id);
198  if(itr!=simhits_ee.end())hit_energy = itr->second;
199  hgcdigi_simenergy_.emplace_back(hit_energy);
200  }
201  }
202 
203  for(const auto& digi : fh_digis)
204  {
205  const HGCalDetId id(digi.id());
206  hgcdigi_id_.emplace_back(id.rawId());
208  hgcdigi_side_.emplace_back(id.zside());
209  hgcdigi_layer_.emplace_back(triggerTools_.layerWithOffset(id));
210  hgcdigi_wafer_.emplace_back(id.wafer());
211  hgcdigi_wafertype_.emplace_back(id.waferType());
212  hgcdigi_cell_.emplace_back(id.cell());
213  GlobalPoint cellpos = triggerGeometry_->fhGeometry()->getPosition(id.rawId());
214  hgcdigi_eta_.emplace_back(cellpos.eta());
215  hgcdigi_phi_.emplace_back(cellpos.phi());
216  hgcdigi_z_.emplace_back(cellpos.z());
217  hgcdigi_data_.emplace_back(digi[kIntimeSample].data());
218  int is_adc=0;
219  if (!(digi[kIntimeSample].mode())) is_adc =1;
220  hgcdigi_isadc_.emplace_back(is_adc);
221  if (is_Simhit_comp_) {
222  double hit_energy=0;
223  auto itr = simhits_fh.find(id);
224  if(itr!=simhits_fh.end())hit_energy = itr->second;
225  hgcdigi_simenergy_.emplace_back(hit_energy);
226  }
227  }
228 
229  for(const auto& digi : bh_digis)
230  {
231  const HcalDetId id(digi.id());
232  bhdigi_id_.emplace_back(id.rawId());
233  bhdigi_subdet_.emplace_back(id.subdetId());
234  bhdigi_side_.emplace_back(id.zside());
235  bhdigi_layer_.emplace_back(triggerTools_.layerWithOffset(id));
236  bhdigi_ieta_.emplace_back(id.ieta());
237  bhdigi_iphi_.emplace_back(id.iphi());
238  GlobalPoint cellpos = triggerGeometry_->bhGeometry()->getPosition(id.rawId());
239  bhdigi_eta_.emplace_back(cellpos.eta());
240  bhdigi_phi_.emplace_back(cellpos.phi());
241  bhdigi_z_.emplace_back(cellpos.z());
242  bhdigi_data_.emplace_back(digi[kIntimeSample].data());
243  if (is_Simhit_comp_) {
244  double hit_energy=0;
245  auto itr = simhits_bh.find(id);
246  if(itr!=simhits_bh.end())hit_energy = itr->second;
247  bhdigi_simenergy_.emplace_back(hit_energy);
248  }
249  }
250 }
251 
252 void
254 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)
255 {
256 
258  e.getByToken(SimHits_inputee_,ee_simhits_h);
259  const edm::PCaloHitContainer& ee_simhits = *ee_simhits_h;
261  e.getByToken(SimHits_inputfh_,fh_simhits_h);
262  const edm::PCaloHitContainer& fh_simhits = *fh_simhits_h;
264  e.getByToken(SimHits_inputbh_,bh_simhits_h);
265  const edm::PCaloHitContainer& bh_simhits = *bh_simhits_h;
266 
267  //EE
268  int layer=0,cell=0, sec=0, subsec=0, zp=0,subdet=0;
269  ForwardSubdetector mysubdet;
270 
271  for( const auto& simhit : ee_simhits ) {
272  HGCalTestNumbering::unpackHexagonIndex(simhit.id(), subdet, zp, layer, sec, subsec, cell);
273  mysubdet = (ForwardSubdetector)(subdet);
274  std::pair<int,int> recoLayerCell = triggerGeometry_->eeTopology().dddConstants().simToReco(cell,layer,sec,triggerGeometry_->eeTopology().detectorType());
275  cell = recoLayerCell.first;
276  layer = recoLayerCell.second;
277  if (layer<0 || cell<0) {
278  continue;
279  }
280  auto itr_insert = simhits_ee.emplace(HGCalDetId(mysubdet,zp,layer,subsec,sec,cell), 0.);
281  itr_insert.first->second += simhit.energy();
282  }
283 
284  // FH
285  layer=0; cell=0; sec=0; subsec=0; zp=0; subdet=0;
286 
287  for( const auto& simhit : fh_simhits ) {
288  HGCalTestNumbering::unpackHexagonIndex(simhit.id(), subdet, zp, layer, sec, subsec, cell);
289  mysubdet = (ForwardSubdetector)(subdet);
290  std::pair<int,int> recoLayerCell = triggerGeometry_->fhTopology().dddConstants().simToReco(cell,layer,sec,triggerGeometry_->fhTopology().detectorType());
291  cell = recoLayerCell.first;
292  layer = recoLayerCell.second;
293  if (layer<0 || cell<0) {
294  continue;
295  }
296  auto itr_insert = simhits_fh.emplace(HGCalDetId(mysubdet,zp,layer,subsec,sec,cell), 0.);
297  itr_insert.first->second += simhit.energy();
298  }
299  // BH
300  for( const auto& simhit : bh_simhits ) {
302  if (id.subdetId()!=HcalEndcap) continue;
303  auto itr_insert = simhits_bh.emplace(id, 0.);
304  itr_insert.first->second += simhit.energy();
305  }
306 }
307 
308 
309 void
312 {
313  hgcdigi_n_ = 0;
314  hgcdigi_id_.clear();
315  hgcdigi_subdet_.clear();
316  hgcdigi_side_.clear();
317  hgcdigi_layer_.clear();
318  hgcdigi_wafer_.clear();
319  hgcdigi_wafertype_.clear();
320  hgcdigi_cell_.clear();
321  hgcdigi_eta_.clear();
322  hgcdigi_phi_.clear();
323  hgcdigi_z_.clear();
324  hgcdigi_data_.clear();
325  hgcdigi_isadc_.clear();
326  if (is_Simhit_comp_) hgcdigi_simenergy_.clear();
327 
328  bhdigi_n_ = 0;
329  bhdigi_id_.clear();
330  bhdigi_subdet_.clear();
331  bhdigi_side_.clear();
332  bhdigi_layer_.clear();
333  bhdigi_ieta_.clear();
334  bhdigi_iphi_.clear();
335  bhdigi_eta_.clear();
336  bhdigi_phi_.clear();
337  bhdigi_z_.clear();
338  bhdigi_data_.clear();
339  if (is_Simhit_comp_) bhdigi_simenergy_.clear();
340 }
T getParameter(std::string const &) const
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:167
std::vector< PCaloHit > PCaloHitContainer
void eventSetup(const edm::EventSetup &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
const HGCalGeometry * eeGeometry() const
const HGCalTopology & eeTopology() const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
edm::ESHandle< HGCalTriggerGeometryBase > triggerGeometry_
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)
GlobalPoint getPosition(const DetId &id) const
bool detectorType() const
unsigned layerWithOffset(const DetId &) const
int zside(DetId const &)
ForwardSubdetector
void fill(const edm::Event &e, const edm::EventSetup &es) final
std::vector< uint32_t > hgcdigi_data_
void initialize(TTree &, const edm::ParameterSet &, edm::ConsumesCollector &&) final
std::pair< int, int > simToReco(int cell, int layer, int mod, bool half) const
const HGCalGeometry * fhGeometry() const
const HcalTopology & bhTopology() const
T z() const
Definition: PV3DBase.h:64
HGCalTriggerNtupleHGCDigis(const edm::ParameterSet &conf)
GlobalPoint getPosition(const DetId &id) const
const HcalGeometry * bhGeometry() const
std::vector< uint32_t > bhdigi_data_
const HGCalDDDConstants & dddConstants() const
const T & get() const
Definition: EventSetup.h:59
T eta() const
Definition: PV3DBase.h:76
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
size_type size() const
#define DEFINE_EDM_PLUGIN(factory, type, name)
DetId relabel(const uint32_t testId) const
Definition: tree.py:1
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
const HGCalTopology & fhTopology() const