CMS 3D CMS Logo

CaloTruthAccumulator.cc
Go to the documentation of this file.
1 // S Zenz/L Gray, May 2016
2 // Loosely based on TrackingTruthAccumulator (M Grimes)
3 
6 
12 
20 
26 
31 
35 
38 
39 #include <iterator>
40 
42  messageCategory_("CaloTruthAccumulator"),
43  maximumPreviousBunchCrossing_( config.getParameter<unsigned int>("maximumPreviousBunchCrossing") ),
44  maximumSubsequentBunchCrossing_( config.getParameter<unsigned int>("maximumSubsequentBunchCrossing") ),
45  simTrackLabel_( config.getParameter<edm::InputTag>("simTrackCollection") ),
46  simVertexLabel_( config.getParameter<edm::InputTag>("simVertexCollection") ),
47  collectionTags_( ),
48  genParticleLabel_( config.getParameter<edm::InputTag>("genParticleCollection") ),
49  hepMCproductLabel_( config.getParameter<edm::InputTag>("HepMCProductLabel") ),
50  minEnergy_( config.getParameter<double>("MinEnergy") ),
51  maxPseudoRapidity_( config.getParameter<double>("MaxPseudoRapidity") )
52 {
54 
55  mixMod.produces<SimClusterCollection>("MergedCaloTruth");
56  mixMod.produces<CaloParticleCollection>("MergedCaloTruth");
57 
58  iC.consumes<std::vector<SimTrack> >(simTrackLabel_);
59  iC.consumes<std::vector<SimVertex> >(simVertexLabel_);
60  iC.consumes<std::vector<reco::GenParticle> >(genParticleLabel_);
61  iC.consumes<std::vector<int> >(genParticleLabel_);
62  iC.consumes<std::vector<int> >(hepMCproductLabel_);
63 
64  // Fill the collection tags
65  const edm::ParameterSet& simHitCollectionConfig=config.getParameterSet("simHitCollections");
66  std::vector<std::string> parameterNames=simHitCollectionConfig.getParameterNames();
67 
68  for( const auto& parameterName : parameterNames )
69  {
70  std::vector<edm::InputTag> tags=simHitCollectionConfig.getParameter<std::vector<edm::InputTag> >(parameterName);
71  collectionTags_.insert(collectionTags_.end(), tags.begin(), tags.end());
72  }
73 
74  for( const auto& collectionTag : collectionTags_ ) {
75  iC.consumes<std::vector<PCaloHit> >(collectionTag);
76  }
77 }
78 
81  iSetup.get<CaloGeometryRecord>().get(geom);
82  const HGCalGeometry *eegeom, *fhgeom;
83  const HcalGeometry *bhgeom;
84 
85  eegeom = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(DetId::Forward,HGCEE));
86  fhgeom = dynamic_cast<const HGCalGeometry*>(geom->getSubdetectorGeometry(DetId::Forward,HGCHEF));
87  bhgeom = dynamic_cast<const HcalGeometry*>(geom->getSubdetectorGeometry(DetId::Hcal,HcalEndcap));
88 
89  hgtopo_[0] = &(eegeom->topology());
90  hgtopo_[1] = &(fhgeom->topology());
91 
92  for( unsigned i = 0; i < 2; ++i ) {
93  hgddd_[i] = &(hgtopo_[i]->dddConstants());
94  }
95 
96  hcddd_ = bhgeom->topology().dddConstants();
97 
98  caloStartZ = hgddd_[0]->waferZ(1,false)*10.0; // get the location of the first plane of silicon, put in mm
99 }
100 
102 {
105 
106  m_detIdToCluster.clear();
107  m_detIdToTotalSimEnergy.clear();
108 }
109 
114 
116  // Call the templated version that does the same for both signal and pileup events
117 
119  event.getByLabel(hepMCproductLabel_, hepmc);
120 
121  edm::LogInfo(messageCategory_) << " CaloTruthAccumulator::accumulate (signal)";
122  accumulateEvent( event, setup, hepmc );
123 }
124 
125 void CaloTruthAccumulator::accumulate( PileUpEventPrincipal const& event, edm::EventSetup const& setup, edm::StreamID const& ) {
126  // If this bunch crossing is outside the user configured limit, don't do anything.
127  if( event.bunchCrossing()>=-static_cast<int>(maximumPreviousBunchCrossing_) && event.bunchCrossing()<=static_cast<int>(maximumSubsequentBunchCrossing_) )
128  {
129  //edm::LogInfo(messageCategory_) << "Analysing pileup event for bunch crossing " << event.bunchCrossing();
130 
131  //simply create empty handle as we do not have a HepMCProduct in PU anyway
133  edm::LogInfo(messageCategory_) << " CaloTruthAccumulator::accumulate (pileup) bunchCrossing=" << event.bunchCrossing();
134  accumulateEvent( event, setup, hepmc );
135  }
136  else edm::LogInfo(messageCategory_) << "Skipping pileup event for bunch crossing " << event.bunchCrossing();
137 }
138 
140  edm::LogInfo("CaloTruthAccumulator") << "Adding " << output_.pSimClusters->size()
141  << " SimParticles and " << output_.pCaloParticles->size()
142  << " CaloParticles to the event.";
143 
144  // now we need to normalize the hits and energies into hits and fractions
145  // (since we have looped over all pileup events)
146 
147  for( auto& sc : *(output_.pSimClusters) ) {
148  auto hitsAndEnergies = sc.hits_and_fractions();
150  for( auto& hAndE : hitsAndEnergies ) {
151  const float totalenergy = m_detIdToTotalSimEnergy[hAndE.first];
152  float fraction = 0.;
153  if(totalenergy>0) fraction = hAndE.second/totalenergy;
154  else edm::LogWarning("CaloTruthAccumulator") << "TotalSimEnergy for hit " << hAndE.first << " is 0! The fraction for this hit cannot be computed.";
155  sc.addRecHitAndFraction(hAndE.first,fraction);
156  }
157  }
158 
159  // save the SimCluster orphan handle so we can fill the calo particles
160  auto scHandle = event.put( std::move(output_.pSimClusters), "MergedCaloTruth" );
161 
162  // now fill the calo particles
163  for( unsigned i = 0; i < output_.pCaloParticles->size(); ++i ) {
164  auto& cp = (*output_.pCaloParticles)[i];
165  for( unsigned j = m_caloParticles.sc_start_[i]; j < m_caloParticles.sc_stop_[i]; ++j ) {
166  edm::Ref<SimClusterCollection> ref(scHandle,j);
167  cp.addSimCluster(ref);
168  }
169  }
170 
171  event.put( std::move(output_.pCaloParticles), "MergedCaloTruth" );
172 
174 
175  std::unordered_map<Index_t,float>().swap(m_detIdToTotalSimEnergy);
176 
177  std::unordered_map<Barcode_t,Index_t>().swap(m_genParticleBarcodeToIndex);
178  std::unordered_map<Barcode_t,Index_t>().swap(m_simTrackBarcodeToIndex);
179  std::unordered_map<Barcode_t,Index_t>().swap(m_genBarcodeToSimTrackIndex);
180  std::unordered_map<Barcode_t,Index_t>().swap(m_simVertexBarcodeToIndex);
181 
182  std::unordered_multimap<Index_t,Index_t>().swap(m_detIdToCluster);
183  std::unordered_multimap<Barcode_t,Index_t>().swap(m_simHitBarcodeToIndex);
184  std::unordered_multimap<Barcode_t,Barcode_t>().swap(m_simVertexBarcodeToSimTrackBarcode);
185  std::unordered_map<Barcode_t,Barcode_t>().swap(m_simTrackBarcodeToSimVertexParentBarcode);
186  std::unordered_multimap<Barcode_t,Index_t>().swap(m_simTrackToSimVertex);
187  std::unordered_multimap<Barcode_t,Index_t>().swap(m_simVertexToSimTrackParent);
188  std::vector<Barcode_t>().swap(m_simVertexBarcodes);
189  std::unordered_map<Index_t,float>().swap(m_detIdToTotalSimEnergy);
190 }
191 
192 template<class T>
194  const edm::EventSetup& setup,
195  const edm::Handle< edm::HepMCProduct >& hepMCproduct) {
196  //
197  // Get the collections
198  //
199 
201  edm::Handle< std::vector<int> > hGenParticleIndices;
202 
203  event.getByLabel( simTrackLabel_, hSimTracks );
204  event.getByLabel( simVertexLabel_, hSimVertices );
205 
206  event.getByLabel( genParticleLabel_, hGenParticles );
207  event.getByLabel( genParticleLabel_, hGenParticleIndices );
208 
209  std::vector<std::pair<DetId,const PCaloHit*> > simHitPointers;
210  fillSimHits( simHitPointers, event, setup );
211 
212  // Clear maps from previous event fill them for this one
213  m_simHitBarcodeToIndex.clear();
215  for (unsigned int i = 0 ; i < simHitPointers.size(); ++i) {
216  m_simHitBarcodeToIndex.emplace(simHitPointers[i].second->geantTrackId(),i);
217  }
219  if( hGenParticles.isValid() && hGenParticleIndices.isValid() ) {
220  for (unsigned int i = 0 ; i < hGenParticles->size() ; ++i) {
221  m_genParticleBarcodeToIndex.emplace(hGenParticleIndices->at(i),i);
222  }
223  }
227  m_simTrackBarcodeToIndex.clear();
228  for (unsigned int i = 0 ; i < hSimTracks->size() ; i++) {
229  if( !hSimTracks->at(i).noGenpart() ) {
230  m_genBarcodeToSimTrackIndex.emplace(hSimTracks->at(i).genpartIndex(), i);
231  }
232  if( !hSimTracks->at(i).noVertex() ) {
233  m_simVertexBarcodeToSimTrackBarcode.emplace(hSimTracks->at(i).vertIndex(), hSimTracks->at(i).trackId());
234  m_simTrackBarcodeToSimVertexParentBarcode.emplace(hSimTracks->at(i).trackId(), hSimTracks->at(i).vertIndex());
235  }
236  m_simTrackBarcodeToIndex.emplace(hSimTracks->at(i).trackId(), i);
237  }
238  m_simVertexBarcodes.clear();
240  m_simTrackToSimVertex.clear();
242  for (unsigned int i = 0 ; i < hSimVertices->size() ; i++) {
243  m_simVertexBarcodes.push_back(i);
244  m_simVertexBarcodeToIndex.emplace(hSimVertices->at(i).vertexId(), i);
245  if (!hSimVertices->at(i).noParent()) {
246  m_simTrackToSimVertex.emplace(hSimVertices->at(i).parentIndex(), i);
247  m_simVertexToSimTrackParent.emplace( hSimVertices->at(i).vertexId(), hSimVertices->at(i).parentIndex() );
248  }
249  }
250 
251  std::vector<Index_t> tracksToBecomeClustersInitial;
252  std::vector<Barcode_t> descendantTracks;
253  SimClusterCollection simClustersForGenParts;
254  std::vector<std::unique_ptr<SimHitInfoPerSimTrack_t> > hitInfoList;
255  std::vector<std::vector<uint32_t> > simClusterPrimitives;
256  std::unordered_multimap<Index_t,Index_t> genPartsToSimClusters;
257  const auto& simTracks = *hSimTracks;
258  // loop over
259  for (unsigned int i = 0 ; i < simTracks.size() ; ++i) {
260  if ( simTracks[i].momentum().E() < minEnergy_ || std::abs(simTracks[i].momentum().Eta()) >= maxPseudoRapidity_ ) continue;
261  if ( simTracks[i].noGenpart() ) continue;
262  auto temp = CaloTruthAccumulator::descendantSimClusters( simTracks[i].trackId(),simHitPointers );
263  if( !temp.empty() ) {
264  output_.pCaloParticles->emplace_back(simTracks[i]);
265  m_caloParticles.sc_start_.push_back(output_.pSimClusters->size());
266  auto mbegin = std::make_move_iterator(temp.begin());
267  auto mend = std::make_move_iterator(temp.end());
268  output_.pSimClusters->insert(output_.pSimClusters->end(), mbegin, mend);
269  m_caloParticles.sc_stop_.push_back(output_.pSimClusters->size());
270  }
271  }
272 }
273 
274 std::vector<Barcode_t> CaloTruthAccumulator::descendantTrackBarcodes( Barcode_t barcode ) {
275  std::vector<Barcode_t> result;
276  if (m_simTrackToSimVertex.count(barcode)) {
277  auto vertex_range = m_simTrackToSimVertex.equal_range(barcode);
278  for ( auto vertex_iter = vertex_range.first ; vertex_iter != vertex_range.second ; vertex_iter++ ) {
279  Index_t decayVertexIndex = vertex_iter->second;
280  Barcode_t decayVertexBarcode = m_simVertexBarcodes[decayVertexIndex];
281  auto track_range = m_simVertexBarcodeToSimTrackBarcode.equal_range( decayVertexBarcode );
282  for ( auto track_iter = track_range.first ; track_iter != track_range.second ; track_iter++ ) {
283  result.push_back( track_iter->second );
284  std::vector<Barcode_t> daughter_result = CaloTruthAccumulator::descendantTrackBarcodes( track_iter->second );
285  result.insert(result.end(),daughter_result.begin(),daughter_result.end());
286  }
287  }
288  }
289  return result;
290 }
291 
292 SimClusterCollection CaloTruthAccumulator::descendantSimClusters( Barcode_t barcode, const std::vector<std::pair<DetId,const PCaloHit*> >& hits ) {
294  const auto& simTracks = *hSimTracks;
295  if ( CaloTruthAccumulator::consideredBarcode( barcode ) ) {
296  LogDebug("CaloParticles") << "SCZ DEBUG Ignoring descendantSimClusters call because this particle is already marked used: " << barcode << std::endl;
297  //return result;
298  }
299 
300  std::unique_ptr<SimHitInfoPerSimTrack_t> hit_info = std::move(CaloTruthAccumulator::attachedSimHitInfo(barcode,hits, true, false, false));
301  //std::unique_ptr<SimHitInfoPerSimTrack_t> inclusive_hit_info = std::move(CaloTruthAccumulator::allAttachedSimHitInfo(barcode,hits, false) );
302 
303  const auto& simTrack = simTracks[m_simTrackBarcodeToIndex[barcode]];
305  const auto& vtx = hSimVertices->at(m_simVertexBarcodeToIndex[vtxBarcode]);
306  const bool isInCalo = (std::abs(vtx.position().z()) > caloStartZ - 30.0); // add a buffer region in front of the calo face
307 
308  if (!hit_info->empty()) {
309  // define the sim cluster starting from the earliest particle that has hits in the calorimeter
310  // grab everything that descends from it
311  std::unique_ptr<SimHitInfoPerSimTrack_t> marked_hit_info;
312 
313 
314  if( isInCalo ) {
315  marked_hit_info = std::move( CaloTruthAccumulator::allAttachedSimHitInfo(barcode,hits,true) );
316  } else {
317  marked_hit_info = std::move( CaloTruthAccumulator::attachedSimHitInfo(barcode,hits,true,false,true) );
318  }
319 
320  if( !marked_hit_info->empty() ) {
321  result.emplace_back(simTrack);
322  auto& simcluster = result.back();
323 
324  std::unordered_map<uint32_t,float> acc_energy;
325 
326  for( const auto& hit_and_energy : *marked_hit_info ) {
327  const uint32_t id = hit_and_energy.first.rawId();
328  if( acc_energy.count(id) ) acc_energy[id] += hit_and_energy.second;
329  else acc_energy[id] = hit_and_energy.second;
330  }
331 
332  for( const auto& hit_and_energy : acc_energy ) {
333  simcluster.addRecHitAndFraction(hit_and_energy.first,hit_and_energy.second);
334  }
335  }
336  }
337 
338  if ( m_simTrackToSimVertex.count(barcode) ) {
339  auto vertex_range = m_simTrackToSimVertex.equal_range(barcode);
340  for ( auto vertex_iter = vertex_range.first ; vertex_iter != vertex_range.second ; vertex_iter++ ) {
341  Index_t decayVertexIndex = vertex_iter->second;
342  Barcode_t decayVertexBarcode = m_simVertexBarcodes[decayVertexIndex];
343  auto track_range = m_simVertexBarcodeToSimTrackBarcode.equal_range( decayVertexBarcode );
344  for ( auto track_iter = track_range.first ; track_iter != track_range.second ; track_iter++ ) {
345 
346  auto daughter_result = CaloTruthAccumulator::descendantSimClusters(track_iter->second,hits);
347  result.insert(result.end(),daughter_result.begin(),daughter_result.end());
348  }
349  }
350  }
351 
352  return result;
353 }
354 
355 std::unique_ptr<SimHitInfoPerSimTrack_t> CaloTruthAccumulator::attachedSimHitInfo( Barcode_t barcode , const std::vector<std::pair<DetId,const PCaloHit*> >& hits,
356  bool includeOwn , bool includeOther, bool markUsed ) {
357  const auto& simTracks = *hSimTracks;
358  std::unique_ptr<SimHitInfoPerSimTrack_t> result(new SimHitInfoPerSimTrack_t);
359 
360  const auto& simTrack = simTracks[m_simTrackBarcodeToIndex[barcode]];
361 
362  if ( markUsed ) {
363  if ( CaloTruthAccumulator::consideredBarcode( barcode ) ) {
364  return result;
365  }
367  }
368  if (includeOwn) {
369  auto range = m_simHitBarcodeToIndex.equal_range( barcode );
370  unsigned n = 0;
371  for ( auto iter = range.first ; iter != range.second ; iter++ ) {
372  const auto& the_hit = hits[iter->second];
373  result->emplace_back(the_hit.first,the_hit.second->energy());
374  ++n;
375  }
376  }
377 
378  // need to sim to the next sim track if we explicitly ask or
379  // if we are in the calorimeter next (no interaction)
380  // or if this is a continuation of the same particle
381  if (m_simTrackToSimVertex.count(barcode)) {
382  auto vertex_range = m_simTrackToSimVertex.equal_range(barcode);
383  for ( auto vertex_iter = vertex_range.first ; vertex_iter != vertex_range.second ; vertex_iter++ ) {
384  Index_t decayVertexIndex = vertex_iter->second;
385  const auto& nextVtx = (*hSimVertices)[decayVertexIndex];
386  const bool nextInCalo = (std::abs(nextVtx.position().z()) > caloStartZ*0.1 - 20.0); // add a buffer region in front of the calo face
387 
388  Barcode_t decayVertexBarcode = m_simVertexBarcodes[decayVertexIndex];
389  auto track_range = m_simVertexBarcodeToSimTrackBarcode.equal_range( decayVertexBarcode );
390  for ( auto track_iter = track_range.first ; track_iter != track_range.second ; track_iter++ ) {
391  if( !barcodeLogicWarningAlready_ && track_iter->second < barcode ) {
393  edm::LogWarning(messageCategory_) << " Daughter particle has a lower barcode than parent. This may screw up the logic!" << std::endl;
394  }
395  const auto& daughter = simTracks[m_simTrackBarcodeToIndex[track_iter->second]];
396 
397  if( includeOther || nextInCalo ) {
398  std::unique_ptr<SimHitInfoPerSimTrack_t> daughter_result = std::move(CaloTruthAccumulator::allAttachedSimHitInfo(track_iter->second,hits,markUsed));
399  result->insert(result->end(),daughter_result->begin(),daughter_result->end());
400  } else if ( daughter.type() == simTrack.type() ) {
401  std::unique_ptr<SimHitInfoPerSimTrack_t> daughter_result = std::move(CaloTruthAccumulator::attachedSimHitInfo(track_iter->second,hits,includeOwn, includeOther, markUsed));
402  result->insert(result->end(),daughter_result->begin(),daughter_result->end());
403  }
404  }
405  }
406  }
407  return result;
408 }
409 
410 std::unique_ptr<SimHitInfoPerSimTrack_t>
412  const std::vector<std::pair<DetId,const PCaloHit*> >& hits,
413  bool markUsed) {
414  return CaloTruthAccumulator::attachedSimHitInfo(barcode,hits,false,true,markUsed);
415 }
416 
417 std::unique_ptr<SimHitInfoPerSimTrack_t>
419  const std::vector<std::pair<DetId,const PCaloHit*> >& hits,
420  bool markUsed ) {
421  return CaloTruthAccumulator::attachedSimHitInfo(barcode,hits,true,true,markUsed);
422 }
423 
424 template<class T> void CaloTruthAccumulator::fillSimHits( std::vector<std::pair<DetId, const PCaloHit*> >& returnValue, const T& event, const edm::EventSetup& setup ) {
425  // loop over the collections
426  for( const auto& collectionTag : collectionTags_ ) {
428  const bool isHcal = ( collectionTag.instance().find("HcalHits") != std::string::npos );
429  event.getByLabel( collectionTag, hSimHits );
430  for( const auto& simHit : *hSimHits ) {
431  DetId id(0);
432  const uint32_t simId = simHit.id();
433  if( isHcal ) {
435  if(hid.subdet()==HcalEndcap) id = hid;
436  } else {
437  int subdet, layer, cell, sec, subsec, zp;
438  HGCalTestNumbering::unpackHexagonIndex(simId, subdet, zp, layer, sec, subsec, cell);
439  const HGCalDDDConstants* ddd = hgddd_[subdet-3];
440  std::pair<int,int> recoLayerCell = ddd->simToReco(cell,layer,sec,
441  hgtopo_[subdet-3]->detectorType());
442  cell = recoLayerCell.first;
443  layer = recoLayerCell.second;
444  // skip simhits with bad barcodes or non-existant layers
445  if( layer == -1 || simHit.geantTrackId() == 0 ) continue;
446  id = HGCalDetId((ForwardSubdetector)subdet,zp,layer,subsec,sec,cell);
447  }
448 
449  if( DetId(0) == id ) continue;
450 
451  uint32_t detId = id.rawId();
452  returnValue.emplace_back(id, &simHit);
453 
454  if( m_detIdToTotalSimEnergy.count(detId) ) m_detIdToTotalSimEnergy[detId] += simHit.energy();
455  else m_detIdToTotalSimEnergy[detId] = simHit.energy();
456  }
457  } // end of loop over InputTags
458 }
459 
460 // Register with the framework
#define LogDebug(id)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:45
void setConsideredBarcode(Barcode_t barcode)
calo_particles m_caloParticles
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:161
std::unordered_multimap< Index_t, Index_t > m_detIdToCluster
edm::Handle< std::vector< SimVertex > > hSimVertices
std::unordered_multimap< Barcode_t, Index_t > m_simTrackToSimVertex
void initializeEvent(const edm::Event &event, const edm::EventSetup &setup) override
std::vector< Barcode_t > descendantTrackBarcodes(Barcode_t barcode)
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
std::unordered_multimap< Barcode_t, Index_t > m_simHitBarcodeToIndex
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:210
const unsigned int maximumSubsequentBunchCrossing_
void beginLuminosityBlock(edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) override
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
std::unordered_map< Barcode_t, Index_t > m_genParticleBarcodeToIndex
std::unordered_map< Barcode_t, Index_t > m_genBarcodeToSimTrackIndex
std::vector< edm::InputTag > collectionTags_
Definition: config.py:1
std::unique_ptr< SimHitInfoPerSimTrack_t > attachedSimHitInfo(Barcode_t st, const std::vector< std::pair< DetId, const PCaloHit * > > &hits, bool includeOwn=true, bool includeOther=false, bool markUsed=false)
ForwardSubdetector
OutputCollections output_
void addRecHitAndFraction(uint32_t hit, float fraction)
add rechit with fraction
Definition: SimCluster.h:204
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
const HcalTopology & topology() const
Definition: HcalGeometry.h:120
U second(std::pair< T, U > const &p)
edm::InputTag hepMCproductLabel_
Needed to add HepMC::GenVertex to SimVertex.
const std::string messageCategory_
The message category used to send messages to MessageLogger.
std::pair< int, int > simToReco(int cell, int layer, int mod, bool half) const
void accumulate(const edm::Event &event, const edm::EventSetup &setup) override
std::unique_ptr< SimHitInfoPerSimTrack_t > allAttachedSimHitInfo(Barcode_t st, const std::vector< std::pair< DetId, const PCaloHit * > > &hits, bool markUsed=false)
simTrack
per collection params
const HGCalTopology & topology() const
Definition: HGCalGeometry.h:96
void fillSimHits(std::vector< std::pair< DetId, const PCaloHit * > > &returnValue, const T &event, const edm::EventSetup &setup)
Fills the supplied vector with pointers to the SimHits, checking for bad modules if required...
void clearHitsAndFractions()
clear the hits and fractions list
Definition: SimCluster.h:219
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const edm::InputTag simVertexLabel_
void addSimCluster(const SimClusterRef &ref)
Definition: CaloParticle.h:62
std::vector< std::string > getParameterNames() const
bool isValid() const
Definition: HandleBase.h:74
bool consideredBarcode(Barcode_t barcode)
std::unordered_map< Barcode_t, Barcode_t > m_simTrackBarcodeToSimVertexParentBarcode
const HGCalDDDConstants * hgddd_[2]
std::vector< Barcode_t > m_simVertexBarcodes
double waferZ(int layer, bool reco) const
const edm::InputTag simTrackLabel_
If bremsstrahlung merging, whether to also add the unmerged collection to the event or not...
CaloTruthAccumulator(const edm::ParameterSet &config, edm::stream::EDProducerBase &mixMod, edm::ConsumesCollector &iC)
std::unordered_map< Barcode_t, Index_t > m_simTrackBarcodeToIndex
edm::Handle< std::vector< SimTrack > > hSimTracks
std::unordered_map< Barcode_t, Index_t > m_simVertexBarcodeToIndex
Definition: DetId.h:18
void accumulateEvent(const T &event, const edm::EventSetup &setup, const edm::Handle< edm::HepMCProduct > &hepMCproduct)
Both forms of accumulate() delegate to this templated method.
std::unique_ptr< CaloParticleCollection > pCaloParticles
const HGCalDDDConstants & dddConstants() const
ParameterSet const & getParameterSet(std::string const &) const
const T & get() const
Definition: EventSetup.h:55
unsigned Index_t
std::unique_ptr< SimHitInfoPerSimTrack_t > descendantOnlySimHitInfo(Barcode_t st, const std::vector< std::pair< DetId, const PCaloHit * > > &hits, bool markUsed=false)
const HcalDDDRecConstants * hcddd_
void finalizeEvent(edm::Event &event, const edm::EventSetup &setup) override
std::vector< SimHitInfoPerRecoDetId_t > SimHitInfoPerSimTrack_t
int Barcode_t
std::set< Barcode_t > m_simTracksConsideredForSimClusters
HLT enums.
#define DEFINE_DIGI_ACCUMULATOR(type)
std::vector< CaloParticle > CaloParticleCollection
SimClusterCollection descendantSimClusters(Barcode_t barcode, const std::vector< std::pair< DetId, const PCaloHit * > > &hits)
std::unique_ptr< SimClusterCollection > pSimClusters
const HGCalTopology * hgtopo_[2]
DetId relabel(const uint32_t testId) const
std::unordered_multimap< Barcode_t, Index_t > m_simVertexToSimTrackParent
std::vector< SimCluster > SimClusterCollection
Definition: SimClusterFwd.h:8
std::unordered_map< Index_t, float > m_detIdToTotalSimEnergy
long double T
std::unordered_multimap< Barcode_t, Barcode_t > m_simVertexBarcodeToSimTrackBarcode
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
Definition: event.py:1
const unsigned int maximumPreviousBunchCrossing_
maximum distance for HepMC::GenVertex to be added to SimVertex