CMS 3D CMS Logo

HGCDigitizer.cc
Go to the documentation of this file.
20 
21 #include <algorithm>
22 #include <boost/foreach.hpp>
24 
25 using namespace hgc_digi;
26 
27 namespace {
28 
29  constexpr std::array<double,3> occupancyGuesses = { { 0.5,0.2,0.2 } };
30 
31 
32  float getPositionDistance(const HGCalGeometry* geom, const DetId& id) {
33  return geom->getPosition(id).mag();
34  }
35 
36  float getPositionDistance(const HcalGeometry* geom, const DetId& id) {
37  return geom->getGeometry(id)->getPosition().mag();
38  }
39 
40  int getCellThickness(const HGCalGeometry* geom, const DetId& detid ) {
41  const auto& topo = geom->topology();
42  const auto& dddConst = topo.dddConstants();
43  uint32_t id(detid.rawId());
44  HGCalDetId hid(id);
45  int wafer = HGCalDetId(id).wafer();
46  int waferTypeL = dddConst.waferTypeL(wafer);
47  return waferTypeL;
48  }
49 
50  int getCellThickness(const HcalGeometry* geom, const DetId& detid ) {
51  return 1;
52  }
53 
54  void getValidDetIds(const HGCalGeometry* geom, std::unordered_set<DetId>& valid) {
55  const std::vector<DetId>& ids = geom->getValidDetIds();
56  valid.reserve(ids.size());
57  valid.insert(ids.begin(),ids.end());
58  }
59 
60  void getValidDetIds(const HcalGeometry* geom, std::unordered_set<DetId>& valid) {
61  const std::vector<DetId>& ids = geom->getValidDetIds();
62  for( const auto& id : ids ) {
63  if( HcalEndcap == id.subdetId() &&
64  DetId::Hcal == id.det() )
65  valid.emplace(id);
66  }
67  valid.reserve(valid.size());
68  }
69 
70  DetId simToReco(const HcalGeometry* geom, unsigned simid) {
71  DetId result(0);
72  const auto& topo = geom->topology();
73  const auto* dddConst = topo.dddConstants();
74  HcalDetId id = HcalHitRelabeller::relabel(simid,dddConst);
75 
76  if (id.subdet()==int(HcalEndcap)) {
77  result = id;
78  }
79 
80  return result;
81  }
82 
83  DetId simToReco(const HGCalGeometry* geom, unsigned simId) {
84  DetId result(0);
85  const auto& topo = geom->topology();
86  const auto& dddConst = topo.dddConstants();
87 
88  int subdet(DetId(simId).subdetId()), layer, cell, sec, subsec, zp;
89 
90  const bool isSqr = (dddConst.geomMode() == HGCalGeometryMode::Square);
91  if (isSqr) {
92  HGCalTestNumbering::unpackSquareIndex(simId, zp, layer, sec, subsec, cell);
93  } else {
94  HGCalTestNumbering::unpackHexagonIndex(simId, subdet, zp, layer, sec, subsec, cell);
95  //sec is wafer and subsec is celltyp
96  }
97  //skip this hit if after ganging it is not valid
98  std::pair<int,int> recoLayerCell=dddConst.simToReco(cell,layer,sec,topo.detectorType());
99  cell = recoLayerCell.first;
100  layer = recoLayerCell.second;
101  if (layer<0 || cell<0) {
102  return result;
103  }
104 
105  //assign the RECO DetId
106  result = HGCalDetId((ForwardSubdetector)subdet,zp,layer,subsec,sec,cell);
107 
108  return result;
109  }
110 
111  float getCCE(const HGCalGeometry* geom,
112  const DetId& detid,
113  const std::vector<float>&cces) {
114  if( cces.empty() ) return 1.f;
115  const auto& topo = geom->topology();
116  const auto& dddConst = topo.dddConstants();
117  uint32_t id(detid.rawId());
118  HGCalDetId hid(id);
119  int wafer = HGCalDetId(id).wafer();
120  int waferTypeL = dddConst.waferTypeL(wafer);
121  return cces[waferTypeL-1];
122  }
123 
124  float getCCE(const HcalGeometry* geom,
125  const DetId& id,
126  const std::vector<float>&cces) {
127  return 1.f;
128  }
129 
130 }
131 
132 //
135  simHitAccumulator_( new HGCSimHitDataAccumulator() ),
136  mySubDet_(ForwardSubdetector::ForwardEmpty),
137  refSpeed_(0.1*CLHEP::c_light), //[CLHEP::c_light]=mm/ns convert to cm/ns
138  averageOccupancies_(occupancyGuesses),
139  nEvents_(1)
140 {
141  //configure from cfg
142  hitCollection_ = ps.getParameter< std::string >("hitCollection");
143  digiCollection_ = ps.getParameter< std::string >("digiCollection");
144  maxSimHitsAccTime_ = ps.getParameter< uint32_t >("maxSimHitsAccTime");
145  bxTime_ = ps.getParameter< double >("bxTime");
146  digitizationType_ = ps.getParameter< uint32_t >("digitizationType");
147  verbosity_ = ps.getUntrackedParameter< uint32_t >("verbosity",0);
148  tofDelay_ = ps.getParameter< double >("tofDelay");
149 
150  std::unordered_set<DetId>().swap(validIds_);
151 
152  iC.consumes<std::vector<PCaloHit> >(edm::InputTag("g4SimHits",hitCollection_));
153  const auto& myCfg_ = ps.getParameter<edm::ParameterSet>("digiCfg");
154 
155  if( myCfg_.existsAs<std::vector<double> >( "chargeCollectionEfficiencies" ) ) {
156  cce_.clear();
157  const auto& temp = myCfg_.getParameter<std::vector<double> >("chargeCollectionEfficiencies");
158  for( double cce : temp ) {
159  cce_.push_back(cce);
160  }
161  } else {
162  std::vector<float>().swap(cce_);
163  }
164 
165  if(hitCollection_.find("HitsEE")!=std::string::npos) {
167  theHGCEEDigitizer_=std::unique_ptr<HGCEEDigitizer>(new HGCEEDigitizer(ps) );
168  }
169  if(hitCollection_.find("HitsHEfront")!=std::string::npos)
170  {
172  theHGCHEfrontDigitizer_=std::unique_ptr<HGCHEfrontDigitizer>(new HGCHEfrontDigitizer(ps) );
173  }
174  if(hitCollection_.find("HcalHits")!=std::string::npos)
175  {
177  theHGCHEbackDigitizer_=std::unique_ptr<HGCHEbackDigitizer>(new HGCHEbackDigitizer(ps) );
178  }
179 }
180 
181 //
183 {
184  // reserve memory for a full detector
186  switch(mySubDet_) {
188  idx = 0;
189  break;
191  idx = 1;
192  break;
194  idx = 2;
195  break;
196  default:
197  break;
198  }
199  simHitAccumulator_->reserve( averageOccupancies_[idx]*validIds_.size() );
200 }
201 
202 //
203 void HGCDigitizer::finalizeEvent(edm::Event& e, edm::EventSetup const& es, CLHEP::HepRandomEngine* hre)
204 {
205  hitRefs_bx0.clear();
206 
207  const CaloSubdetectorGeometry* theGeom = ( nullptr == gHGCal_ ?
208  static_cast<const CaloSubdetectorGeometry*>(gHcal_) :
209  static_cast<const CaloSubdetectorGeometry*>(gHGCal_) );
210 
211  ++nEvents_;
213  switch(mySubDet_) {
215  idx = 0;
216  break;
218  idx = 1;
219  break;
221  idx = 2;
222  break;
223  default:
224  break;
225  }
226  // release memory for unfilled parts of hash table
227  if( validIds_.size()*averageOccupancies_[idx] > simHitAccumulator_->size() ) {
228  simHitAccumulator_->reserve(simHitAccumulator_->size());
229  }
230  //update occupancy guess
231  const double thisOcc = simHitAccumulator_->size()/((double)validIds_.size());
233 
234  if( producesEEDigis() )
235  {
236  std::unique_ptr<HGCEEDigiCollection> digiResult(new HGCEEDigiCollection() );
237  theHGCEEDigitizer_->run(digiResult,*simHitAccumulator_,theGeom,validIds_,digitizationType_, hre);
238  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " EE hits";
239  e.put(std::move(digiResult),digiCollection());
240  }
241  if( producesHEfrontDigis())
242  {
243  std::unique_ptr<HGCHEDigiCollection> digiResult(new HGCHEDigiCollection() );
245  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " HE front hits";
246  e.put(std::move(digiResult),digiCollection());
247  }
248  if( producesHEbackDigis() )
249  {
250  std::unique_ptr<HGCBHDigiCollection> digiResult(new HGCBHDigiCollection() );
252  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " HE back hits";
253  e.put(std::move(digiResult),digiCollection());
254  }
255 
257 }
258 
259 //
260 void HGCDigitizer::accumulate(edm::Event const& e, edm::EventSetup const& eventSetup, CLHEP::HepRandomEngine* hre) {
261 
262  //get inputs
264  e.getByLabel(edm::InputTag("g4SimHits",hitCollection_),hits);
265  if( !hits.isValid() ){
266  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of g4SimHits";
267  return;
268  }
269 
270  //accumulate in-time the main event
271  if( nullptr != gHGCal_ ) {
272  accumulate(hits, 0, gHGCal_, hre);
273  } else if( nullptr != gHcal_ ) {
274  accumulate(hits, 0, gHcal_, hre);
275  } else {
276  throw cms::Exception("BadConfiguration")
277  << "HGCDigitizer is not producing EE, FH, or BH digis!";
278  }
279 }
280 
281 //
282 void HGCDigitizer::accumulate(PileUpEventPrincipal const& e, edm::EventSetup const& eventSetup, CLHEP::HepRandomEngine* hre) {
283 
284  //get inputs
286  e.getByLabel(edm::InputTag("g4SimHits",hitCollection_),hits);
287  if( !hits.isValid() ){
288  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of g4SimHits";
289  return;
290  }
291 
292  //accumulate for the simulated bunch crossing
293  if( nullptr != gHGCal_ ) {
294  accumulate(hits, e.bunchCrossing(), gHGCal_, hre);
295  } else if ( nullptr != gHcal_ ) {
296  accumulate(hits, e.bunchCrossing(), gHcal_, hre);
297  } else {
298  throw cms::Exception("BadConfiguration")
299  << "HGCDigitizer is not producing EE, FH, or BH digis!";
300  }
301 }
302 
303 //
304 template<typename GEOM>
306  int bxCrossing,
307  const GEOM* geom,
308  CLHEP::HepRandomEngine* hre) {
309  if( nullptr == geom ) return;
310 
311 
312 
313  //configuration to apply for the computation of time-of-flight
314  bool weightToAbyEnergy(false);
315  std::array<float, 3> tdcForToAOnset{ {0.f, 0.f, 0.f} };
316  float keV2fC(0.f);
317  switch( mySubDet_ ) {
319  weightToAbyEnergy = theHGCEEDigitizer_->toaModeByEnergy();
320  tdcForToAOnset = theHGCEEDigitizer_->tdcForToAOnset();
321  keV2fC = theHGCEEDigitizer_->keV2fC();
322  break;
324  weightToAbyEnergy = theHGCHEfrontDigitizer_->toaModeByEnergy();
325  tdcForToAOnset = theHGCHEfrontDigitizer_->tdcForToAOnset();
326  keV2fC = theHGCHEfrontDigitizer_->keV2fC();
327  break;
329  weightToAbyEnergy = theHGCHEbackDigitizer_->toaModeByEnergy();
330  tdcForToAOnset = theHGCHEbackDigitizer_->tdcForToAOnset();
331  keV2fC = theHGCHEbackDigitizer_->keV2fC();
332  break;
333  default:
334  break;
335  }
336 
337  //create list of tuples (pos in container, RECO DetId, time) to be sorted first
338  int nchits=(int)hits->size();
339  std::vector< HGCCaloHitTuple_t > hitRefs;
340  hitRefs.reserve(nchits);
341  for(int i=0; i<nchits; ++i) {
342  const auto& the_hit = hits->at(i);
343 
344  DetId id = simToReco(geom,the_hit.id());
345 
346  if (verbosity_>0) {
347  if (producesEEDigis())
348  edm::LogInfo("HGCDigitizer") << " i/p " << std::hex << the_hit.id() << " o/p " << id.rawId() << std::dec << std::endl;
349  else
350  edm::LogInfo("HGCDigitizer") << " i/p " << std::hex << the_hit.id() << " o/p " << id.rawId() << std::dec << std::endl;
351  }
352 
353  if( 0 != id.rawId() ) {
354  hitRefs.emplace_back( i, id.rawId(), (float)the_hit.time() );
355  }
356  }
357  std::sort(hitRefs.begin(),hitRefs.end(),this->orderByDetIdThenTime);
358 
359  //loop over sorted hits
360  nchits = hitRefs.size();
361  for(int i=0; i<nchits; ++i) {
362  const int hitidx = std::get<0>(hitRefs[i]);
363  const uint32_t id = std::get<1>(hitRefs[i]);
364 
365  //get the data for this cell, if not available then we skip it
366 
367  if( !validIds_.count(id) ) continue;
368  HGCSimHitDataAccumulator::iterator simHitIt = simHitAccumulator_->emplace(id,HGCCellInfo()).first;
369 
370  if(id==0) continue; // to be ignored at RECO level
371 
372  const float toa = std::get<2>(hitRefs[i]);
373  const PCaloHit &hit=hits->at( hitidx );
374  const float charge = hit.energy()*1e6*keV2fC*getCCE(geom,id,cce_);
375 
376  //distance to the center of the detector
377  const float dist2center( getPositionDistance(geom,id) );
378 
379  //hit time: [time()]=ns [centerDist]=cm [refSpeed_]=cm/ns + delay by 1ns
380  //accumulate in 15 buckets of 25ns (9 pre-samples, 1 in-time, 5 post-samples)
381  const float tof = toa-dist2center/refSpeed_+tofDelay_ ;
382  const int itime= std::floor( tof/bxTime_ ) + 9;
383 
384  //no need to add bx crossing - tof comes already corrected from the mixing module
385  //itime += bxCrossing;
386  //itime += 9;
387 
388  if(itime<0 || itime>14) continue;
389 
390  //check if time index is ok and store energy
391  if(itime >= (int)simHitIt->second.hit_info[0].size() ) continue;
392 
393  (simHitIt->second).hit_info[0][itime] += charge;
394 
395 
396  //working version with pileup only for in-time hits
397  int waferThickness = getCellThickness(geom,id);
398  bool orderChanged = false;
399  if(itime == 9){
400  if(hitRefs_bx0[id].empty()){
401  hitRefs_bx0[id].push_back(std::pair<float, float>(charge, tof));
402  }
403  else if(tof <= hitRefs_bx0[id].back().second){
404  std::vector<std::pair<float, float> >::iterator findPos =
405  std::upper_bound(hitRefs_bx0[id].begin(), hitRefs_bx0[id].end(), std::pair<float, float>(0.f,tof),
406  [](const auto& i, const auto& j){return i.second < j.second;});
407 
408  std::vector<std::pair<float, float> >::iterator insertedPos =
409  hitRefs_bx0[id].insert(findPos, (findPos == hitRefs_bx0[id].begin()) ?
410  std::pair<float, float>(charge,tof) : std::pair<float, float>((findPos-1)->first+charge,tof));
411 
412  for(std::vector<std::pair<float, float> >::iterator step = insertedPos+1; step != hitRefs_bx0[id].end(); ++step){
413  step->first += charge;
414  if(step->first > tdcForToAOnset[waferThickness-1] && step->second != hitRefs_bx0[id].back().second){
415  hitRefs_bx0[id].resize(std::upper_bound(hitRefs_bx0[id].begin(), hitRefs_bx0[id].end(), std::pair<float, float>(0.f,step->second),
416  [](const auto& i, const auto& j){return i.second < j.second;}) - hitRefs_bx0[id].begin());
417  for(auto stepEnd = step+1; stepEnd != hitRefs_bx0[id].end(); ++stepEnd) stepEnd->first += charge;
418  break;
419  }
420  }
421  orderChanged = true;
422  }
423  else{
424  if(hitRefs_bx0[id].back().first <= tdcForToAOnset[waferThickness-1]){
425  hitRefs_bx0[id].push_back(std::pair<float, float>(hitRefs_bx0[id].back().first+charge, tof));
426  }
427  }
428  }
429 
430  float accChargeForToA = hitRefs_bx0[id].empty() ? 0.f : hitRefs_bx0[id].back().first;
431 
432  //time-of-arrival (check how to be used)
433  if(weightToAbyEnergy) (simHitIt->second).hit_info[1][itime] += charge*tof;
434  else if(accChargeForToA > tdcForToAOnset[waferThickness-1] &&
435  ((simHitIt->second).hit_info[1][itime] == 0 || orderChanged == true) ){
436  float fireTDC = hitRefs_bx0[id].back().second;
437  if (hitRefs_bx0[id].size() > 1){
438  float chargeBeforeThr = 0.f;
439  float tofchargeBeforeThr = 0.f;
440  for(const auto& step : hitRefs_bx0[id]){
441  if(step.first + chargeBeforeThr <= tdcForToAOnset[waferThickness-1]){
442  chargeBeforeThr += step.first;
443  tofchargeBeforeThr = step.second;
444  }
445  else break;
446  }
447  float deltaQ = accChargeForToA - chargeBeforeThr;
448  float deltaTOF = fireTDC - tofchargeBeforeThr;
449  fireTDC = (tdcForToAOnset[waferThickness-1] - chargeBeforeThr) * deltaTOF / deltaQ + tofchargeBeforeThr;
450  }
451  (simHitIt->second).hit_info[1][itime] = fireTDC;
452  }
453 
454  }
455  hitRefs.clear();
456 }
457 
458 //
460 {
461  //get geometry
463  es.get<CaloGeometryRecord>().get(geom);
464 
465  gHGCal_ = nullptr;
466  gHcal_ = nullptr;
467 
468  if( producesEEDigis() ) gHGCal_ = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(DetId::Forward, HGCEE));
469  if( producesHEfrontDigis() ) gHGCal_ = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(DetId::Forward, HGCHEF));
470  if( producesHEbackDigis() ) gHcal_ = dynamic_cast<const HcalGeometry*>(geom->getSubdetectorGeometry(DetId::Hcal, HcalEndcap));
471 
472  int nadded(0);
473  //valid ID lists
474  if( nullptr != gHGCal_ ) {
475  getValidDetIds( gHGCal_, validIds_ );
476  } else if( nullptr != gHcal_ ) {
477  getValidDetIds( gHcal_, validIds_ );
478  } else {
479  throw cms::Exception("BadConfiguration")
480  << "HGCDigitizer is not producing EE, FH, or BH digis!";
481  }
482 
483  if (verbosity_ > 0)
484  edm::LogInfo("HGCDigitizer")
485  << "Added " << nadded << ":" << validIds_.size()
486  << " detIds without " << hitCollection_
487  << " in first event processed" << std::endl;
488 }
489 
490 //
492 {
493  std::unordered_set<DetId>().swap(validIds_);
494 }
495 
496 //
498 {
499  for( HGCSimHitDataAccumulator::iterator it = simHitAccumulator_->begin(); it!=simHitAccumulator_->end(); it++)
500  {
501  it->second.hit_info[0].fill(0.);
502  it->second.hit_info[1].fill(0.);
503  }
504 }
505 
506 template void HGCDigitizer::accumulate<HcalGeometry>(edm::Handle<edm::PCaloHitContainer> const &hits, int bxCrossing,const HcalGeometry *geom, CLHEP::HepRandomEngine* hre);
507 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:116
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:45
T getUntrackedParameter(std::string const &, T const &) const
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:161
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
std::vector< float > cce_
Definition: HGCDigitizer.h:118
edm::SortedCollection< HGCEEDataFrame > HGCEEDigiCollection
void finalizeEvent(edm::Event &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *hre)
ForwardSubdetector mySubDet_
Definition: HGCDigitizer.h:103
int digitizationType_
Definition: HGCDigitizer.h:84
std::string hitCollection_
Definition: HGCDigitizer.h:81
void resetSimHitDataAccumulator()
virtual const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
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 initializeEvent(edm::Event const &e, edm::EventSetup const &c)
actions at the start/end of event
GlobalPoint getPosition(const DetId &id) const
ForwardSubdetector
const HGCalGeometry * gHGCal_
Definition: HGCDigitizer.h:99
bool producesEEDigis()
Definition: HGCDigitizer.h:67
#define constexpr
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
std::string digiCollection_
Definition: HGCDigitizer.h:81
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
const HcalTopology & topology() const
Definition: HcalGeometry.h:119
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
U second(std::pair< T, U > const &p)
std::unordered_map< uint32_t, HGCCellInfo > HGCSimHitDataAccumulator
T mag() const
Definition: PV3DBase.h:67
edm::SortedCollection< HGCBHDataFrame > HGCBHDigiCollection
void beginRun(const edm::EventSetup &es)
actions at the start/end of run
const HGCalTopology & topology() const
Definition: HGCalGeometry.h:96
double f[11][100]
std::unordered_set< DetId > validIds_
Definition: HGCDigitizer.h:98
#define end
Definition: vmac.h:37
int wafer() const
get the wafer #
Definition: HGCalDetId.h:42
bool isValid() const
Definition: HandleBase.h:74
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:464
bool producesHEfrontDigis()
Definition: HGCDigitizer.h:68
std::unique_ptr< HGCHEbackDigitizer > theHGCHEbackDigitizer_
Definition: HGCDigitizer.h:94
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:76
bool producesHEbackDigis()
Definition: HGCDigitizer.h:69
const HGCalDDDConstants & dddConstants() const
static void unpackSquareIndex(const uint32_t &idx, int &z, int &lay, int &sec, int &subsec, int &cell)
const T & get() const
Definition: EventSetup.h:55
int maxSimHitsAccTime_
Definition: HGCDigitizer.h:87
std::unique_ptr< HGCHEfrontDigitizer > theHGCHEfrontDigitizer_
Definition: HGCDigitizer.h:95
std::unique_ptr< HGCEEDigitizer > theHGCEEDigitizer_
Definition: HGCDigitizer.h:93
void accumulate(edm::Event const &e, edm::EventSetup const &c, CLHEP::HepRandomEngine *hre)
handle SimHit accumulation
#define begin
Definition: vmac.h:30
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
std::map< uint32_t, std::vector< std::pair< float, float > > > hitRefs_bx0
Definition: HGCDigitizer.h:120
HGCDigitizer(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
step
std::array< double, 3 > averageOccupancies_
Definition: HGCDigitizer.h:115
DetId relabel(const uint32_t testId) const
std::unique_ptr< hgc::HGCSimHitDataAccumulator > simHitAccumulator_
Definition: HGCDigitizer.h:89
static bool orderByDetIdThenTime(const HGCCaloHitTuple_t &a, const HGCCaloHitTuple_t &b)
Definition: HGCDigitizer.h:37
double bxTime_
Definition: HGCDigitizer.h:88
std::string digiCollection()
Definition: HGCDigitizer.h:70
edm::SortedCollection< HGCHEDataFrame > HGCHEDigiCollection
const HcalGeometry * gHcal_
Definition: HGCDigitizer.h:100
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:106