CMS 3D CMS Logo

Basic2DClusterForEachSeed.cc
Go to the documentation of this file.
5 
7 public:
9  : InitialClusteringStepBase(conf, cc) {}
10  ~Basic2DClusterForEachSeed() override = default;
11 
13  const std::vector<bool>&,
14  const std::vector<bool>&,
16  const HcalPFCuts*) override;
17 };
18 
20 
22  const std::vector<bool>& rechitMask,
23  const std::vector<bool>& seedable,
25  const HcalPFCuts*) {
26  auto const& hits = *input;
27 
28  // loop over seeds and make clusters
29  reco::PFCluster cluster;
30  for (unsigned int hit = 0; hit < hits.size(); ++hit) {
31  if (!rechitMask[hit] || !seedable[hit])
32  continue; // if not seed, ignore.
33  cluster.reset();
34 
35  // seed
36  auto refhit = makeRefhit(input, hit);
37  auto rhf = reco::PFRecHitFraction(refhit, 1.0); // entire rechit energy should go to a cluster
38 
39  // add the hit to the cluster
40  cluster.addRecHitFraction(rhf);
41 
42  // extract
43  const auto rh_energy = refhit->energy();
44 
45  // fill cluster information
46  cluster.setSeed(refhit->detId());
47  cluster.setEnergy(rh_energy);
48  cluster.setTime(refhit->time());
49  cluster.setLayer(refhit->layer());
50  cluster.setPosition(math::XYZPoint(refhit->position().x(), refhit->position().y(), refhit->position().z()));
51  cluster.calculatePositionREP();
52  cluster.setDepth(refhit->depth());
53 
54  output.push_back(cluster);
55 
56  } // looping over seeds ends
57 }
void setLayer(PFLayer::Layer layer)
set layer
Definition: PFCluster.cc:49
Basic2DClusterForEachSeed(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
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
void setPosition(const math::XYZPoint &p)
Definition: CaloCluster.h:140
void setTime(float time, float timeError=0)
Definition: PFCluster.h:84
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
void setEnergy(double energy)
Definition: CaloCluster.h:136
~Basic2DClusterForEachSeed() override=default
static std::string const input
Definition: EdmProvDump.cc:50
void setSeed(const DetId &id)
Definition: CaloCluster.h:146
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:95
void buildClusters(const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &, const std::vector< bool > &, reco::PFClusterCollection &, const HcalPFCuts *) override
void reset()
resets clusters parameters
Definition: PFCluster.cc:17
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
void setDepth(double depth)
Definition: PFCluster.h:89
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: output.py:1