CMS 3D CMS Logo

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