CMS 3D CMS Logo

PreMixingCaloParticleWorker.cc
Go to the documentation of this file.
7 
13 
16 
18 public:
20  ~PreMixingCaloParticleWorker() override = default;
21 
22  void initializeEvent(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
23  void addSignals(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
24  void addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &iSetup) override;
25  void put(edm::Event &iEvent,
26  edm::EventSetup const &iSetup,
27  std::vector<PileupSummaryInfo> const &ps,
28  int bunchSpacing) override;
29 
30 private:
31  using EnergyMap = std::vector<std::pair<unsigned, float>>;
32 
33  void add(const SimClusterCollection &clusters, const CaloParticleCollection &particles, const EnergyMap &energyMap);
34 
38 
41 
42  std::unordered_map<unsigned, float> totalEnergy_;
43 
44  std::unique_ptr<SimClusterCollection> newClusters_;
45  std::unique_ptr<CaloParticleCollection> newParticles_;
47 };
48 
52  : sigClusterToken_(iC.consumes<SimClusterCollection>(ps.getParameter<edm::InputTag>("labelSig"))),
53  sigParticleToken_(iC.consumes<CaloParticleCollection>(ps.getParameter<edm::InputTag>("labelSig"))),
54  sigEnergyToken_(iC.consumes<EnergyMap>(ps.getParameter<edm::InputTag>("labelSig"))),
55  particlePileInputTag_(ps.getParameter<edm::InputTag>("pileInputTag")),
56  particleCollectionDM_(ps.getParameter<std::string>("collectionDM")) {
59 }
60 
62  newClusters_ = std::make_unique<SimClusterCollection>();
63  newParticles_ = std::make_unique<CaloParticleCollection>();
64 
65  // need RefProds in order to re-key the CaloParticle->SimCluster refs
66  // TODO: try to remove const_cast, requires making Event non-const in
67  // BMixingModule::initializeEvent
68  clusterRef_ = const_cast<edm::Event &>(iEvent).getRefBeforePut<SimClusterCollection>(particleCollectionDM_);
69 }
70 
73  iEvent.getByToken(sigClusterToken_, clusters);
74 
76  iEvent.getByToken(sigParticleToken_, particles);
77 
79  iEvent.getByToken(sigEnergyToken_, energy);
80 
81  if (clusters.isValid() && particles.isValid() && energy.isValid()) {
82  add(*clusters, *particles, *energy);
83  }
84 }
85 
88  pep.getByLabel(particlePileInputTag_, clusters);
89 
91  pep.getByLabel(particlePileInputTag_, particles);
92 
94  pep.getByLabel(particlePileInputTag_, energy);
95 
96  if (clusters.isValid() && particles.isValid() && energy.isValid()) {
97  add(*clusters, *particles, *energy);
98  }
99 }
100 
103  const EnergyMap &energy) {
104  const size_t startingIndex = newClusters_->size();
105 
106  // Copy SimClusters
107  newClusters_->reserve(newClusters_->size() + clusters.size());
108  std::copy(clusters.begin(), clusters.end(), std::back_inserter(*newClusters_));
109 
110  // Copy CaloParticles
111  newParticles_->reserve(newParticles_->size() + particles.size());
112  for (const auto &p : particles) {
113  newParticles_->push_back(p);
114  auto &particle = newParticles_->back();
115 
116  // re-key the refs to SimClusters
117  particle.clearSimClusters();
118  for (const auto &ref : p.simClusters()) {
119  particle.addSimCluster(SimClusterRef(clusterRef_, startingIndex + ref.index()));
120  }
121  }
122 
123  // Add energies
124  for (const auto elem : energy) {
125  totalEnergy_[elem.first] += elem.second;
126  }
127 }
128 
130  edm::EventSetup const &iSetup,
131  std::vector<PileupSummaryInfo> const &ps,
132  int bunchSpacing) {
133  for (auto &sc : *newClusters_) {
134  auto hitsAndEnergies = sc.hits_and_fractions();
136  for (auto &hAndE : hitsAndEnergies) {
137  const float totalenergy = totalEnergy_[hAndE.first];
138  float fraction = 0.;
139  if (totalenergy > 0)
140  fraction = hAndE.second / totalenergy;
141  else
142  edm::LogWarning("PreMixingParticleWorker")
143  << "TotalSimEnergy for hit " << hAndE.first << " is 0! The fraction for this hit cannot be computed.";
144  sc.addRecHitAndFraction(hAndE.first, fraction);
145  }
146  }
147 
148  // clear memory
149  std::unordered_map<unsigned, float>{}.swap(totalEnergy_);
150 
151  iEvent.put(std::move(newClusters_), particleCollectionDM_);
153 }
154 
BranchAliasSetterT< ProductType > produces()
declare what type of product will make and with which optional label
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
void addSignals(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
def copy(args, dbName)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::Ref< SimClusterCollection > SimClusterRef
Definition: SimClusterFwd.h:10
void addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &iSetup) override
edm::EDGetTokenT< CaloParticleCollection > sigParticleToken_
edm::EDGetTokenT< EnergyMap > sigEnergyToken_
void put(edm::Event &iEvent, edm::EventSetup const &iSetup, std::vector< PileupSummaryInfo > const &ps, int bunchSpacing) override
void initializeEvent(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
void addRecHitAndFraction(uint32_t hit, float fraction)
add rechit with fraction
Definition: SimCluster.h:175
int iEvent
Definition: GenABIO.cc:224
~PreMixingCaloParticleWorker() override=default
std::unordered_map< unsigned, float > totalEnergy_
void clearHitsAndFractions()
clear the hits and fractions list
Definition: SimCluster.h:190
std::unique_ptr< SimClusterCollection > newClusters_
edm::EDGetTokenT< SimClusterCollection > sigClusterToken_
bool isValid() const
Definition: HandleBase.h:74
PreMixingCaloParticleWorker(const edm::ParameterSet &ps, edm::ProducerBase &producer, edm::ConsumesCollector &&iC)
def elem(elemtype, innerHTML='', html_class='', kwargs)
Definition: HTMLExport.py:19
void add(const SimClusterCollection &clusters, const CaloParticleCollection &particles, const EnergyMap &energyMap)
HLT enums.
std::vector< CaloParticle > CaloParticleCollection
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
std::vector< std::pair< unsigned, float >> EnergyMap
std::vector< SimCluster > SimClusterCollection
Definition: SimClusterFwd.h:8
#define DEFINE_PREMIXING_WORKER(TYPE)
std::unique_ptr< CaloParticleCollection > newParticles_
def move(src, dest)
Definition: eostools.py:511
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:181