CMS 3D CMS Logo

GenericSimClusterMapper.cc
Go to the documentation of this file.
3 
5 
8 
9 #ifdef PFLOW_DEBUG
10 #define LOGVERB(x) edm::LogVerbatim(x)
11 #define LOGWARN(x) edm::LogWarning(x)
12 #define LOGERR(x) edm::LogError(x)
13 #define LOGDRESSED(x) edm::LogInfo(x)
14 #else
15 #define LOGVERB(x) LogTrace(x)
16 #define LOGWARN(x) edm::LogWarning(x)
17 #define LOGERR(x) edm::LogError(x)
18 #define LOGDRESSED(x) LogDebug(x)
19 #endif
20 
22 
24  const std::vector<bool>& rechitMask,
25  const std::vector<bool>& seedable,
27  const SimClusterCollection& simClusters = *_simClusterH;
28  auto const& hits = *input;
29 
30  // for quick indexing back to hit energy
31  std::unordered_map<uint32_t, size_t> detIdToIndex(hits.size());
32  for (uint32_t i = 0; i < hits.size(); ++i) {
33  detIdToIndex[hits[i].detId()] = i;
34  auto ref = makeRefhit(input, i);
35  }
36 
37  for (const auto& sc : simClusters) {
38  output.emplace_back();
39  reco::PFCluster& back = output.back();
41  double energy = 0.0, highest_energy = 0.0;
42  auto hitsAndFractions = sc.hits_and_fractions();
43  for (const auto& hAndF : hitsAndFractions) {
44  auto itr = detIdToIndex.find(hAndF.first);
45  if (itr == detIdToIndex.end())
46  continue; // hit wasn't saved in reco
47  auto ref = makeRefhit(input, itr->second);
48  const double hit_energy = hAndF.second * ref->energy();
49  energy += hit_energy;
50  back.addRecHitFraction(reco::PFRecHitFraction(ref, hAndF.second));
51  if (hit_energy > highest_energy || highest_energy == 0.0) {
52  highest_energy = hit_energy;
53  seed = ref;
54  }
55  }
56  if (!back.hitsAndFractions().empty()) {
57  back.setSeed(seed->detId());
58  back.setEnergy(energy);
59  back.setCorrectedEnergy(energy);
60  } else {
61  back.setSeed(-1);
62  back.setEnergy(0.f);
63  }
64  }
65 }
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:46
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
void setEnergy(double energy)
Definition: CaloCluster.h:135
bool ev
void updateEvent(const edm::Event &) final
static std::string const input
Definition: EdmProvDump.cc:48
void setSeed(const DetId &id)
Definition: CaloCluster.h:145
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:209
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:136
double f[11][100]
edm::EDGetTokenT< SimClusterCollection > _simClusterToken
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:77
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
std::vector< SimCluster > SimClusterCollection
Definition: SimClusterFwd.h:8
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