CMS 3D CMS Logo

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