CMS 3D CMS Logo

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