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  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  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  void getValidDetIds(const HcalGeometry* geom, std::unordered_set<DetId>& valid) {
46  const std::vector<DetId>& ids = geom->getValidDetIds();
47  for( const auto& id : ids ) {
48  if( HcalEndcap == id.subdetId() &&
49  DetId::Hcal == id.det() )
50  valid.emplace(id);
51  }
52  valid.reserve(valid.size());
53  }
54 
55  DetId simToReco(const HcalGeometry* geom, unsigned simid) {
56  DetId result(0);
57  const auto& topo = geom->topology();
58  const auto* dddConst = topo.dddConstants();
59  HcalDetId id = HcalHitRelabeller::relabel(simid,dddConst);
60 
61  if (id.subdet()==int(HcalEndcap)) {
62  result = id;
63  }
64 
65  return result;
66  }
67 
68  DetId simToReco(const HGCalGeometry* geom, unsigned simId) {
69  DetId result(0);
70  const auto& topo = geom->topology();
71  const auto& dddConst = topo.dddConstants();
72 
73  int subdet(DetId(simId).subdetId()), layer, cell, sec, subsec, zp;
74 
75  const bool isSqr = (dddConst.geomMode() == HGCalGeometryMode::Square);
76  if (isSqr) {
77  HGCalTestNumbering::unpackSquareIndex(simId, zp, layer, sec, subsec, cell);
78  } else {
79  HGCalTestNumbering::unpackHexagonIndex(simId, subdet, zp, layer, sec, subsec, cell);
80  //sec is wafer and subsec is celltyp
81  }
82  //skip this hit if after ganging it is not valid
83  std::pair<int,int> recoLayerCell=dddConst.simToReco(cell,layer,sec,topo.detectorType());
84  cell = recoLayerCell.first;
85  layer = recoLayerCell.second;
86  if (layer<0 || cell<0) {
87  return result;
88  }
89 
90  //assign the RECO DetId
91  result = HGCalDetId((ForwardSubdetector)subdet,zp,layer,subsec,sec,cell);
92 
93  return result;
94  }
95 
96  float getCCE(const HGCalGeometry* geom,
97  const DetId& detid,
98  const std::vector<float>&cces) {
99  if( cces.empty() ) return 1.f;
100  const auto& topo = geom->topology();
101  const auto& dddConst = topo.dddConstants();
102  uint32_t id(detid.rawId());
103  HGCalDetId hid(id);
104  int wafer = HGCalDetId(id).wafer();
105  int waferTypeL = dddConst.waferTypeL(wafer);
106  return cces[waferTypeL-1];
107  }
108 
109  float getCCE(const HcalGeometry* geom,
110  const DetId& id,
111  const std::vector<float>&cces) {
112  return 1.f;
113  }
114 
115 }
116 
117 //
120  simHitAccumulator_( new HGCSimHitDataAccumulator() ),
121  mySubDet_(ForwardSubdetector::ForwardEmpty),
122  refSpeed_(0.1*CLHEP::c_light), //[CLHEP::c_light]=mm/ns convert to cm/ns
123  averageOccupancies_(occupancyGuesses),
124  nEvents_(1)
125 {
126  //configure from cfg
127  hitCollection_ = ps.getParameter< std::string >("hitCollection");
128  digiCollection_ = ps.getParameter< std::string >("digiCollection");
129  maxSimHitsAccTime_ = ps.getParameter< uint32_t >("maxSimHitsAccTime");
130  bxTime_ = ps.getParameter< double >("bxTime");
131  digitizationType_ = ps.getParameter< uint32_t >("digitizationType");
132  verbosity_ = ps.getUntrackedParameter< uint32_t >("verbosity",0);
133  tofDelay_ = ps.getParameter< double >("tofDelay");
134 
135  std::unordered_set<DetId>().swap(validIds_);
136 
137  iC.consumes<std::vector<PCaloHit> >(edm::InputTag("g4SimHits",hitCollection_));
138  const auto& myCfg_ = ps.getParameter<edm::ParameterSet>("digiCfg");
139 
140  if( myCfg_.existsAs<std::vector<double> >( "chargeCollectionEfficiencies" ) ) {
141  cce_.clear();
142  const auto& temp = myCfg_.getParameter<std::vector<double> >("chargeCollectionEfficiencies");
143  for( double cce : temp ) {
144  cce_.push_back(cce);
145  }
146  } else {
147  std::vector<float>().swap(cce_);
148  }
149 
150  if(hitCollection_.find("HitsEE")!=std::string::npos) {
152  theHGCEEDigitizer_=std::unique_ptr<HGCEEDigitizer>(new HGCEEDigitizer(ps) );
153  }
154  if(hitCollection_.find("HitsHEfront")!=std::string::npos)
155  {
157  theHGCHEfrontDigitizer_=std::unique_ptr<HGCHEfrontDigitizer>(new HGCHEfrontDigitizer(ps) );
158  }
159  if(hitCollection_.find("HcalHits")!=std::string::npos)
160  {
162  theHGCHEbackDigitizer_=std::unique_ptr<HGCHEbackDigitizer>(new HGCHEbackDigitizer(ps) );
163  }
164 }
165 
166 //
168 {
169  // reserve memory for a full detector
171  switch(mySubDet_) {
173  idx = 0;
174  break;
176  idx = 1;
177  break;
179  idx = 2;
180  break;
181  default:
182  break;
183  }
184  simHitAccumulator_->reserve( averageOccupancies_[idx]*validIds_.size() );
185 }
186 
187 //
188 void HGCDigitizer::finalizeEvent(edm::Event& e, edm::EventSetup const& es, CLHEP::HepRandomEngine* hre)
189 {
190 
191  const CaloSubdetectorGeometry* theGeom = ( nullptr == gHGCal_ ?
192  static_cast<const CaloSubdetectorGeometry*>(gHcal_) :
193  static_cast<const CaloSubdetectorGeometry*>(gHGCal_) );
194 
195  ++nEvents_;
197  switch(mySubDet_) {
199  idx = 0;
200  break;
202  idx = 1;
203  break;
205  idx = 2;
206  break;
207  default:
208  break;
209  }
210  // release memory for unfilled parts of hash table
211  if( validIds_.size()*averageOccupancies_[idx] > simHitAccumulator_->size() ) {
212  simHitAccumulator_->reserve(simHitAccumulator_->size());
213  }
214  //update occupancy guess
215  const double thisOcc = simHitAccumulator_->size()/((double)validIds_.size());
217 
218  if( producesEEDigis() )
219  {
220  std::unique_ptr<HGCEEDigiCollection> digiResult(new HGCEEDigiCollection() );
221  theHGCEEDigitizer_->run(digiResult,*simHitAccumulator_,theGeom,validIds_,digitizationType_, hre);
222  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " EE hits";
223  e.put(std::move(digiResult),digiCollection());
224  }
225  if( producesHEfrontDigis())
226  {
227  std::unique_ptr<HGCHEDigiCollection> digiResult(new HGCHEDigiCollection() );
229  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " HE front hits";
230  e.put(std::move(digiResult),digiCollection());
231  }
232  if( producesHEbackDigis() )
233  {
234  std::unique_ptr<HGCBHDigiCollection> digiResult(new HGCBHDigiCollection() );
236  edm::LogInfo("HGCDigitizer") << " @ finalize event - produced " << digiResult->size() << " HE back hits";
237  e.put(std::move(digiResult),digiCollection());
238  }
239 
241 }
242 
243 //
244 void HGCDigitizer::accumulate(edm::Event const& e, edm::EventSetup const& eventSetup, CLHEP::HepRandomEngine* hre) {
245 
246  //get inputs
248  e.getByLabel(edm::InputTag("g4SimHits",hitCollection_),hits);
249  if( !hits.isValid() ){
250  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of g4SimHits";
251  return;
252  }
253 
254  //accumulate in-time the main event
255  if( nullptr != gHGCal_ ) {
256  accumulate(hits, 0, gHGCal_, hre);
257  } else if( nullptr != gHcal_ ) {
258  accumulate(hits, 0, gHcal_, hre);
259  } else {
260  throw cms::Exception("BadConfiguration")
261  << "HGCDigitizer is not producing EE, FH, or BH digis!";
262  }
263 }
264 
265 //
266 void HGCDigitizer::accumulate(PileUpEventPrincipal const& e, edm::EventSetup const& eventSetup, CLHEP::HepRandomEngine* hre) {
267 
268  //get inputs
270  e.getByLabel(edm::InputTag("g4SimHits",hitCollection_),hits);
271  if( !hits.isValid() ){
272  edm::LogError("HGCDigitizer") << " @ accumulate : can't find " << hitCollection_ << " collection of g4SimHits";
273  return;
274  }
275 
276  //accumulate for the simulated bunch crossing
277  if( nullptr != gHGCal_ ) {
278  accumulate(hits, e.bunchCrossing(), gHGCal_, hre);
279  } else if ( nullptr != gHcal_ ) {
280  accumulate(hits, e.bunchCrossing(), gHcal_, hre);
281  } else {
282  throw cms::Exception("BadConfiguration")
283  << "HGCDigitizer is not producing EE, FH, or BH digis!";
284  }
285 }
286 
287 //
288 template<typename GEOM>
290  int bxCrossing,
291  const GEOM* geom,
292  CLHEP::HepRandomEngine* hre) {
293  if( nullptr == geom ) return;
294 
295 
296 
297  //configuration to apply for the computation of time-of-flight
298  bool weightToAbyEnergy(false);
299  float tdcOnset(0.f),keV2fC(0.f);
300  switch( mySubDet_ ) {
302  weightToAbyEnergy = theHGCEEDigitizer_->toaModeByEnergy();
303  tdcOnset = theHGCEEDigitizer_->tdcOnset();
304  keV2fC = theHGCEEDigitizer_->keV2fC();
305  break;
307  weightToAbyEnergy = theHGCHEfrontDigitizer_->toaModeByEnergy();
308  tdcOnset = theHGCHEfrontDigitizer_->tdcOnset();
309  keV2fC = theHGCHEfrontDigitizer_->keV2fC();
310  break;
312  weightToAbyEnergy = theHGCHEbackDigitizer_->toaModeByEnergy();
313  tdcOnset = theHGCHEbackDigitizer_->tdcOnset();
314  keV2fC = theHGCHEbackDigitizer_->keV2fC();
315  break;
316  default:
317  break;
318  }
319 
320  //create list of tuples (pos in container, RECO DetId, time) to be sorted first
321  int nchits=(int)hits->size();
322  std::vector< HGCCaloHitTuple_t > hitRefs;
323  hitRefs.reserve(nchits);
324  for(int i=0; i<nchits; ++i) {
325  const auto& the_hit = hits->at(i);
326 
327  DetId id = simToReco(geom,the_hit.id());
328 
329  if (verbosity_>0) {
330  if (producesEEDigis())
331  edm::LogInfo("HGCDigitizer") << " i/p " << std::hex << the_hit.id() << " o/p " << id.rawId() << std::dec << std::endl;
332  else
333  edm::LogInfo("HGCDigitizer") << " i/p " << std::hex << the_hit.id() << " o/p " << id.rawId() << std::dec << std::endl;
334  }
335 
336  if( 0 != id.rawId() ) {
337  hitRefs.emplace_back( i, id.rawId(), (float)the_hit.time() );
338  }
339  }
340  std::sort(hitRefs.begin(),hitRefs.end(),this->orderByDetIdThenTime);
341 
342  //loop over sorted hits
343  nchits = hitRefs.size();
344  for(int i=0; i<nchits; ++i) {
345  const int hitidx = std::get<0>(hitRefs[i]);
346  const uint32_t id = std::get<1>(hitRefs[i]);
347 
348  //get the data for this cell, if not available then we skip it
349 
350  if( !validIds_.count(id) ) continue;
351  HGCSimHitDataAccumulator::iterator simHitIt = simHitAccumulator_->emplace(id,HGCCellInfo()).first;
352 
353  if(id==0) continue; // to be ignored at RECO level
354 
355  const float toa = std::get<2>(hitRefs[i]);
356  const PCaloHit &hit=hits->at( hitidx );
357  const float charge = hit.energy()*1e6*keV2fC*getCCE(geom,id,cce_);
358 
359  //distance to the center of the detector
360  const float dist2center( getPositionDistance(geom,id) );
361 
362  //hit time: [time()]=ns [centerDist]=cm [refSpeed_]=cm/ns + delay by 1ns
363  //accumulate in 15 buckets of 25ns (9 pre-samples, 1 in-time, 5 post-samples)
364  const float tof = toa-dist2center/refSpeed_+tofDelay_ ;
365  const int itime= std::floor( tof/bxTime_ ) + 9;
366 
367  //no need to add bx crossing - tof comes already corrected from the mixing module
368  //itime += bxCrossing;
369  //itime += 9;
370 
371  if(itime<0 || itime>14) continue;
372 
373  //check if time index is ok and store energy
374  if(itime >= (int)simHitIt->second.hit_info[0].size() ) continue;
375 
376  (simHitIt->second).hit_info[0][itime] += charge;
377  float accCharge=(simHitIt->second).hit_info[0][itime];
378 
379  //time-of-arrival (check how to be used)
380  if(weightToAbyEnergy) (simHitIt->second).hit_info[1][itime] += charge*tof;
381  else if((simHitIt->second).hit_info[1][itime]==0) {
382  if( accCharge>tdcOnset)
383  {
384  //extrapolate linear using previous simhit if it concerns to the same DetId
385  float fireTDC=tof;
386  if(i>0)
387  {
388  uint32_t prev_id = std::get<1>(hitRefs[i-1]);
389  if(prev_id==id)
390  {
391  float prev_toa = std::get<2>(hitRefs[i-1]);
392  float prev_tof(prev_toa-dist2center/refSpeed_+tofDelay_);
393  //float prev_charge = std::get<3>(hitRefs[i-1]);
394  float deltaQ2TDCOnset = tdcOnset-((simHitIt->second).hit_info[0][itime]-charge);
395  float deltaQ = charge;
396  float deltaT = (tof-prev_tof);
397  fireTDC = deltaT*(deltaQ2TDCOnset/deltaQ)+prev_tof;
398  }
399  }
400 
401  (simHitIt->second).hit_info[1][itime]=fireTDC;
402  }
403  }
404  }
405  hitRefs.clear();
406 }
407 
408 //
410 {
411  //get geometry
413  es.get<CaloGeometryRecord>().get(geom);
414 
415  gHGCal_ = nullptr;
416  gHcal_ = nullptr;
417 
418  if( producesEEDigis() ) gHGCal_ = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(DetId::Forward, HGCEE));
419  if( producesHEfrontDigis() ) gHGCal_ = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(DetId::Forward, HGCHEF));
420  if( producesHEbackDigis() ) gHcal_ = dynamic_cast<const HcalGeometry*>(geom->getSubdetectorGeometry(DetId::Hcal, HcalEndcap));
421 
422  int nadded(0);
423  //valid ID lists
424  if( nullptr != gHGCal_ ) {
425  getValidDetIds( gHGCal_, validIds_ );
426  } else if( nullptr != gHcal_ ) {
427  getValidDetIds( gHcal_, validIds_ );
428  } else {
429  throw cms::Exception("BadConfiguration")
430  << "HGCDigitizer is not producing EE, FH, or BH digis!";
431  }
432 
433  if (verbosity_ > 0)
434  edm::LogInfo("HGCDigitizer")
435  << "Added " << nadded << ":" << validIds_.size()
436  << " detIds without " << hitCollection_
437  << " in first event processed" << std::endl;
438 }
439 
440 //
442 {
443  std::unordered_set<DetId>().swap(validIds_);
444 }
445 
446 //
448 {
449  for( HGCSimHitDataAccumulator::iterator it = simHitAccumulator_->begin(); it!=simHitAccumulator_->end(); it++)
450  {
451  it->second.hit_info[0].fill(0.);
452  it->second.hit_info[1].fill(0.);
453  }
454 }
455 
456 template void HGCDigitizer::accumulate<HcalGeometry>(edm::Handle<edm::PCaloHitContainer> const &hits, int bxCrossing,const HcalGeometry *geom, CLHEP::HepRandomEngine* hre);
457 template void HGCDigitizer::accumulate<HGCalGeometry>(edm::Handle<edm::PCaloHitContainer> const &hits, int bxCrossing,const HGCalGeometry *geom, CLHEP::HepRandomEngine* hre);
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:122
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 HcalTopology & topology() const
Definition: HcalGeometry.h:117
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
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
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:416
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
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.h:108
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
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
HGCDigitizer(const edm::ParameterSet &ps, edm::ConsumesCollector &iC)
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