CMS 3D CMS Logo

GenericSimClusterMapper.cc
Go to the documentation of this file.
10 
13 
14 public:
16  : InitialClusteringStepBase(conf, cc) {
17  _simClusterToken = cc.consumes<SimClusterCollection>(conf.getParameter<edm::InputTag>("simClusterSrc"));
18  }
19  ~GenericSimClusterMapper() override = default;
20  GenericSimClusterMapper(const B2DGT&) = delete;
21  B2DGT& operator=(const B2DGT&) = delete;
22 
23  void updateEvent(const edm::Event&) final;
24 
26  const std::vector<bool>&,
27  const std::vector<bool>&,
29  const HcalPFCuts*) override;
30 
31 private:
34 };
35 
37 
38 #ifdef PFLOW_DEBUG
39 #define LOGVERB(x) edm::LogVerbatim(x)
40 #define LOGWARN(x) edm::LogWarning(x)
41 #define LOGERR(x) edm::LogError(x)
42 #define LOGDRESSED(x) edm::LogInfo(x)
43 #else
44 #define LOGVERB(x) LogTrace(x)
45 #define LOGWARN(x) edm::LogWarning(x)
46 #define LOGERR(x) edm::LogError(x)
47 #define LOGDRESSED(x) LogDebug(x)
48 #endif
49 
51 
53  const std::vector<bool>& rechitMask,
54  const std::vector<bool>& seedable,
56  const HcalPFCuts* hcalCuts) {
58  auto const& hits = *input;
59 
60  // for quick indexing back to hit energy
61  std::unordered_map<uint32_t, size_t> detIdToIndex(hits.size());
62  for (uint32_t i = 0; i < hits.size(); ++i) {
63  detIdToIndex[hits[i].detId()] = i;
64  auto ref = makeRefhit(input, i);
65  }
66 
67  for (const auto& sc : simClusters) {
68  output.emplace_back();
69  reco::PFCluster& back = output.back();
71  double energy = 0.0, highest_energy = 0.0;
72  auto hitsAndFractions = sc.hits_and_fractions();
73  for (const auto& hAndF : hitsAndFractions) {
74  auto itr = detIdToIndex.find(hAndF.first);
75  if (itr == detIdToIndex.end())
76  continue; // hit wasn't saved in reco
77  auto ref = makeRefhit(input, itr->second);
78  const double hit_energy = hAndF.second * ref->energy();
79  energy += hit_energy;
80  back.addRecHitFraction(reco::PFRecHitFraction(ref, hAndF.second));
81  if (hit_energy > highest_energy || highest_energy == 0.0) {
82  highest_energy = hit_energy;
83  seed = ref;
84  }
85  }
86  if (!back.hitsAndFractions().empty()) {
87  back.setSeed(seed->detId());
88  back.setEnergy(energy);
90  } else {
91  back.setSeed(-1);
92  back.setEnergy(0.f);
93  }
94  }
95 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
reco::PFRecHitRef makeRefhit(const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
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
GenericSimClusterMapper(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
void updateEvent(const edm::Event &) final
static std::string const input
Definition: EdmProvDump.cc:50
void setSeed(const DetId &id)
Definition: CaloCluster.h:146
edm::Handle< SimClusterCollection > _simClusterH
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:137
double f[11][100]
edm::EDGetTokenT< SimClusterCollection > _simClusterToken
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)
Definition: output.py:1
std::vector< SimCluster > SimClusterCollection
~GenericSimClusterMapper() override=default
void buildClusters(const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &, const std::vector< bool > &, reco::PFClusterCollection &, const HcalPFCuts *) override