CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GenericSimClusterMapper.cc
Go to the documentation of this file.
8 
11 
12 public:
14  : InitialClusteringStepBase(conf, cc) {
16  }
17  ~GenericSimClusterMapper() override = default;
18  GenericSimClusterMapper(const B2DGT&) = delete;
19  B2DGT& operator=(const B2DGT&) = delete;
20 
21  void updateEvent(const edm::Event&) final;
22 
24  const std::vector<bool>&,
25  const std::vector<bool>&,
27 
28 private:
31 };
32 
34 
35 #ifdef PFLOW_DEBUG
36 #define LOGVERB(x) edm::LogVerbatim(x)
37 #define LOGWARN(x) edm::LogWarning(x)
38 #define LOGERR(x) edm::LogError(x)
39 #define LOGDRESSED(x) edm::LogInfo(x)
40 #else
41 #define LOGVERB(x) LogTrace(x)
42 #define LOGWARN(x) edm::LogWarning(x)
43 #define LOGERR(x) edm::LogError(x)
44 #define LOGDRESSED(x) LogDebug(x)
45 #endif
46 
47 void GenericSimClusterMapper::updateEvent(const edm::Event& ev) { ev.getByToken(_simClusterToken, _simClusterH); }
48 
50  const std::vector<bool>& rechitMask,
51  const std::vector<bool>& seedable,
53  const SimClusterCollection& simClusters = *_simClusterH;
54  auto const& hits = *input;
55 
56  // for quick indexing back to hit energy
57  std::unordered_map<uint32_t, size_t> detIdToIndex(hits.size());
58  for (uint32_t i = 0; i < hits.size(); ++i) {
59  detIdToIndex[hits[i].detId()] = i;
60  auto ref = makeRefhit(input, i);
61  }
62 
63  for (const auto& sc : simClusters) {
64  output.emplace_back();
65  reco::PFCluster& back = output.back();
67  double energy = 0.0, highest_energy = 0.0;
68  auto hitsAndFractions = sc.hits_and_fractions();
69  for (const auto& hAndF : hitsAndFractions) {
70  auto itr = detIdToIndex.find(hAndF.first);
71  if (itr == detIdToIndex.end())
72  continue; // hit wasn't saved in reco
73  auto ref = makeRefhit(input, itr->second);
74  const double hit_energy = hAndF.second * ref->energy();
75  energy += hit_energy;
76  back.addRecHitFraction(reco::PFRecHitFraction(ref, hAndF.second));
77  if (hit_energy > highest_energy || highest_energy == 0.0) {
78  highest_energy = hit_energy;
79  seed = ref;
80  }
81  }
82  if (!back.hitsAndFractions().empty()) {
83  back.setSeed(seed->detId());
84  back.setEnergy(energy);
85  back.setCorrectedEnergy(energy);
86  } else {
87  back.setSeed(-1);
88  back.setEnergy(0.f);
89  }
90  }
91 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
GenericSimClusterMapper B2DGT
B2DGT & operator=(const B2DGT &)=delete
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
void setEnergy(double energy)
Definition: CaloCluster.h:136
std::vector< PFRecHit > PFRecHitCollection
collection of PFRecHit objects
Definition: PFRecHitFwd.h:9
GenericSimClusterMapper(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
void updateEvent(const edm::Event &) final
static std::string const input
Definition: EdmProvDump.cc:47
void setSeed(const DetId &id)
Definition: CaloCluster.h:146
reco::PFRecHitRef makeRefhit(const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
edm::Handle< SimClusterCollection > _simClusterH
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
void buildClusters(const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &, const std::vector< bool > &, reco::PFClusterCollection &) override
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:137
#define override(base_class)
edm::EDGetTokenT< SimClusterCollection > _simClusterToken
std::vector< l1t::PFCluster > PFClusterCollection
Definition: PFCluster.h:73
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< SimCluster > SimClusterCollection
Definition: SimClusterFwd.h:8
~GenericSimClusterMapper() override=default