CMS 3D CMS Logo

HGCDigitizer.cc
Go to the documentation of this file.
19 #include <algorithm>
20 #include <boost/foreach.hpp>
22 
23 //#define EDM_ML_DEBUG
24 using namespace std;
25 using namespace hgc_digi;
26 
27 typedef std::unordered_map<uint32_t, std::vector<std::pair<float, float>>> IdHit_Map;
28 typedef std::tuple<float, float, float> hit_timeStamp;
29 typedef std::unordered_map<uint32_t, std::vector<hit_timeStamp>> hitRec_container;
30 
31 namespace {
32 
33  constexpr std::array<double, 4> occupancyGuesses = {{0.5, 0.2, 0.2, 0.8}};
34 
35  int getCellThickness(const HGCalGeometry* geom, const DetId& detid) {
36  const auto& dddConst = geom->topology().dddConstants();
37  return (1 + dddConst.waferType(detid));
38  }
39 
40  void getValidDetIds(const HGCalGeometry* geom, std::unordered_set<DetId>& valid) {
41  const std::vector<DetId>& ids = geom->getValidDetIds();
42  valid.reserve(ids.size());
43  valid.insert(ids.begin(), ids.end());
44  }
45 
46  DetId simToReco(const HGCalGeometry* geom, unsigned simId) {
47  DetId result(0);
48  const auto& topo = geom->topology();
49  const auto& dddConst = topo.dddConstants();
50 
51  if (dddConst.waferHexagon8() || dddConst.tileTrapezoid()) {
52  result = DetId(simId);
53  } else {
54  int subdet(DetId(simId).subdetId()), layer, cell, sec, subsec, zp;
55  HGCalTestNumbering::unpackHexagonIndex(simId, subdet, zp, layer, sec, subsec, cell);
56  //sec is wafer and subsec is celltyp
57  //skip this hit if after ganging it is not valid
58  auto recoLayerCell = dddConst.simToReco(cell, layer, sec, topo.detectorType());
59  cell = recoLayerCell.first;
60  layer = recoLayerCell.second;
61  if (layer < 0 || cell < 0) {
62  return result;
63  } else {
64  //assign the RECO DetId
65  result = HGCalDetId((ForwardSubdetector)subdet, zp, layer, subsec, sec, cell);
66  }
67  }
68  return result;
69  }
70 
71  void saveSimHitAccumulator_forPreMix(PHGCSimAccumulator& simResult,
72  const hgc::HGCPUSimHitDataAccumulator& simData,
73  const std::unordered_set<DetId>& validIds,
74  const float minCharge,
75  const float maxCharge) {
76  constexpr auto nEnergies = std::tuple_size<decltype(hgc_digi::HGCCellHitInfo().PUhit_info)>::value;
77  static_assert(nEnergies <= PHGCSimAccumulator::SimHitCollection::energyMask + 1,
78  "PHGCSimAccumulator bit pattern needs to updated");
80  "PHGCSimAccumulator bit pattern needs to updated");
81  const float minPackChargeLog = minCharge > 0.f ? std::log(minCharge) : -2;
82  const float maxPackChargeLog = std::log(maxCharge);
84 
85  simResult.reserve(simData.size());
86 
87  for (const auto& id : validIds) {
88  auto found = simData.find(id);
89  if (found == simData.end())
90  continue;
91 
92  const hgc_digi::PUSimHitData& accCharge_across_bxs = found->second.PUhit_info[0];
93  const hgc_digi::PUSimHitData& timing_across_bxs = found->second.PUhit_info[1];
94  for (size_t iSample = 0; iSample < hgc_digi::nSamples; ++iSample) {
95  const std::vector<float>& accCharge_inthis_bx = accCharge_across_bxs[iSample];
96  const std::vector<float>& timing_inthis_bx = timing_across_bxs[iSample];
97  std::vector<unsigned short> vc, vt;
98  size_t nhits = accCharge_inthis_bx.size();
99 
100  for (size_t ihit = 0; ihit < nhits; ++ihit) {
101  if (accCharge_inthis_bx[ihit] > minCharge) {
102  unsigned short c =
103  logintpack::pack16log(accCharge_inthis_bx[ihit], minPackChargeLog, maxPackChargeLog, base);
104  unsigned short t = logintpack::pack16log(timing_inthis_bx[ihit], minPackChargeLog, maxPackChargeLog, base);
105  vc.emplace_back(c);
106  vt.emplace_back(t);
107  }
108  }
109  simResult.emplace_back(id.rawId(), iSample, vc, vt);
110  }
111  }
112  simResult.shrink_to_fit();
113  }
114 
115  void loadSimHitAccumulator_forPreMix(hgc::HGCSimHitDataAccumulator& simData,
117  const HGCalGeometry* geom,
118  IdHit_Map& hitRefs_bx0,
119  const PHGCSimAccumulator& simAccumulator,
120  const float minCharge,
121  const float maxCharge,
122  bool setIfZero,
123  const std::array<float, 3>& tdcForToAOnset,
124  const bool minbiasFlag,
125  std::unordered_map<uint32_t, bool>& hitOrder_monitor,
126  const unsigned int thisBx) {
127  const float minPackChargeLog = minCharge > 0.f ? std::log(minCharge) : -2;
128  const float maxPackChargeLog = std::log(maxCharge);
130  for (const auto& detIdIndexHitInfo : simAccumulator) {
131  unsigned int detId = detIdIndexHitInfo.detId();
132 
133  auto simIt = simData.emplace(detId, HGCCellInfo()).first;
134  size_t nhits = detIdIndexHitInfo.nhits();
135 
136  hitOrder_monitor[detId] = false;
137  if (nhits > 0) {
138  unsigned short iSample = detIdIndexHitInfo.sampleIndex();
139 
140  const auto& unsigned_charge_array = detIdIndexHitInfo.chargeArray();
141  const auto& unsigned_time_array = detIdIndexHitInfo.timeArray();
142 
143  float p_charge, p_time;
144  unsigned short unsigned_charge, unsigned_time;
145 
146  for (size_t ihit = 0; ihit < nhits; ++ihit) {
147  unsigned_charge = (unsigned_charge_array[ihit] & PHGCSimAccumulator::SimHitCollection::dataMask);
148  unsigned_time = (unsigned_time_array[ihit] & PHGCSimAccumulator::SimHitCollection::dataMask);
149  p_time = logintpack::unpack16log(unsigned_time, minPackChargeLog, maxPackChargeLog, base);
150  p_charge = logintpack::unpack16log(unsigned_charge, minPackChargeLog, maxPackChargeLog, base);
151 
152  (simIt->second).hit_info[0][iSample] += p_charge;
153  if (iSample == (unsigned short)thisBx) {
154  if (hitRefs_bx0[detId].empty()) {
155  hitRefs_bx0[detId].emplace_back(p_charge, p_time);
156  } else {
157  if (p_time < hitRefs_bx0[detId].back().second) {
158  auto findPos = std::upper_bound(hitRefs_bx0[detId].begin(),
159  hitRefs_bx0[detId].end(),
160  std::make_pair(0.f, p_time),
161  [](const auto& i, const auto& j) { return i.second < j.second; });
162 
163  auto insertedPos = findPos;
164  if (findPos == hitRefs_bx0[detId].begin()) {
165  insertedPos = hitRefs_bx0[detId].emplace(insertedPos, p_charge, p_time);
166  } else {
167  auto prevPos = findPos - 1;
168  if (prevPos->second == p_time) {
169  prevPos->first = prevPos->first + p_charge;
170  insertedPos = prevPos;
171  } else if (prevPos->second < p_time) {
172  insertedPos = hitRefs_bx0[detId].emplace(findPos, (prevPos)->first + p_charge, p_time);
173  }
174  }
175 
176  for (auto step = insertedPos; step != hitRefs_bx0[detId].end(); ++step) {
177  if (step != insertedPos)
178  step->first += p_charge;
179  }
180 
181  hitOrder_monitor[detId] = true;
182 
183  } else if (p_time == hitRefs_bx0[detId].back().second) {
184  hitRefs_bx0[detId].back().first += p_charge;
185  } else if (p_time > hitRefs_bx0[detId].back().second) {
186  hitRefs_bx0[detId].emplace_back(hitRefs_bx0[detId].back().first + p_charge, p_time);
187  }
188  }
189  }
190  }
191  }
192  }
193 
194  if (minbiasFlag) {
195  for (const auto& hitmapIterator : hitRefs_bx0) {
196  unsigned int detectorId = hitmapIterator.first;
197  auto simIt = simData.emplace(detectorId, HGCCellInfo()).first;
198  const bool orderChanged = hitOrder_monitor[detectorId];
199  int waferThickness = getCellThickness(geom, detectorId);
200  float cell_threshold = tdcForToAOnset[waferThickness - 1];
201  const auto& hitRec = hitmapIterator.second;
202  float accChargeForToA(0.f), fireTDC(0.f);
203  const auto aboveThrPos = std::upper_bound(
204  hitRec.begin(), hitRec.end(), std::make_pair(cell_threshold, 0.f), [](const auto& i, const auto& j) {
205  return i.first < j.first;
206  });
207 
208  if (aboveThrPos == hitRec.end()) {
209  accChargeForToA = hitRec.back().first;
210  fireTDC = 0.f;
211  } else if (hitRec.end() - aboveThrPos > 0 || orderChanged) {
212  accChargeForToA = aboveThrPos->first;
213  fireTDC = aboveThrPos->second;
214  if (aboveThrPos - hitRec.begin() >= 1) {
215  const auto& belowThrPos = aboveThrPos - 1;
216  float chargeBeforeThr = belowThrPos->first;
217  float timeBeforeThr = belowThrPos->second;
218  float deltaQ = accChargeForToA - chargeBeforeThr;
219  float deltaTOF = fireTDC - timeBeforeThr;
220  fireTDC = (cell_threshold - chargeBeforeThr) * deltaTOF / deltaQ + timeBeforeThr;
221  }
222  }
223  (simIt->second).hit_info[1][9] = fireTDC;
224  }
225  }
226  }
227 } //namespace
228 
230  : simHitAccumulator_(new HGCSimHitDataAccumulator()),
231  pusimHitAccumulator_(new HGCPUSimHitDataAccumulator()),
232  digiCollection_(ps.getParameter<std::string>("digiCollection")),
233  digitizationType_(ps.getParameter<uint32_t>("digitizationType")),
234  premixStage1_(ps.getParameter<bool>("premixStage1")),
235  premixStage1MinCharge_(ps.getParameter<double>("premixStage1MinCharge")),
236  premixStage1MaxCharge_(ps.getParameter<double>("premixStage1MaxCharge")),
237  maxSimHitsAccTime_(ps.getParameter<uint32_t>("maxSimHitsAccTime")),
238  bxTime_(ps.getParameter<double>("bxTime")),
239  hitsProducer_(ps.getParameter<std::string>("hitsProducer")),
240  hitCollection_(ps.getParameter<std::string>("hitCollection")),
241  hitToken_(iC.consumes<std::vector<PCaloHit>>(edm::InputTag(hitsProducer_, hitCollection_))),
242  geomToken_(iC.esConsumes()),
243  verbosity_(ps.getUntrackedParameter<uint32_t>("verbosity", 0)),
244  tofDelay_(ps.getParameter<double>("tofDelay")),
245  averageOccupancies_(occupancyGuesses),
246  nEvents_(1) {
247  //configure from cfg
248 
249  const auto& myCfg_ = ps.getParameter<edm::ParameterSet>("digiCfg");
250 
251  if (myCfg_.existsAs<edm::ParameterSet>("chargeCollectionEfficiencies")) {
252  cce_.clear();
253  const auto& temp = myCfg_.getParameter<edm::ParameterSet>("chargeCollectionEfficiencies")
254  .getParameter<std::vector<double>>("values");
255  for (double cce : temp) {
256  cce_.emplace_back(cce);
257  }
258  } else {
259  std::vector<float>().swap(cce_);
260  }
261 
262  auto const& pluginName = ps.getParameter<std::string>("digitizer");
264 }
265 
266 //
268  if (geomWatcher_.check(es)) {
269  std::unordered_set<DetId>().swap(validIds_);
270 
271  //get geometry
272  CaloGeometry const& geom = es.getData(geomToken_);
273 
274  gHGCal_ =
275  dynamic_cast<const HGCalGeometry*>(geom.getSubdetectorGeometry(theDigitizer_->det(), theDigitizer_->subdet()));
276 
277  int nadded(0);
278  //valid ID lists
279  if (nullptr != gHGCal_) {
280  getValidDetIds(gHGCal_, validIds_);
281  } else {
282  throw cms::Exception("BadConfiguration") << "HGCDigitizer is not producing EE, FH, or BH digis!";
283  }
284 
285  if (verbosity_ > 0)
286  edm::LogInfo("HGCDigitizer") << "Added " << nadded << ":" << validIds_.size() << " detIds without "
287  << hitCollection_ << " in first event processed" << std::endl;
288  }
289 
290  // reserve memory for a full detector
291  unsigned idx = getType();
294 }
295 
296 //
297 void HGCDigitizer::finalizeEvent(edm::Event& e, edm::EventSetup const& es, CLHEP::HepRandomEngine* hre) {
298  hitRefs_bx0.clear();
299  PhitRefs_bx0.clear();
300  hitOrder_monitor.clear();
301 
302  const CaloSubdetectorGeometry* theGeom = static_cast<const CaloSubdetectorGeometry*>(gHGCal_);
303 
304  ++nEvents_;
305 
306  unsigned idx = getType();
307  // release memory for unfilled parts of hash table
308  if (validIds_.size() * averageOccupancies_[idx] > simHitAccumulator_->size()) {
309  simHitAccumulator_->reserve(simHitAccumulator_->size());
310  pusimHitAccumulator_->reserve(simHitAccumulator_->size());
311  }
312  //update occupancy guess
313  const double thisOcc = simHitAccumulator_->size() / ((double)validIds_.size());
315 
316  if (premixStage1_) {
317  auto simRecord = std::make_unique<PHGCSimAccumulator>();
318 
319  if (!pusimHitAccumulator_->empty()) {
320  saveSimHitAccumulator_forPreMix(
322  }
323 
324  e.put(std::move(simRecord), digiCollection());
325 
326  } else {
327  auto digiResult = std::make_unique<HGCalDigiCollection>();
328  theDigitizer_->run(digiResult, *simHitAccumulator_, theGeom, validIds_, digitizationType_, hre);
329  edm::LogVerbatim("HGCDigitizer") << "HGCDigitizer:: finalize event - produced " << digiResult->size()
330  << " hits in det/subdet " << theDigitizer_->det() << "/"
331  << theDigitizer_->subdet();
332 #ifdef EDM_ML_DEBUG
333  checkPosition(&(*digiResult));
334 #endif
335  e.put(std::move(digiResult), digiCollection());
336  }
337 
340 }
343  CLHEP::HepRandomEngine* hre) {
344  //get inputs
345 
347  if (!hits.isValid()) {
348  edm::LogError("HGCDigitizer") << " @ accumulate_minbias : can't find " << hitCollection_ << " collection of "
349  << hitsProducer_;
350  return;
351  }
352 
353  //accumulate in-time the main event
354  if (nullptr != gHGCal_) {
356  } else {
357  throw cms::Exception("BadConfiguration") << "HGCDigitizer is not producing EE, FH, or BH digis!";
358  }
359 }
360 
361 //
362 void HGCDigitizer::accumulate(edm::Event const& e, edm::EventSetup const& eventSetup, CLHEP::HepRandomEngine* hre) {
363  //get inputs
365  if (!hits.isValid()) {
366  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of "
367  << hitsProducer_;
368  return;
369  }
370 
371  //accumulate in-time the main event
372  if (nullptr != gHGCal_) {
373  accumulate(hits, 0, gHGCal_, hre);
374  } else {
375  throw cms::Exception("BadConfiguration") << "HGCDigitizer is not producing EE, FH, or BH digis!";
376  }
377 }
378 
379 //
382  CLHEP::HepRandomEngine* hre) {
385  e.getByLabel(hitTag, hits);
386 
387  if (!hits.isValid()) {
388  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of "
389  << hitsProducer_;
390  return;
391  }
392 
393  if (nullptr != gHGCal_) {
394  accumulate_forPreMix(hits, e.bunchCrossing(), gHGCal_, hre);
395  } else {
396  throw cms::Exception("BadConfiguration") << "HGCDigitizer is not producing EE, FH, or BH digis!";
397  }
398 }
399 
400 //
403  CLHEP::HepRandomEngine* hre) {
404  //get inputs
407  e.getByLabel(hitTag, hits);
408 
409  if (!hits.isValid()) {
410  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of "
411  << hitsProducer_;
412  return;
413  }
414 
415  //accumulate for the simulated bunch crossing
416  if (nullptr != gHGCal_) {
417  accumulate(hits, e.bunchCrossing(), gHGCal_, hre);
418  } else {
419  throw cms::Exception("BadConfiguration") << "HGCDigitizer is not producing EE, FH, or BH digis!";
420  }
421 }
422 
423 //
425  int bxCrossing,
426  const HGCalGeometry* geom,
427  CLHEP::HepRandomEngine* hre) {
428  if (nullptr == geom)
429  return;
430 
431  auto keV2fC = theDigitizer_->keV2fC();
432  auto tdcForToAOnset = theDigitizer_->tdcForToAOnset();
433 
434  int nchits = (int)hits->size();
435  int count_thisbx = 0;
436  std::vector<HGCCaloHitTuple_t> hitRefs;
437  hitRefs.reserve(nchits);
438  for (int i = 0; i < nchits; ++i) {
439  const auto& the_hit = hits->at(i);
440  DetId id = simToReco(geom, the_hit.id());
441  // to be written the verbosity block
442  if (id.rawId() != 0) {
443  hitRefs.emplace_back(i, id.rawId(), (float)the_hit.time());
444  }
445  }
446  std::sort(hitRefs.begin(), hitRefs.end(), this->orderByDetIdThenTime);
447 
448  nchits = hitRefs.size();
449  for (int i = 0; i < nchits; ++i) {
450  const int hitidx = std::get<0>(hitRefs[i]);
451  const uint32_t id = std::get<1>(hitRefs[i]);
452  if (!validIds_.count(id))
453  continue;
454 
455  if (id == 0)
456  continue;
457 
458  const float toa = std::get<2>(hitRefs[i]);
459  const PCaloHit& hit = hits->at(hitidx);
460  const float charge = hit.energy() * 1e6 * keV2fC; // * getCCE(geom, id, cce_);
461 
462  const float tof = toa + tofDelay_;
463  const int itime = std::floor(tof / bxTime_) + 9;
464 
465  if (itime < 0 || itime > (int)maxBx_)
466  continue;
467 
468  if (itime >= (int)(maxBx_ + 1))
469  continue;
470 
471  int waferThickness = getCellThickness(geom, id);
472  if (itime == (int)thisBx_) {
473  ++count_thisbx;
474  if (PhitRefs_bx0[id].empty()) {
475  PhitRefs_bx0[id].emplace_back(charge, charge, tof);
476  } else if (tof > std::get<2>(PhitRefs_bx0[id].back())) {
477  PhitRefs_bx0[id].emplace_back(charge, charge + std::get<1>(PhitRefs_bx0[id].back()), tof);
478  } else if (tof == std::get<2>(PhitRefs_bx0[id].back())) {
479  std::get<0>(PhitRefs_bx0[id].back()) += charge;
480  std::get<1>(PhitRefs_bx0[id].back()) += charge;
481  } else {
482  //find position to insert new entry preserving time sorting
483  auto findPos = std::upper_bound(PhitRefs_bx0[id].begin(),
484  PhitRefs_bx0[id].end(),
485  hit_timeStamp(charge, 0.f, tof),
486  [](const auto& i, const auto& j) { return std::get<2>(i) <= std::get<2>(j); });
487 
488  auto insertedPos = findPos;
489 
490  if (tof == std::get<2>(*(findPos - 1))) {
491  std::get<0>(*(findPos - 1)) += charge;
492  std::get<1>(*(findPos - 1)) += charge;
493 
494  } else {
495  insertedPos = PhitRefs_bx0[id].insert(findPos,
496  (findPos == PhitRefs_bx0[id].begin())
497  ? hit_timeStamp(charge, charge, tof)
498  : hit_timeStamp(charge, charge + std::get<1>(*(findPos - 1)), tof));
499  }
500  //cumulate the charge of new entry for all elements that follow in the sorted list
501  //and resize list accounting for cases when the inserted element itself crosses the threshold
502 
503  for (auto step = insertedPos; step != PhitRefs_bx0[id].end(); ++step) {
504  if (step != insertedPos)
505  std::get<1>(*(step)) += charge;
506 
507  // resize the list stopping with the first timeStamp with cumulative charge above threshold
508  if (std::get<1>(*step) > tdcForToAOnset[waferThickness - 1] &&
509  std::get<2>(*step) != std::get<2>(PhitRefs_bx0[id].back())) {
510  PhitRefs_bx0[id].resize(
511  std::upper_bound(PhitRefs_bx0[id].begin(),
512  PhitRefs_bx0[id].end(),
513  hit_timeStamp(charge, 0.f, std::get<2>(*step)),
514  [](const auto& i, const auto& j) { return std::get<2>(i) < std::get<2>(j); }) -
515  PhitRefs_bx0[id].begin());
516  for (auto stepEnd = step + 1; stepEnd != PhitRefs_bx0[id].end(); ++stepEnd)
517  std::get<1>(*stepEnd) += charge;
518  break;
519  }
520  }
521  }
522  }
523  }
524 
525  for (const auto& hitCollection : PhitRefs_bx0) {
526  const uint32_t detectorId = hitCollection.first;
527  auto simHitIt = pusimHitAccumulator_->emplace(detectorId, HGCCellHitInfo()).first;
528 
529  for (const auto& hit_timestamp : PhitRefs_bx0[detectorId]) {
530  (simHitIt->second).PUhit_info[1][thisBx_].push_back(std::get<2>(hit_timestamp));
531  (simHitIt->second).PUhit_info[0][thisBx_].push_back(std::get<0>(hit_timestamp));
532  }
533  }
534 
535  if (nchits == 0) {
536  HGCPUSimHitDataAccumulator::iterator simHitIt = pusimHitAccumulator_->emplace(0, HGCCellHitInfo()).first;
537  (simHitIt->second).PUhit_info[1][9].push_back(0.0);
538  (simHitIt->second).PUhit_info[0][9].push_back(0.0);
539  }
540  hitRefs.clear();
541  PhitRefs_bx0.clear();
542 }
543 
544 //
546  int bxCrossing,
547  const HGCalGeometry* geom,
548  CLHEP::HepRandomEngine* hre) {
549  if (nullptr == geom)
550  return;
551 
552  //configuration to apply for the computation of time-of-flight
553  auto weightToAbyEnergy = theDigitizer_->toaModeByEnergy();
554  auto tdcForToAOnset = theDigitizer_->tdcForToAOnset();
555  auto keV2fC = theDigitizer_->keV2fC();
556 
557  //create list of tuples (pos in container, RECO DetId, time) to be sorted first
558  int nchits = (int)hits->size();
559 
560  std::vector<HGCCaloHitTuple_t> hitRefs;
561  hitRefs.reserve(nchits);
562  for (int i = 0; i < nchits; ++i) {
563  const auto& the_hit = hits->at(i);
564  DetId id = simToReco(geom, the_hit.id());
565 
566  if (verbosity_ > 0) {
567  edm::LogVerbatim("HGCDigitizer") << "HGCDigitizer::i/p " << std::hex << the_hit.id() << " o/p " << id.rawId()
568  << std::dec;
569  }
570 
571  if (0 != id.rawId()) {
572  hitRefs.emplace_back(i, id.rawId(), (float)the_hit.time());
573  }
574  }
575 
576  std::sort(hitRefs.begin(), hitRefs.end(), this->orderByDetIdThenTime);
577  //loop over sorted hits
578  nchits = hitRefs.size();
579  for (int i = 0; i < nchits; ++i) {
580  const int hitidx = std::get<0>(hitRefs[i]);
581  const uint32_t id = std::get<1>(hitRefs[i]);
582 
583  //get the data for this cell, if not available then we skip it
584 
585  if (!validIds_.count(id))
586  continue;
587  HGCSimHitDataAccumulator::iterator simHitIt = simHitAccumulator_->emplace(id, HGCCellInfo()).first;
588 
589  if (id == 0)
590  continue; // to be ignored at RECO level
591 
592  const float toa = std::get<2>(hitRefs[i]);
593  const PCaloHit& hit = hits->at(hitidx);
594  const float charge = hit.energy() * 1e6 * keV2fC;
595 
596  //hit time: [time()]=ns + delay
597  //accumulate in 15 buckets of 25ns (9 pre-samples, 1 in-time, 5 post-samples)
598  const float tof = toa + tofDelay_;
599  const int itime = std::floor(tof / bxTime_) + 9;
600 
601  //no need to add bx crossing - tof comes already corrected from the mixing module
602  //itime += bxCrossing;
603  //itime += 9;
604 
605  if (itime < 0 || itime > (int)maxBx_)
606  continue;
607 
608  //check if time index is ok and store energy
609  if (itime >= (int)simHitIt->second.hit_info[0].size())
610  continue;
611 
612  (simHitIt->second).hit_info[0][itime] += charge;
613 
614  //for time-of-arrival: save the time-sorted list of timestamps with cumulative charge just above threshold
615  //working version with pileup only for in-time hits
616  int waferThickness = getCellThickness(geom, id);
617  bool orderChanged = false;
618  if (itime == (int)thisBx_) {
619  //if start empty => just add charge and time
620  if (hitRefs_bx0[id].empty()) {
621  hitRefs_bx0[id].emplace_back(charge, tof);
622 
623  } else if (tof <= hitRefs_bx0[id].back().second) {
624  //find position to insert new entry preserving time sorting
625  std::vector<std::pair<float, float>>::iterator findPos =
626  std::upper_bound(hitRefs_bx0[id].begin(),
627  hitRefs_bx0[id].end(),
628  std::pair<float, float>(0.f, tof),
629  [](const auto& i, const auto& j) { return i.second <= j.second; });
630 
631  std::vector<std::pair<float, float>>::iterator insertedPos = findPos;
632  if (findPos->second == tof) {
633  //just merge timestamps with exact timing
634  findPos->first += charge;
635  } else {
636  //insert new element cumulating the charge
637  insertedPos = hitRefs_bx0[id].insert(findPos,
638  (findPos == hitRefs_bx0[id].begin())
639  ? std::pair<float, float>(charge, tof)
640  : std::pair<float, float>((findPos - 1)->first + charge, tof));
641  }
642 
643  //cumulate the charge of new entry for all elements that follow in the sorted list
644  //and resize list accounting for cases when the inserted element itself crosses the threshold
645  for (std::vector<std::pair<float, float>>::iterator step = insertedPos; step != hitRefs_bx0[id].end(); ++step) {
646  if (step != insertedPos)
647  step->first += charge;
648  // resize the list stopping with the first timeStamp with cumulative charge above threshold
649  if (step->first > tdcForToAOnset[waferThickness - 1] && step->second != hitRefs_bx0[id].back().second) {
650  hitRefs_bx0[id].resize(std::upper_bound(hitRefs_bx0[id].begin(),
651  hitRefs_bx0[id].end(),
652  std::pair<float, float>(0.f, step->second),
653  [](const auto& i, const auto& j) { return i.second < j.second; }) -
654  hitRefs_bx0[id].begin());
655  for (auto stepEnd = step + 1; stepEnd != hitRefs_bx0[id].end(); ++stepEnd)
656  stepEnd->first += charge;
657  break;
658  }
659  }
660 
661  orderChanged = true;
662  } else {
663  //add new entry at the end of the list
664  if (hitRefs_bx0[id].back().first <= tdcForToAOnset[waferThickness - 1]) {
665  hitRefs_bx0[id].emplace_back(hitRefs_bx0[id].back().first + charge, tof);
666  }
667  }
668  }
669  float accChargeForToA = hitRefs_bx0[id].empty() ? 0.f : hitRefs_bx0[id].back().first;
670  //now compute the firing ToA through the interpolation of the consecutive time-stamps at threshold
671  if (weightToAbyEnergy)
672  (simHitIt->second).hit_info[1][itime] += charge * tof;
673  else if (accChargeForToA > tdcForToAOnset[waferThickness - 1] &&
674  ((simHitIt->second).hit_info[1][itime] == 0 || orderChanged == true)) {
675  float fireTDC = hitRefs_bx0[id].back().second;
676  if (hitRefs_bx0[id].size() > 1) {
677  float chargeBeforeThr = (hitRefs_bx0[id].end() - 2)->first;
678  float tofchargeBeforeThr = (hitRefs_bx0[id].end() - 2)->second;
679 
680  float deltaQ = accChargeForToA - chargeBeforeThr;
681  float deltaTOF = fireTDC - tofchargeBeforeThr;
682  fireTDC = (tdcForToAOnset[waferThickness - 1] - chargeBeforeThr) * deltaTOF / deltaQ + tofchargeBeforeThr;
683  }
684  (simHitIt->second).hit_info[1][itime] = fireTDC;
685  }
686  }
687  hitRefs.clear();
688 }
689 void HGCDigitizer::accumulate_forPreMix(const PHGCSimAccumulator& simAccumulator, const bool minbiasFlag) {
690  //configuration to apply for the computation of time-of-flight
691  auto weightToAbyEnergy = theDigitizer_->toaModeByEnergy();
692  auto tdcForToAOnset = theDigitizer_->tdcForToAOnset();
693 
694  if (nullptr != gHGCal_) {
695  loadSimHitAccumulator_forPreMix(*simHitAccumulator_,
697  gHGCal_,
698  hitRefs_bx0,
699  simAccumulator,
702  !weightToAbyEnergy,
703  tdcForToAOnset,
704  minbiasFlag,
706  thisBx_);
707  }
708 }
709 
710 //
712  for (HGCSimHitDataAccumulator::iterator it = simHitAccumulator_->begin(); it != simHitAccumulator_->end(); it++) {
713  it->second.hit_info[0].fill(0.);
714  it->second.hit_info[1].fill(0.);
715  }
716 }
717 
718 uint32_t HGCDigitizer::getType() const {
720  switch (theDigitizer_->det()) {
721  case DetId::HGCalEE:
722  idx = 0;
723  break;
724  case DetId::HGCalHSi:
725  idx = 1;
726  break;
727  case DetId::HGCalHSc:
728  idx = 2;
729  break;
730  case DetId::Forward:
731  idx = 3;
732  break;
733  default:
734  break;
735  }
736  return idx;
737 }
738 
740  const double tol(0.5);
741  if (nullptr != gHGCal_) {
742  for (const auto& digi : *(digis)) {
743  const DetId& id = digi.id();
744  const GlobalPoint& global = gHGCal_->getPosition(id);
745  double r = global.perp();
746  double z = std::abs(global.z());
747  std::pair<double, double> zrange = gHGCal_->topology().dddConstants().rangeZ(true);
748  std::pair<double, double> rrange = gHGCal_->topology().dddConstants().rangeR(z, true);
749  bool ok = ((r >= rrange.first) && (r <= rrange.second) && (z >= zrange.first) && (z <= zrange.second));
750  std::string ck = (((r < rrange.first - tol) || (r > rrange.second + tol) || (z < zrange.first - tol) ||
751  (z > zrange.second + tol))
752  ? "***** ERROR *****"
753  : "");
754  bool val = gHGCal_->topology().valid(id);
755  if ((!ok) || (!val)) {
756  if (id.det() == DetId::HGCalEE || id.det() == DetId::HGCalHSi) {
757  edm::LogVerbatim("HGCDigitizer") << "Check " << HGCSiliconDetId(id) << " " << global << " R " << r << ":"
758  << rrange.first << ":" << rrange.second << " Z " << z << ":" << zrange.first
759  << ":" << zrange.second << " Flag " << ok << ":" << val << " " << ck;
760  } else if (id.det() == DetId::HGCalHSc) {
761  edm::LogVerbatim("HGCDigitizer") << "Check " << HGCScintillatorDetId(id) << " " << global << " R " << r << ":"
762  << rrange.first << ":" << rrange.second << " Z " << z << ":" << zrange.first
763  << ":" << zrange.second << " Flag " << ok << ":" << val << " " << ck;
764  } else if ((id.det() == DetId::Forward) && (id.subdetId() == static_cast<int>(HFNose))) {
765  edm::LogVerbatim("HGCDigitizer") << "Check " << HFNoseDetId(id) << " " << global << " R " << r << ":"
766  << rrange.first << ":" << rrange.second << " Z " << z << ":" << zrange.first
767  << ":" << zrange.second << " Flag " << ok << ":" << val << " " << ck;
768  } else {
769  edm::LogVerbatim("HGCDigitizer")
770  << "Check " << std::hex << id.rawId() << std::dec << " " << id.det() << ":" << id.subdetId() << " "
771  << global << " R " << r << ":" << rrange.first << ":" << rrange.second << " Z " << z << ":"
772  << zrange.first << ":" << zrange.second << " Flag " << ok << ":" << val << " " << ck;
773  }
774  }
775  }
776  }
777 }
const std::string hitsProducer_
Definition: HGCDigitizer.h:101
size
Write out results.
Log< level::Info, true > LogVerbatim
uint32_t nEvents_
Definition: HGCDigitizer.h:126
void emplace_back(unsigned int detId, unsigned short sampleIndex, const std::vector< unsigned short > &accCharge, const std::vector< unsigned short > &timing)
void checkPosition(const HGCalDigiCollection *digis) const
std::pair< double, double > rangeZ(bool reco) const
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< float > cce_
Definition: HGCDigitizer.h:131
T perp() const
Definition: PV3DBase.h:69
void finalizeEvent(edm::Event &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *hre)
std::unique_ptr< HGCDigitizerBase > theDigitizer_
Definition: HGCDigitizer.h:110
void resetSimHitDataAccumulator()
T z() const
Definition: PV3DBase.h:61
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geomToken_
Definition: HGCDigitizer.h:113
bool valid(const DetId &id) const override
Is this a valid cell id.
base
Main Program
Definition: newFWLiteAna.py:92
std::unordered_map< uint32_t, HGCCellInfo > HGCSimHitDataAccumulator
void initializeEvent(edm::Event const &e, edm::EventSetup const &c)
actions at the start/end of event
static constexpr unsigned dataMask
double unpack16log(int16_t i, double lmin, double lmax, uint16_t base=32768)
Definition: liblogintpack.h:62
Log< level::Error, false > LogError
const double premixStage1MaxCharge_
Definition: HGCDigitizer.h:95
ForwardSubdetector
const HGCalGeometry * gHGCal_
Definition: HGCDigitizer.h:116
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
uint32_t getType() const
U second(std::pair< T, U > const &p)
const double premixStage1MinCharge_
Definition: HGCDigitizer.h:93
static const unsigned int maxBx_
Definition: HGCDigitizer.h:129
edm::ESWatcher< CaloGeometryRecord > geomWatcher_
Definition: HGCDigitizer.h:114
std::unordered_map< uint32_t, std::vector< hit_timeStamp > > hitRec_container
Definition: HGCDigitizer.cc:29
const uint32_t verbosity_
Definition: HGCDigitizer.h:119
static constexpr unsigned sampleOffset
std::unordered_map< uint32_t, std::vector< std::tuple< float, float, float > > > PhitRefs_bx0
Definition: HGCDigitizer.h:133
std::unordered_map< uint32_t, bool > hitOrder_monitor
Definition: HGCDigitizer.h:134
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::array< HGCSimDataCollection, nSamples > PUSimHitData
double f[11][100]
std::unordered_set< DetId > validIds_
Definition: HGCDigitizer.h:115
bool getData(T &iHolder) const
Definition: EventSetup.h:122
constexpr size_t nSamples
const HGCalTopology & topology() const
static const unsigned int thisBx_
Definition: HGCDigitizer.h:130
std::unordered_map< uint32_t, std::vector< std::pair< float, float > > > IdHit_Map
Definition: HGCDigitizer.cc:27
std::pair< double, double > rangeR(double z, bool reco) const
Log< level::Info, false > LogInfo
Definition: DetId.h:17
void reserve(size_t size)
const bool premixStage1_
Definition: HGCDigitizer.h:90
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
deadvectors [0] push_back({0.0175431, 0.538005, 6.80997, 13.29})
const edm::EDGetTokenT< std::vector< PCaloHit > > hitToken_
Definition: HGCDigitizer.h:103
GlobalPoint getPosition(const DetId &id, bool debug=false) const
std::array< double, 4 > averageOccupancies_
Definition: HGCDigitizer.h:125
void accumulate(edm::Event const &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *hre)
handle SimHit accumulation
HLT enums.
const float tofDelay_
Definition: HGCDigitizer.h:122
std::unordered_map< uint32_t, std::vector< std::pair< float, float > > > hitRefs_bx0
Definition: HGCDigitizer.h:132
const double bxTime_
Definition: HGCDigitizer.h:99
HGCDigitizer(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
std::tuple< float, float, float > hit_timeStamp
Definition: HGCDigitizer.cc:28
step
Definition: StallMonitor.cc:98
#define get
std::unique_ptr< hgc::HGCSimHitDataAccumulator > simHitAccumulator_
Definition: HGCDigitizer.h:81
static bool orderByDetIdThenTime(const HGCCaloHitTuple_t &a, const HGCCaloHitTuple_t &b)
Definition: HGCDigitizer.h:36
void accumulate_forPreMix(edm::Event const &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *hre)
static constexpr unsigned energyMask
int16_t pack16log(double x, double lmin, double lmax, uint16_t base=32768)
Definition: liblogintpack.h:27
std::unique_ptr< hgc::HGCPUSimHitDataAccumulator > pusimHitAccumulator_
Definition: HGCDigitizer.h:82
std::string digiCollection()
Definition: HGCDigitizer.h:76
const HGCalDDDConstants & dddConstants() const
Definition: HGCalTopology.h:98
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
def move(src, dest)
Definition: eostools.py:511
const int digitizationType_
Definition: HGCDigitizer.h:87
const std::string hitCollection_
Definition: HGCDigitizer.h:102
static constexpr unsigned sampleMask
std::unordered_map< uint32_t, HGCCellHitInfo > HGCPUSimHitDataAccumulator