CMS 3D CMS Logo

HGCDigitizer.cc
Go to the documentation of this file.
19 
20 #include <algorithm>
21 #include <boost/foreach.hpp>
23 
24 using namespace hgc_digi;
25 
26 namespace {
27 
28  constexpr std::array<double,3> occupancyGuesses = { { 0.5,0.2,0.2 } };
29 
30 
31  float getPositionDistance(const HGCalGeometry* geom, const DetId& id) {
32  return geom->getPosition(id).mag();
33  }
34 
35  float getPositionDistance(const HcalGeometry* geom, const DetId& id) {
36  return geom->getGeometry(id)->getPosition().mag();
37  }
38 
39  int getCellThickness(const HGCalGeometry* geom, const DetId& detid ) {
40  const auto& dddConst = geom->topology().dddConstants();
41  return dddConst.waferType(detid);
42  }
43 
44  int getCellThickness(const HcalGeometry* geom, const DetId& detid ) {
45  return 1;
46  }
47 
48  void getValidDetIds(const HGCalGeometry* geom, std::unordered_set<DetId>& valid) {
49  const std::vector<DetId>& ids = geom->getValidDetIds();
50  valid.reserve(ids.size());
51  valid.insert(ids.begin(),ids.end());
52  }
53 
54  void getValidDetIds(const HcalGeometry* geom, std::unordered_set<DetId>& valid) {
55  const std::vector<DetId>& ids = geom->getValidDetIds();
56  for( const auto& id : ids ) {
57  if( HcalEndcap == id.subdetId() &&
58  DetId::Hcal == id.det() )
59  valid.emplace(id);
60  }
61  valid.reserve(valid.size());
62  }
63 
64  DetId simToReco(const HcalGeometry* geom, unsigned simid) {
65  DetId result(0);
66  const auto& topo = geom->topology();
67  const auto* dddConst = topo.dddConstants();
68  HcalDetId id = HcalHitRelabeller::relabel(simid,dddConst);
69 
70  if (id.subdet()==int(HcalEndcap)) {
71  result = id;
72  }
73 
74  return result;
75  }
76 
77  DetId simToReco(const HGCalGeometry* geom, unsigned simId) {
78  DetId result(0);
79  const auto& topo = geom->topology();
80  const auto& dddConst = topo.dddConstants();
81 
82  if ((dddConst.geomMode() == HGCalGeometryMode::Hexagon8) ||
83  (dddConst.geomMode() == HGCalGeometryMode::Hexagon8Full) ||
84  (dddConst.geomMode() == HGCalGeometryMode::Trapezoid)) {
85  result = DetId(simId);
86  } else {
87  int subdet(DetId(simId).subdetId()), layer, cell, sec, subsec, zp;
88  HGCalTestNumbering::unpackHexagonIndex(simId, subdet, zp, layer, sec, subsec, cell);
89  //sec is wafer and subsec is celltyp
90  //skip this hit if after ganging it is not valid
91  auto recoLayerCell = dddConst.simToReco(cell,layer,sec,topo.detectorType());
92  cell = recoLayerCell.first;
93  layer = recoLayerCell.second;
94  if (layer<0 || cell<0) {
95  return result;
96  } else {
97  //assign the RECO DetId
98  result = HGCalDetId((ForwardSubdetector)subdet,zp,layer,subsec,sec,cell);
99  }
100  }
101  return result;
102  }
103 
104  float getCCE(const HGCalGeometry* geom,
105  const DetId& detid,
106  const std::vector<float>&cces) {
107  if ( cces.empty() ) return 1.f;
108  if ((detid.det() == DetId::Hcal) || (detid.det() == DetId::HGCalHSc)) {
109  return 1.f;
110  } else {
111  const auto& topo = geom->topology();
112  const auto& dddConst = topo.dddConstants();
113  int waferTypeL = dddConst.waferType(detid);
114  return cces[waferTypeL-1];
115  }
116  }
117 
118  float getCCE(const HcalGeometry* geom,
119  const DetId& id,
120  const std::vector<float>&cces) {
121  return 1.f;
122  }
123 
124  // Dumps the internals of the SimHit accumulator to the digis for premixing
125  void saveSimHitAccumulator(PHGCSimAccumulator& simResult, const hgc::HGCSimHitDataAccumulator& simData, const std::unordered_set<DetId>& validIds, const float minCharge, const float maxCharge) {
126  constexpr auto nEnergies = std::tuple_size<decltype(hgc_digi::HGCCellInfo().hit_info)>::value;
127  static_assert(nEnergies <= PHGCSimAccumulator::Data::energyMask+1, "PHGCSimAccumulator bit pattern needs to updated");
128  static_assert(hgc_digi::nSamples <= PHGCSimAccumulator::Data::sampleMask+1, "PHGCSimAccumulator bit pattern needs to updated");
129 
130  const float minPackChargeLog = minCharge > 0.f ? std::log(minCharge) : -2;
131  const float maxPackChargeLog = std::log(maxCharge);
133 
134  simResult.reserve(simData.size());
135  // mimicing the digitization
136  for(const auto& id: validIds) {
137  auto found = simData.find(id);
138  if(found == simData.end())
139  continue;
140  // store only non-zero
141  for(size_t iEn = 0; iEn < nEnergies; ++iEn) {
142  const auto& samples = found->second.hit_info[iEn];
143  for(size_t iSample = 0; iSample < hgc_digi::nSamples; ++iSample) {
144  if(samples[iSample] > minCharge) {
145  const auto packed = logintpack::pack16log(samples[iSample], minPackChargeLog, maxPackChargeLog, base);
146  simResult.emplace_back(id.rawId(), iEn, iSample, packed);
147  }
148  }
149  }
150  }
151  simResult.shrink_to_fit();
152  }
153 
154  // Loads the internals of the SimHit accumulator from the digis for premixing
155  void loadSimHitAccumulator(hgc::HGCSimHitDataAccumulator& simData, const PHGCSimAccumulator& simAccumulator, const float minCharge, const float maxCharge, bool setIfZero) {
156  const float minPackChargeLog = minCharge > 0.f ? std::log(minCharge) : -2;
157  const float maxPackChargeLog = std::log(maxCharge);
159 
160  for(const auto& detIdIndexHitInfo: simAccumulator) {
161  auto simIt = simData.emplace(detIdIndexHitInfo.detId(), HGCCellInfo()).first;
162  auto& hit_info = simIt->second.hit_info;
163 
164  size_t iEn = detIdIndexHitInfo.energyIndex();
165  size_t iSample = detIdIndexHitInfo.sampleIndex();
166 
167  float value = logintpack::unpack16log(detIdIndexHitInfo.data(), minPackChargeLog, maxPackChargeLog, base);
168 
169  if(iEn == 0 || !setIfZero) {
170  hit_info[iEn][iSample] += value;
171  }
172  else if(hit_info[iEn][iSample] == 0) {
173  hit_info[iEn][iSample] = value;
174  }
175  }
176  }
177 }
178 
179 //
182  simHitAccumulator_( new HGCSimHitDataAccumulator() ),
183  myDet_(DetId::Forward),
184  mySubDet_(ForwardSubdetector::ForwardEmpty),
185  refSpeed_(0.1*CLHEP::c_light), //[CLHEP::c_light]=mm/ns convert to cm/ns
186  averageOccupancies_(occupancyGuesses),
187  nEvents_(1)
188 {
189  //configure from cfg
190  hitCollection_ = ps.getParameter< std::string >("hitCollection");
191  digiCollection_ = ps.getParameter< std::string >("digiCollection");
192  maxSimHitsAccTime_ = ps.getParameter< uint32_t >("maxSimHitsAccTime");
193  bxTime_ = ps.getParameter< double >("bxTime");
194  geometryType_ = ps.getParameter< uint32_t >("geometryType");
195  digitizationType_ = ps.getParameter< uint32_t >("digitizationType");
196  verbosity_ = ps.getUntrackedParameter< uint32_t >("verbosity",0);
197  tofDelay_ = ps.getParameter< double >("tofDelay");
198  premixStage1_ = ps.getParameter<bool>("premixStage1");
199  premixStage1MinCharge_ = ps.getParameter<double>("premixStage1MinCharge");
200  premixStage1MaxCharge_ = ps.getParameter<double>("premixStage1MaxCharge");
201 
202  std::unordered_set<DetId>().swap(validIds_);
203 
204  iC.consumes<std::vector<PCaloHit> >(edm::InputTag("g4SimHits",hitCollection_));
205  const auto& myCfg_ = ps.getParameter<edm::ParameterSet>("digiCfg");
206 
207  if( myCfg_.existsAs<edm::ParameterSet>("chargeCollectionEfficiencies")) {
208  cce_.clear();
209  const auto& temp = myCfg_.getParameter<edm::ParameterSet>("chargeCollectionEfficiencies").getParameter<std::vector<double>>("values");
210  for( double cce : temp ) {
211  cce_.emplace_back(cce);
212  }
213  } else {
214  std::vector<float>().swap(cce_);
215  }
216 
217  if(hitCollection_.find("HitsEE")!=std::string::npos) {
218  if (geometryType_ == 0) {
220  } else {
222  }
223  theHGCEEDigitizer_=std::make_unique<HGCEEDigitizer>(ps);
224  }
225  if(hitCollection_.find("HitsHEfront")!=std::string::npos) {
226  if (geometryType_ == 0) {
228  } else {
230  }
231  theHGCHEfrontDigitizer_=std::make_unique<HGCHEfrontDigitizer>(ps);
232  }
233  if(hitCollection_.find("HcalHits")!=std::string::npos and geometryType_ == 0) {
235  theHGCHEbackDigitizer_=std::make_unique<HGCHEbackDigitizer>(ps);
236  }
237  if(hitCollection_.find("HitsHEback")!=std::string::npos and geometryType_ == 1) {
239  theHGCHEbackDigitizer_=std::make_unique<HGCHEbackDigitizer>(ps);
240  }
241 }
242 
243 //
245 {
246  // reserve memory for a full detector
247  unsigned idx = getType();
248  simHitAccumulator_->reserve( averageOccupancies_[idx]*validIds_.size() );
249 }
250 
251 //
252 void HGCDigitizer::finalizeEvent(edm::Event& e, edm::EventSetup const& es, CLHEP::HepRandomEngine* hre)
253 {
254  hitRefs_bx0.clear();
255 
256  const CaloSubdetectorGeometry* theGeom = ( nullptr == gHGCal_ ?
257  static_cast<const CaloSubdetectorGeometry*>(gHcal_) :
258  static_cast<const CaloSubdetectorGeometry*>(gHGCal_) );
259 
260  ++nEvents_;
261  unsigned idx = getType();
262  // release memory for unfilled parts of hash table
263  if( validIds_.size()*averageOccupancies_[idx] > simHitAccumulator_->size() ) {
264  simHitAccumulator_->reserve(simHitAccumulator_->size());
265  }
266  //update occupancy guess
267  const double thisOcc = simHitAccumulator_->size()/((double)validIds_.size());
269 
270  if(premixStage1_) {
271  std::unique_ptr<PHGCSimAccumulator> simResult;
272  if(!simHitAccumulator_->empty()) {
273  simResult = std::make_unique<PHGCSimAccumulator>(simHitAccumulator_->begin()->first);
274  saveSimHitAccumulator(*simResult, *simHitAccumulator_, validIds_, premixStage1MinCharge_, premixStage1MaxCharge_);
275  }
276  e.put(std::move(simResult), digiCollection());
277  } else {
278  if( producesEEDigis() ) {
279  auto digiResult = std::make_unique<HGCalDigiCollection>();
280  theHGCEEDigitizer_->run(digiResult,*simHitAccumulator_,theGeom,validIds_,digitizationType_, hre);
281  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " EE hits";
282  e.put(std::move(digiResult),digiCollection());
283  }
284  if( producesHEfrontDigis()) {
285  auto digiResult = std::make_unique<HGCalDigiCollection>();
287  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " HE front hits";
288  e.put(std::move(digiResult),digiCollection());
289  }
290  if( producesHEbackDigis() ) {
291  auto digiResult = std::make_unique<HGCalDigiCollection>();
293  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " HE back hits";
294  e.put(std::move(digiResult),digiCollection());
295  }
296  }
297 
299 }
300 
301 //
302 void HGCDigitizer::accumulate(edm::Event const& e, edm::EventSetup const& eventSetup, CLHEP::HepRandomEngine* hre) {
303 
304  //get inputs
306  e.getByLabel(edm::InputTag("g4SimHits",hitCollection_),hits);
307  if( !hits.isValid() ){
308  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of g4SimHits";
309  return;
310  }
311 
312  //accumulate in-time the main event
313  if( nullptr != gHGCal_ ) {
314  accumulate(hits, 0, gHGCal_, hre);
315  } else if( nullptr != gHcal_ ) {
316  accumulate(hits, 0, gHcal_, hre);
317  } else {
318  throw cms::Exception("BadConfiguration")
319  << "HGCDigitizer is not producing EE, FH, or BH digis!";
320  }
321 }
322 
323 //
324 void HGCDigitizer::accumulate(PileUpEventPrincipal const& e, edm::EventSetup const& eventSetup, CLHEP::HepRandomEngine* hre) {
325 
326  //get inputs
328  e.getByLabel(edm::InputTag("g4SimHits",hitCollection_),hits);
329  if( !hits.isValid() ){
330  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of g4SimHits";
331  return;
332  }
333 
334  //accumulate for the simulated bunch crossing
335  if( nullptr != gHGCal_ ) {
336  accumulate(hits, e.bunchCrossing(), gHGCal_, hre);
337  } else if ( nullptr != gHcal_ ) {
338  accumulate(hits, e.bunchCrossing(), gHcal_, hre);
339  } else {
340  throw cms::Exception("BadConfiguration")
341  << "HGCDigitizer is not producing EE, FH, or BH digis!";
342  }
343 }
344 
345 //
346 template<typename GEOM>
348  int bxCrossing,
349  const GEOM* geom,
350  CLHEP::HepRandomEngine* hre) {
351  if( nullptr == geom ) return;
352 
353 
354 
355  //configuration to apply for the computation of time-of-flight
356  std::array<float, 3> tdcForToAOnset{ {0.f, 0.f, 0.f} };
357  float keV2fC(0.f);
358  bool weightToAbyEnergy= getWeight(tdcForToAOnset, keV2fC);
359 
360  //create list of tuples (pos in container, RECO DetId, time) to be sorted first
361  int nchits=(int)hits->size();
362  std::vector< HGCCaloHitTuple_t > hitRefs;
363  hitRefs.reserve(nchits);
364  for(int i=0; i<nchits; ++i) {
365  const auto& the_hit = hits->at(i);
366 
367  DetId id = simToReco(geom,the_hit.id());
368 
369  if (verbosity_>0) {
370  if (producesEEDigis())
371  edm::LogInfo("HGCDigitizer") << " i/p " << std::hex << the_hit.id() << " o/p " << id.rawId() << std::dec << std::endl;
372  else
373  edm::LogInfo("HGCDigitizer") << " i/p " << std::hex << the_hit.id() << " o/p " << id.rawId() << std::dec << std::endl;
374  }
375 
376  if( 0 != id.rawId() ) {
377  hitRefs.emplace_back( i, id.rawId(), (float)the_hit.time() );
378  }
379  }
380  std::sort(hitRefs.begin(),hitRefs.end(),this->orderByDetIdThenTime);
381 
382  //loop over sorted hits
383  nchits = hitRefs.size();
384  for(int i=0; i<nchits; ++i) {
385  const int hitidx = std::get<0>(hitRefs[i]);
386  const uint32_t id = std::get<1>(hitRefs[i]);
387 
388  //get the data for this cell, if not available then we skip it
389 
390  if( !validIds_.count(id) ) continue;
391  HGCSimHitDataAccumulator::iterator simHitIt = simHitAccumulator_->emplace(id,HGCCellInfo()).first;
392 
393  if(id==0) continue; // to be ignored at RECO level
394 
395  const float toa = std::get<2>(hitRefs[i]);
396  const PCaloHit &hit=hits->at( hitidx );
397  const float charge = hit.energy()*1e6*keV2fC*getCCE(geom,id,cce_);
398 
399  //distance to the center of the detector
400  const float dist2center( getPositionDistance(geom,id) );
401 
402  //hit time: [time()]=ns [centerDist]=cm [refSpeed_]=cm/ns + delay by 1ns
403  //accumulate in 15 buckets of 25ns (9 pre-samples, 1 in-time, 5 post-samples)
404  const float tof = toa-dist2center/refSpeed_+tofDelay_ ;
405  const int itime= std::floor( tof/bxTime_ ) + 9;
406 
407  //no need to add bx crossing - tof comes already corrected from the mixing module
408  //itime += bxCrossing;
409  //itime += 9;
410 
411  if(itime<0 || itime>14) continue;
412 
413  //check if time index is ok and store energy
414  if(itime >= (int)simHitIt->second.hit_info[0].size() ) continue;
415 
416  (simHitIt->second).hit_info[0][itime] += charge;
417 
418 
419  //working version with pileup only for in-time hits
420  int waferThickness = getCellThickness(geom,id);
421  bool orderChanged = false;
422  if(itime == 9){
423  if(hitRefs_bx0[id].empty()){
424  hitRefs_bx0[id].emplace_back(charge, tof);
425  }
426  else if(tof <= hitRefs_bx0[id].back().second){
427  std::vector<std::pair<float, float> >::iterator findPos =
428  std::upper_bound(hitRefs_bx0[id].begin(), hitRefs_bx0[id].end(), std::pair<float, float>(0.f,tof),
429  [](const auto& i, const auto& j){return i.second < j.second;});
430 
431  std::vector<std::pair<float, float> >::iterator insertedPos =
432  hitRefs_bx0[id].insert(findPos, (findPos == hitRefs_bx0[id].begin()) ?
433  std::pair<float, float>(charge,tof) : std::pair<float, float>((findPos-1)->first+charge,tof));
434 
435  for(std::vector<std::pair<float, float> >::iterator step = insertedPos+1; step != hitRefs_bx0[id].end(); ++step){
436  step->first += charge;
437  if(step->first > tdcForToAOnset[waferThickness-1] && step->second != hitRefs_bx0[id].back().second){
438  hitRefs_bx0[id].resize(std::upper_bound(hitRefs_bx0[id].begin(), hitRefs_bx0[id].end(), std::pair<float, float>(0.f,step->second),
439  [](const auto& i, const auto& j){return i.second < j.second;}) - hitRefs_bx0[id].begin());
440  for(auto stepEnd = step+1; stepEnd != hitRefs_bx0[id].end(); ++stepEnd) stepEnd->first += charge;
441  break;
442  }
443  }
444  orderChanged = true;
445  }
446  else{
447  if(hitRefs_bx0[id].back().first <= tdcForToAOnset[waferThickness-1]){
448  hitRefs_bx0[id].emplace_back(hitRefs_bx0[id].back().first+charge, tof);
449  }
450  }
451  }
452 
453  float accChargeForToA = hitRefs_bx0[id].empty() ? 0.f : hitRefs_bx0[id].back().first;
454 
455  //time-of-arrival (check how to be used)
456  if(weightToAbyEnergy) (simHitIt->second).hit_info[1][itime] += charge*tof;
457  else if(accChargeForToA > tdcForToAOnset[waferThickness-1] &&
458  ((simHitIt->second).hit_info[1][itime] == 0 || orderChanged == true) ){
459  float fireTDC = hitRefs_bx0[id].back().second;
460  if (hitRefs_bx0[id].size() > 1){
461  float chargeBeforeThr = 0.f;
462  float tofchargeBeforeThr = 0.f;
463  for(const auto& step : hitRefs_bx0[id]){
464  if(step.first + chargeBeforeThr <= tdcForToAOnset[waferThickness-1]){
465  chargeBeforeThr += step.first;
466  tofchargeBeforeThr = step.second;
467  }
468  else break;
469  }
470  float deltaQ = accChargeForToA - chargeBeforeThr;
471  float deltaTOF = fireTDC - tofchargeBeforeThr;
472  fireTDC = (tdcForToAOnset[waferThickness-1] - chargeBeforeThr) * deltaTOF / deltaQ + tofchargeBeforeThr;
473  }
474  (simHitIt->second).hit_info[1][itime] = fireTDC;
475  }
476 
477  }
478  hitRefs.clear();
479 }
480 
481 void HGCDigitizer::accumulate(const PHGCSimAccumulator& simAccumulator) {
482  //configuration to apply for the computation of time-of-flight
483  std::array<float, 3> tdcForToAOnset{ {0.f, 0.f, 0.f} };
484  float keV2fC(0.f);
485  bool weightToAbyEnergy= getWeight(tdcForToAOnset, keV2fC);
486  loadSimHitAccumulator(*simHitAccumulator_, simAccumulator, premixStage1MinCharge_, premixStage1MaxCharge_, !weightToAbyEnergy);
487 }
488 
489 //
491  //get geometry
493  es.get<CaloGeometryRecord>().get(geom);
494 
495  gHGCal_ = nullptr;
496  gHcal_ = nullptr;
497 
498  if( producesEEDigis() ) gHGCal_ = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(myDet_,mySubDet_));
499  if( producesHEfrontDigis() ) gHGCal_ = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(myDet_,mySubDet_));
500 
501  if( producesHEbackDigis() ) {
502  if (geometryType_ == 0) {
503  gHcal_ = dynamic_cast<const HcalGeometry*>(geom->getSubdetectorGeometry(DetId::Hcal, HcalEndcap));
504  } else {
505  gHGCal_ = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(myDet_,mySubDet_));
506  }
507  }
508 
509  int nadded(0);
510  //valid ID lists
511  if( nullptr != gHGCal_ ) {
512  getValidDetIds( gHGCal_, validIds_ );
513  } else if( nullptr != gHcal_ ) {
514  getValidDetIds( gHcal_, validIds_ );
515  } else {
516  throw cms::Exception("BadConfiguration")
517  << "HGCDigitizer is not producing EE, FH, or BH digis!";
518  }
519 
520  if (verbosity_ > 0)
521  edm::LogInfo("HGCDigitizer")
522  << "Added " << nadded << ":" << validIds_.size()
523  << " detIds without " << hitCollection_
524  << " in first event processed" << std::endl;
525 }
526 
527 //
529 {
530  std::unordered_set<DetId>().swap(validIds_);
531 }
532 
533 //
535 {
536  for( HGCSimHitDataAccumulator::iterator it = simHitAccumulator_->begin(); it!=simHitAccumulator_->end(); it++)
537  {
538  it->second.hit_info[0].fill(0.);
539  it->second.hit_info[1].fill(0.);
540  }
541 }
542 
543 uint32_t HGCDigitizer::getType() const {
544 
546  if (geometryType_ == 0) {
547  switch(mySubDet_) {
549  idx = 0;
550  break;
552  idx = 1;
553  break;
555  idx = 2;
556  break;
557  default:
558  break;
559  }
560  } else {
561  switch(myDet_) {
562  case DetId::HGCalEE:
563  idx = 0;
564  break;
565  case DetId::HGCalHSi:
566  idx = 1;
567  break;
568  case DetId::HGCalHSc:
569  idx = 2;
570  break;
571  default:
572  break;
573  }
574  }
575  return idx;
576 }
577 
578 bool HGCDigitizer::getWeight(std::array<float,3>& tdcForToAOnset,
579  float& keV2fC) const {
580 
581  bool weightToAbyEnergy(false);
582  if (geometryType_ == 0) {
583  switch( mySubDet_ ) {
585  weightToAbyEnergy = theHGCEEDigitizer_->toaModeByEnergy();
586  tdcForToAOnset = theHGCEEDigitizer_->tdcForToAOnset();
587  keV2fC = theHGCEEDigitizer_->keV2fC();
588  break;
590  weightToAbyEnergy = theHGCHEfrontDigitizer_->toaModeByEnergy();
591  tdcForToAOnset = theHGCHEfrontDigitizer_->tdcForToAOnset();
592  keV2fC = theHGCHEfrontDigitizer_->keV2fC();
593  break;
595  weightToAbyEnergy = theHGCHEbackDigitizer_->toaModeByEnergy();
596  tdcForToAOnset = theHGCHEbackDigitizer_->tdcForToAOnset();
597  keV2fC = theHGCHEbackDigitizer_->keV2fC();
598  break;
599  default:
600  break;
601  }
602  } else {
603  switch(myDet_) {
604  case DetId::HGCalEE:
605  weightToAbyEnergy = theHGCEEDigitizer_->toaModeByEnergy();
606  tdcForToAOnset = theHGCEEDigitizer_->tdcForToAOnset();
607  keV2fC = theHGCEEDigitizer_->keV2fC();
608  break;
609  case DetId::HGCalHSi:
610  weightToAbyEnergy = theHGCHEfrontDigitizer_->toaModeByEnergy();
611  tdcForToAOnset = theHGCHEfrontDigitizer_->tdcForToAOnset();
612  keV2fC = theHGCHEfrontDigitizer_->keV2fC();
613  break;
614  case DetId::HGCalHSc:
615  weightToAbyEnergy = theHGCHEbackDigitizer_->toaModeByEnergy();
616  tdcForToAOnset = theHGCHEbackDigitizer_->tdcForToAOnset();
617  keV2fC = theHGCHEbackDigitizer_->keV2fC();
618  break;
619  default:
620  break;
621  }
622  }
623 return weightToAbyEnergy;
624 }
625 
626 template void HGCDigitizer::accumulate<HcalGeometry>(edm::Handle<edm::PCaloHitContainer> const &hits, int bxCrossing,const HcalGeometry *geom, CLHEP::HepRandomEngine* hre);
627 template void HGCDigitizer::accumulate<HGCalGeometry>(edm::Handle<edm::PCaloHitContainer> const &hits, int bxCrossing,const HGCalGeometry *geom, CLHEP::HepRandomEngine* hre);
size
Write out results.
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
uint32_t nEvents_
Definition: HGCDigitizer.h:138
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
T getUntrackedParameter(std::string const &, T const &) const
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:167
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const override
Get the cell geometry of a given detector id. Should return false if not found.
Definition: HcalGeometry.cc:86
bool premixStage1_
Definition: HGCDigitizer.h:100
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
std::vector< float > cce_
Definition: HGCDigitizer.h:140
void finalizeEvent(edm::Event &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *hre)
ForwardSubdetector mySubDet_
Definition: HGCDigitizer.h:125
int digitizationType_
Definition: HGCDigitizer.h:97
std::string hitCollection_
Definition: HGCDigitizer.h:91
void resetSimHitDataAccumulator()
static constexpr unsigned sampleMask
const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const override
Get a list of valid detector ids (for the given subdetector)
Definition: HcalGeometry.cc:76
void emplace_back(unsigned int detId, unsigned short energyIndex, unsigned short sampleIndex, unsigned short data)
void initializeEvent(edm::Event const &e, edm::EventSetup const &c)
actions at the start/end of event
double unpack16log(int16_t i, double lmin, double lmax, uint16_t base=32768)
Definition: liblogintpack.h:55
GlobalPoint getPosition(const DetId &id) const
DetId::Detector myDet_
Definition: HGCDigitizer.h:124
ForwardSubdetector
const HGCalGeometry * gHGCal_
Definition: HGCDigitizer.h:120
bool producesEEDigis()
Definition: HGCDigitizer.h:70
#define constexpr
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
bool getWeight(std::array< float, 3 > &tdcForToAOnset, float &keV2fC) const
std::string digiCollection_
Definition: HGCDigitizer.h:91
const HcalTopology & topology() const
Definition: HcalGeometry.h:117
static constexpr unsigned sampleOffset
static constexpr unsigned energyMask
U second(std::pair< T, U > const &p)
std::unordered_map< uint32_t, HGCCellInfo > HGCSimHitDataAccumulator
T mag() const
Definition: PV3DBase.h:67
void beginRun(const edm::EventSetup &es)
actions at the start/end of run
const HGCalTopology & topology() const
double f[11][100]
std::unordered_set< DetId > validIds_
Definition: HGCDigitizer.h:119
#define end
Definition: vmac.h:39
Definition: value.py:1
base
Make Sure CMSSW is Setup ##.
bool isValid() const
Definition: HandleBase.h:74
constexpr size_t nSamples
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:535
bool producesHEfrontDigis()
Definition: HGCDigitizer.h:72
int waferType(DetId const &id) const
uint32_t getType() const
std::unique_ptr< HGCHEbackDigitizer > theHGCHEbackDigitizer_
Definition: HGCDigitizer.h:115
Definition: DetId.h:18
const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const override
Get a list of valid detector ids (for the given subdetector)
Definition: HGCalGeometry.h:82
void reserve(size_t size)
double premixStage1MinCharge_
Definition: HGCDigitizer.h:103
bool producesHEbackDigis()
Definition: HGCDigitizer.h:74
const HGCalDDDConstants & dddConstants() const
int maxSimHitsAccTime_
Definition: HGCDigitizer.h:108
std::unique_ptr< HGCHEfrontDigitizer > theHGCHEfrontDigitizer_
Definition: HGCDigitizer.h:116
std::unique_ptr< HGCEEDigitizer > theHGCEEDigitizer_
Definition: HGCDigitizer.h:114
void accumulate(edm::Event const &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *hre)
handle SimHit accumulation
#define begin
Definition: vmac.h:32
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
T get() const
Definition: EventSetup.h:63
std::map< uint32_t, std::vector< std::pair< float, float > > > hitRefs_bx0
Definition: HGCDigitizer.h:142
HGCDigitizer(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
step
std::array< double, 3 > averageOccupancies_
Definition: HGCDigitizer.h:137
DetId relabel(const uint32_t testId) const
std::unique_ptr< hgc::HGCSimHitDataAccumulator > simHitAccumulator_
Definition: HGCDigitizer.h:110
static bool orderByDetIdThenTime(const HGCCaloHitTuple_t &a, const HGCCaloHitTuple_t &b)
Definition: HGCDigitizer.h:38
double bxTime_
Definition: HGCDigitizer.h:109
int16_t pack16log(double x, double lmin, double lmax, uint16_t base=32768)
Definition: liblogintpack.h:26
double premixStage1MaxCharge_
Definition: HGCDigitizer.h:105
std::string digiCollection()
Definition: HGCDigitizer.h:76
const HcalGeometry * gHcal_
Definition: HGCDigitizer.h:121
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:510
uint32_t verbosity_
Definition: HGCDigitizer.h:128
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39