CMS 3D CMS Logo

PFClusterFromHGCalTrackster.cc
Go to the documentation of this file.
7 
9 public:
11  : InitialClusteringStepBase(conf, cc) {
12  filterByTracksterPID_ = conf.getParameter<bool>("filterByTracksterPID");
13  filterByTracksterIteration_ = conf.getParameter<bool>("filterByTracksterIteration");
14  pid_threshold_ = conf.getParameter<double>("pid_threshold");
15  filter_on_categories_ = conf.getParameter<std::vector<int> >("filter_on_categories");
16  filter_on_iterations_ = conf.getParameter<std::vector<int> >("filter_on_iterations");
17 
18  tracksterToken_ = cc.consumes<std::vector<ticl::Trackster> >(conf.getParameter<edm::InputTag>("tracksterSrc"));
20  }
21 
25 
26  void updateEvent(const edm::Event&) final;
27 
29  const std::vector<bool>&,
30  const std::vector<bool>&,
31  reco::PFClusterCollection&) override;
32 
33 private:
37  std::vector<int> filter_on_categories_;
38  std::vector<int> filter_on_iterations_;
39 
42 
45 };
46 
48 
50  ev.getByToken(tracksterToken_, trackstersH_);
51  ev.getByToken(clusterToken_, clusterH_);
52 }
53 
55  const std::vector<bool>& rechitMask,
56  const std::vector<bool>& seedable,
58  auto const& hits = *input;
59 
60  const auto& tracksters = *trackstersH_;
61  const auto& clusters = *clusterH_;
62 
63  // for quick indexing back to hit energy
64  std::unordered_map<uint32_t, size_t> detIdToIndex(hits.size());
65  for (uint32_t i = 0; i < hits.size(); ++i) {
66  detIdToIndex[hits[i].detId()] = i;
67  }
68 
69  for (const auto& tst : tracksters) {
70  // Skip empty tracksters
71  if (tst.vertices().empty()) {
72  continue;
73  }
74  // Filter using trackster PID
76  float probTotal = 0.0f;
77  for (int cat : filter_on_categories_) {
78  probTotal += tst.id_probabilities(cat);
79  }
80  if (probTotal < pid_threshold_) {
81  continue;
82  }
83  }
84 
86  std::find(filter_on_iterations_.begin(), filter_on_iterations_.end(), tst.ticlIteration()) ==
87  filter_on_iterations_.end()) {
88  continue;
89  }
90 
91  DetId seed;
92  double energy = 0.0, highest_energy = 0.0;
93  output.emplace_back();
94  reco::PFCluster& back = output.back();
95 
96  std::vector<std::pair<DetId, float> > hitsAndFractions;
97  int iLC = 0;
98  std::for_each(std::begin(tst.vertices()), std::end(tst.vertices()), [&](unsigned int lcId) {
99  const auto fraction = 1.f / tst.vertex_multiplicity(iLC++);
100  for (const auto& cell : clusters[lcId].hitsAndFractions()) {
101  hitsAndFractions.emplace_back(cell.first, cell.second * fraction);
102  }
103  });
104 
105  for (const auto& hAndF : hitsAndFractions) {
106  auto itr = detIdToIndex.find(hAndF.first);
107  if (itr == detIdToIndex.end()) {
108  continue; // hit wasn't saved in reco
109  }
110  auto ref = makeRefhit(input, itr->second);
111  assert(ref->detId() == hAndF.first.rawId());
112  const double hit_energy = hAndF.second * ref->energy();
113  energy += hit_energy;
114  back.addRecHitFraction(reco::PFRecHitFraction(ref, hAndF.second));
115  // TODO: the following logic to identify the seed of a cluster
116  // could be appropriate for the Run2 Ecal Calorimetric
117  // detector, but could be wrong for the HGCal one. This has to
118  // be reviewd.
119  if (hit_energy > highest_energy || highest_energy == 0.0) {
120  highest_energy = hit_energy;
121  seed = ref->detId();
122  }
123  } // end of hitsAndFractions
124 
125  if (!back.hitsAndFractions().empty()) {
126  back.setSeed(seed);
127  back.setEnergy(energy);
129  } else {
130  back.setSeed(0);
131  back.setEnergy(0.f);
132  }
133  } // end of loop over hgcalTracksters (3D)
134 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
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
edm::Handle< reco::CaloClusterCollection > clusterH_
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
void setEnergy(double energy)
Definition: CaloCluster.h:136
void updateEvent(const edm::Event &) final
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
static std::string const input
Definition: EdmProvDump.cc:50
void setSeed(const DetId &id)
Definition: CaloCluster.h:146
PFClusterFromHGCalTrackster & operator=(const PFClusterFromHGCalTrackster &)=delete
edm::EDGetTokenT< reco::CaloClusterCollection > clusterToken_
def cat(path)
Definition: eostools.py:401
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:137
std::vector< CaloCluster > CaloClusterCollection
collection of CaloCluster objects
double f[11][100]
Definition: DetId.h:17
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)
edm::Handle< std::vector< ticl::Trackster > > trackstersH_
edm::EDGetTokenT< std::vector< ticl::Trackster > > tracksterToken_
PFClusterFromHGCalTrackster(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
void buildClusters(const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &, const std::vector< bool > &, reco::PFClusterCollection &) override