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 
24 }
25 
28  const std::vector<bool>& rechitMask,
29  const std::vector<bool>& seedable,
31  const SimClusterCollection& simClusters = *_simClusterH;
32  auto const& hits = *input;
33 
34  // for quick indexing back to hit energy
35  std::unordered_map<uint32_t, size_t> detIdToIndex(hits.size());
36  for( uint32_t i = 0; i < hits.size(); ++i ) {
37  detIdToIndex[hits[i].detId()] = i;
38  auto ref = makeRefhit(input,i);
39  }
40 
41  for( const auto& sc : simClusters ) {
42  output.emplace_back();
43  reco::PFCluster& back = output.back();
45  double energy = 0.0, highest_energy = 0.0;
46  auto hitsAndFractions = sc.hits_and_fractions() ;
47  for( const auto& hAndF : hitsAndFractions ) {
48  auto itr = detIdToIndex.find(hAndF.first);
49  if( itr == detIdToIndex.end() ) continue; // hit wasn't saved in reco
50  auto ref = makeRefhit(input,itr->second);
51  const double hit_energy = hAndF.second * ref->energy();
52  energy += hit_energy;
53  back.addRecHitFraction(reco::PFRecHitFraction(ref, hAndF.second));
54  if( hit_energy > highest_energy || highest_energy == 0.0) {
55  highest_energy = hit_energy;
56  seed = ref;
57  }
58  }
59  if( !back.hitsAndFractions().empty() ) {
60  back.setSeed(seed->detId());
61  back.setEnergy(energy);
62  back.setCorrectedEnergy(energy);
63  } else {
64  back.setSeed(-1);
65  back.setEnergy(0.f);
66  }
67  }
68 }
69 
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
void setEnergy(double energy)
Definition: CaloCluster.h:111
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:195
bool ev
void updateEvent(const edm::Event &) final
static std::string const input
Definition: EdmProvDump.cc:44
void setSeed(const DetId &id)
Definition: CaloCluster.h:121
reco::PFRecHitRef makeRefhit(const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
edm::Handle< SimClusterCollection > _simClusterH
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:112
double f[11][100]
edm::EDGetTokenT< SimClusterCollection > _simClusterToken
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:99
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
std::vector< SimCluster > SimClusterCollection
Definition: SimClusterFwd.h:8