CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
PFClusterFromHGCalMultiCluster Class Reference

#include <PFClusterFromHGCalMultiCluster.h>

Inheritance diagram for PFClusterFromHGCalMultiCluster:
InitialClusteringStepBase

Public Member Functions

void buildClusters (const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &, const std::vector< bool > &, reco::PFClusterCollection &) override
 
PFClusterFromHGCalMultiClusteroperator= (const PFClusterFromHGCalMultiCluster &)=delete
 
 PFClusterFromHGCalMultiCluster (const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
 
 PFClusterFromHGCalMultiCluster (const PFClusterFromHGCalMultiCluster &)=delete
 
void updateEvent (const edm::Event &) final
 
 ~PFClusterFromHGCalMultiCluster () override
 
- Public Member Functions inherited from InitialClusteringStepBase
 InitialClusteringStepBase (const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
 
 InitialClusteringStepBase (const ICSB &)=delete
 
std::ostream & operator<< (std::ostream &o) const
 
ICSBoperator= (const ICSB &)=delete
 
void reset ()
 
virtual void update (const edm::EventSetup &)
 
virtual ~InitialClusteringStepBase ()=default
 

Private Attributes

edm::Handle< std::vector< reco::HGCalMultiCluster > > clusterH_
 
edm::EDGetTokenT< std::vector< reco::HGCalMultiCluster > > clusterToken_
 

Additional Inherited Members

- Protected Types inherited from InitialClusteringStepBase
typedef std::tuple< std::vector< int >, std::vector< double >, std::vector< double > > I3tuple
 
- Protected Member Functions inherited from InitialClusteringStepBase
reco::PFRecHitRef makeRefhit (const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
 
- Protected Attributes inherited from InitialClusteringStepBase
const std::unordered_map< std::string, int > _layerMap
 
unsigned _nClustersFound
 
unsigned _nSeeds
 
std::unordered_map< int, I3tuple_thresholds
 

Detailed Description

Definition at line 8 of file PFClusterFromHGCalMultiCluster.h.

Constructor & Destructor Documentation

PFClusterFromHGCalMultiCluster::PFClusterFromHGCalMultiCluster ( const edm::ParameterSet conf,
edm::ConsumesCollector sumes 
)
inline

Definition at line 10 of file PFClusterFromHGCalMultiCluster.h.

References clusterToken_, edm::ConsumesCollector::consumes(), and edm::ParameterSet::getParameter().

Referenced by ~PFClusterFromHGCalMultiCluster().

11  : InitialClusteringStepBase(conf, sumes) {
13  sumes.consumes<std::vector<reco::HGCalMultiCluster> >(conf.getParameter<edm::InputTag>("clusterSrc"));
14  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
InitialClusteringStepBase(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
edm::EDGetTokenT< std::vector< reco::HGCalMultiCluster > > clusterToken_
PFClusterFromHGCalMultiCluster::~PFClusterFromHGCalMultiCluster ( )
inlineoverride
PFClusterFromHGCalMultiCluster::PFClusterFromHGCalMultiCluster ( const PFClusterFromHGCalMultiCluster )
delete

Member Function Documentation

void PFClusterFromHGCalMultiCluster::buildClusters ( const edm::Handle< reco::PFRecHitCollection > &  input,
const std::vector< bool > &  rechitMask,
const std::vector< bool > &  seedable,
reco::PFClusterCollection output 
)
overridevirtual

Implements InitialClusteringStepBase.

Definition at line 10 of file PFClusterFromHGCalMultiCluster.cc.

References reco::PFCluster::addRecHitFraction(), GetRecoTauVFromDQM_MC_cff::cl, clusterH_, HCALHighEnergyHPDFilter_cfi::energy, f, hfClusterShapes_cfi::hits, reco::CaloCluster::hitsAndFractions(), mps_fire::i, input, InitialClusteringStepBase::makeRefhit(), SurveyInfoScenario_cff::seed, reco::CaloCluster::setCorrectedEnergy(), reco::CaloCluster::setEnergy(), and reco::CaloCluster::setSeed().

Referenced by ~PFClusterFromHGCalMultiCluster().

13  {
14  const auto& hgcalMultiClusters = *clusterH_;
15  auto const& hits = *input;
16 
17  // for quick indexing back to hit energy
18  std::unordered_map<uint32_t, size_t> detIdToIndex(hits.size());
19  for (uint32_t i = 0; i < hits.size(); ++i) {
20  detIdToIndex[hits[i].detId()] = i;
21  }
22 
23  for (const auto& mcl : hgcalMultiClusters) {
24  DetId seed;
25  double energy = 0.0, highest_energy = 0.0;
26  output.emplace_back();
27  reco::PFCluster& back = output.back();
28  for (const auto& cl : mcl) {
29  const auto& hitsAndFractions = cl->hitsAndFractions();
30  for (const auto& hAndF : hitsAndFractions) {
31  auto itr = detIdToIndex.find(hAndF.first);
32  if (itr == detIdToIndex.end()) {
33  continue; // hit wasn't saved in reco
34  }
35  auto ref = makeRefhit(input, itr->second);
36  assert(ref->detId() == hAndF.first.rawId());
37  const double hit_energy = hAndF.second * ref->energy();
38  energy += hit_energy;
39  back.addRecHitFraction(reco::PFRecHitFraction(ref, hAndF.second));
40  // TODO: the following logic to identify the seed of a cluster
41  // could be appropriate for the Run2 Ecal Calorimetric
42  // detector, but could be wrong for the HGCal one. This has to
43  // be reviewd.
44  if (hit_energy > highest_energy || highest_energy == 0.0) {
45  highest_energy = hit_energy;
46  seed = ref->detId();
47  }
48  } // end of hitsAndFractions
49  } // end of loop over clusters (2D/layer)
50  if (energy <= 1) {
51  output.pop_back();
52  continue;
53  }
54  if (!back.hitsAndFractions().empty()) {
55  back.setSeed(seed);
56  back.setEnergy(energy);
57  back.setCorrectedEnergy(energy);
58  } else {
59  back.setSeed(0);
60  back.setEnergy(0.f);
61  }
62  } // end of loop over hgcalMulticlusters (3D)
63 }
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:46
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
void setEnergy(double energy)
Definition: CaloCluster.h:135
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
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:209
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:136
double f[11][100]
Definition: DetId.h:17
edm::Handle< std::vector< reco::HGCalMultiCluster > > clusterH_
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:77
PFClusterFromHGCalMultiCluster& PFClusterFromHGCalMultiCluster::operator= ( const PFClusterFromHGCalMultiCluster )
delete
void PFClusterFromHGCalMultiCluster::updateEvent ( const edm::Event ev)
finalvirtual

Reimplemented from InitialClusteringStepBase.

Definition at line 8 of file PFClusterFromHGCalMultiCluster.cc.

References clusterH_, clusterToken_, and edm::Event::getByToken().

Referenced by ~PFClusterFromHGCalMultiCluster().

bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::Handle< std::vector< reco::HGCalMultiCluster > > clusterH_
edm::EDGetTokenT< std::vector< reco::HGCalMultiCluster > > clusterToken_

Member Data Documentation

edm::Handle<std::vector<reco::HGCalMultiCluster> > PFClusterFromHGCalMultiCluster::clusterH_
private

Definition at line 28 of file PFClusterFromHGCalMultiCluster.h.

Referenced by buildClusters(), and updateEvent().

edm::EDGetTokenT<std::vector<reco::HGCalMultiCluster> > PFClusterFromHGCalMultiCluster::clusterToken_
private

Definition at line 27 of file PFClusterFromHGCalMultiCluster.h.

Referenced by PFClusterFromHGCalMultiCluster(), and updateEvent().