CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
hgcal::ClusterTools Class Reference

#include <ClusterTools.h>

Public Member Functions

 ClusterTools ()
 
 ClusterTools (const edm::ParameterSet &, edm::ConsumesCollector &)
 
float getClusterHadronFraction (const reco::CaloCluster &) const
 
void getEvent (const edm::Event &)
 
void getEventSetup (const edm::EventSetup &)
 
int getLayer (const DetId) const
 
double getMultiClusterEnergy (const reco::HGCalMultiCluster &) const
 
math::XYZPoint getMultiClusterPosition (const reco::HGCalMultiCluster &) const
 
bool getWidths (const reco::CaloCluster &clus, double &sigmaetaeta, double &sigmaphiphi, double &sigmaetaetalog, double &sigmaphiphilog) const
 
 ~ClusterTools ()
 

Private Member Functions

std::vector< size_t > sort_by_z (const reco::HGCalMultiCluster &v) const
 

Private Attributes

const HGCRecHitCollectionbhrh_
 
const edm::EDGetTokenT< HGCRecHitCollectionbhtok
 
const edm::ESGetToken< CaloGeometry, CaloGeometryRecordcaloGeometryToken_
 
const HGCRecHitCollectioneerh_
 
const edm::EDGetTokenT< HGCRecHitCollectioneetok
 
const HGCRecHitCollectionfhrh_
 
const edm::EDGetTokenT< HGCRecHitCollectionfhtok
 
RecHitTools rhtools_
 

Detailed Description

Definition at line 28 of file ClusterTools.h.

Constructor & Destructor Documentation

◆ ClusterTools() [1/2]

ClusterTools::ClusterTools ( )

Definition at line 19 of file ClusterTools.cc.

19 {}

◆ ClusterTools() [2/2]

ClusterTools::ClusterTools ( const edm::ParameterSet conf,
edm::ConsumesCollector sumes 
)

Definition at line 21 of file ClusterTools.cc.

References edm::ConsumesCollector::esConsumes().

22  : eetok(sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCEEInput"))),
23  fhtok(sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCFHInput"))),
24  bhtok(sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCBHInput"))),
25  caloGeometryToken_{sumes.esConsumes()} {}
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_
Definition: ClusterTools.h:63
const edm::EDGetTokenT< HGCRecHitCollection > fhtok
Definition: ClusterTools.h:62
const edm::EDGetTokenT< HGCRecHitCollection > eetok
Definition: ClusterTools.h:62
const edm::EDGetTokenT< HGCRecHitCollection > bhtok
Definition: ClusterTools.h:62

◆ ~ClusterTools()

hgcal::ClusterTools::~ClusterTools ( )
inline

Definition at line 32 of file ClusterTools.h.

32 {}

Member Function Documentation

◆ getClusterHadronFraction()

float ClusterTools::getClusterHadronFraction ( const reco::CaloCluster clus) const

Definition at line 35 of file ClusterTools.cc.

References bhrh_, eerh_, hcalRecHitTable_cff::energy, Exception, f, fhrh_, edm::SortedCollection< T, SORT >::find(), DetId::Forward, HLT_2024v14_cff::fraction, DetId::Hcal, HcalEndcap, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, HGCEE, HGCHEF, hfClusterShapes_cfi::hits, reco::CaloCluster::hitsAndFractions(), and groupFilesInBlocks::temp.

35  {
36  float energy = 0.f, energyHad = 0.f;
37  const auto& hits = clus.hitsAndFractions();
38  for (const auto& hit : hits) {
39  const auto& id = hit.first;
40  const float fraction = hit.second;
41  if (id.det() == DetId::HGCalEE) {
42  energy += eerh_->find(id)->energy() * fraction;
43  } else if (id.det() == DetId::HGCalHSi) {
44  const float temp = fhrh_->find(id)->energy();
45  energy += temp * fraction;
46  energyHad += temp * fraction;
47  } else if (id.det() == DetId::HGCalHSc) {
48  const float temp = bhrh_->find(id)->energy();
49  energy += temp * fraction;
50  energyHad += temp * fraction;
51  } else if (id.det() == DetId::Forward) {
52  switch (id.subdetId()) {
53  case HGCEE:
54  energy += eerh_->find(id)->energy() * fraction;
55  break;
56  case HGCHEF: {
57  const float temp = fhrh_->find(id)->energy();
58  energy += temp * fraction;
59  energyHad += temp * fraction;
60  } break;
61  default:
62  throw cms::Exception("HGCalClusterTools") << " Cluster contains hits that are not from HGCal! " << std::endl;
63  }
64  } else if (id.det() == DetId::Hcal && id.subdetId() == HcalEndcap) {
65  const float temp = bhrh_->find(id)->energy();
66  energy += temp * fraction;
67  energyHad += temp * fraction;
68  } else {
69  throw cms::Exception("HGCalClusterTools") << " Cluster contains hits that are not from HGCal! " << std::endl;
70  }
71  }
72  float fraction = -1.f;
73  if (energy > 0.f) {
74  fraction = energyHad / energy;
75  }
76  return fraction;
77 }
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
double f[11][100]
const HGCRecHitCollection * eerh_
Definition: ClusterTools.h:65
iterator find(key_type k)
const HGCRecHitCollection * bhrh_
Definition: ClusterTools.h:65
const HGCRecHitCollection * fhrh_
Definition: ClusterTools.h:65

◆ getEvent()

void ClusterTools::getEvent ( const edm::Event ev)

Definition at line 27 of file ClusterTools.cc.

References bhrh_, bhtok, eerh_, eetok, makeMEIFBenchmarkPlots::ev, fhrh_, and fhtok.

27  {
28  eerh_ = &ev.get(eetok);
29  fhrh_ = &ev.get(fhtok);
30  bhrh_ = &ev.get(bhtok);
31 }
const HGCRecHitCollection * eerh_
Definition: ClusterTools.h:65
const edm::EDGetTokenT< HGCRecHitCollection > fhtok
Definition: ClusterTools.h:62
const edm::EDGetTokenT< HGCRecHitCollection > eetok
Definition: ClusterTools.h:62
const HGCRecHitCollection * bhrh_
Definition: ClusterTools.h:65
const edm::EDGetTokenT< HGCRecHitCollection > bhtok
Definition: ClusterTools.h:62
const HGCRecHitCollection * fhrh_
Definition: ClusterTools.h:65

◆ getEventSetup()

void ClusterTools::getEventSetup ( const edm::EventSetup es)

Definition at line 33 of file ClusterTools.cc.

References caloGeometryToken_, edm::EventSetup::getData(), rhtools_, and hgcal::RecHitTools::setGeometry().

T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_
Definition: ClusterTools.h:63
RecHitTools rhtools_
Definition: ClusterTools.h:61
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:79

◆ getLayer()

int ClusterTools::getLayer ( const DetId  detid) const

Definition at line 108 of file ClusterTools.cc.

References hgcal::RecHitTools::getLayerWithOffset(), and rhtools_.

Referenced by getWidths().

108 { return rhtools_.getLayerWithOffset(detid); }
RecHitTools rhtools_
Definition: ClusterTools.h:61
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381

◆ getMultiClusterEnergy()

double ClusterTools::getMultiClusterEnergy ( const reco::HGCalMultiCluster clu) const

Definition at line 110 of file ClusterTools.cc.

References reco::HGCalMultiCluster::clusters().

Referenced by getMultiClusterPosition().

110  {
111  double acc = 0.0;
112  for (const auto& ptr : clu.clusters()) {
113  acc += ptr->energy();
114  }
115  return acc;
116 }
const edm::PtrVector< reco::BasicCluster > & clusters() const

◆ getMultiClusterPosition()

math::XYZPoint ClusterTools::getMultiClusterPosition ( const reco::HGCalMultiCluster clu) const

Definition at line 79 of file ClusterTools.cc.

References reco::HGCalMultiCluster::clusters(), getMultiClusterEnergy(), and mps_merge::weight.

79  {
80  if (clu.clusters().empty())
81  return math::XYZPoint();
82 
83  double acc_x = 0.0;
84  double acc_y = 0.0;
85  double acc_z = 0.0;
86  double totweight = 0.;
87  double mcenergy = getMultiClusterEnergy(clu);
88  for (const auto& ptr : clu.clusters()) {
89  if (mcenergy != 0) {
90  if (ptr->energy() < .01 * mcenergy)
91  continue; //cutoff < 1% layer contribution
92  }
93  const double weight = ptr->energy(); // weigth each corrdinate only by the total energy of the layer cluster
94  acc_x += ptr->x() * weight;
95  acc_y += ptr->y() * weight;
96  acc_z += ptr->z() * weight;
97  totweight += weight;
98  }
99  if (totweight != 0) {
100  acc_x /= totweight;
101  acc_y /= totweight;
102  acc_z /= totweight;
103  }
104  // return x/y/z in absolute coordinates
105  return math::XYZPoint(acc_x, acc_y, acc_z);
106 }
double getMultiClusterEnergy(const reco::HGCalMultiCluster &) const
Definition: weight.py:1
const edm::PtrVector< reco::BasicCluster > & clusters() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12

◆ getWidths()

bool ClusterTools::getWidths ( const reco::CaloCluster clus,
double &  sigmaetaeta,
double &  sigmaphiphi,
double &  sigmaetaetalog,
double &  sigmaphiphilog 
) const

Definition at line 118 of file ClusterTools.cc.

References eerh_, CaloRecHit::energy(), reco::CaloCluster::energy(), PV3DBase< T, PVType, FrameType >::eta(), edm::SortedCollection< T, SORT >::find(), dqmdumpme::first, DetId::Forward, getLayer(), hgcal::RecHitTools::getPosition(), DetId::HGCalEE, HGCEE, reco::CaloCluster::hitsAndFractions(), hgcal::RecHitTools::lastLayerEE(), dqm-mbProfile::log, SiStripPI::max, PV3DBase< T, PVType, FrameType >::phi(), reco::CaloCluster::position(), position, rhtools_, edm::second(), mathSSE::sqrt(), and mps_merge::weight.

122  {
123  if (getLayer(clus.hitsAndFractions()[0].first) > (int)rhtools_.lastLayerEE())
124  return false;
125  const math::XYZPoint& position(clus.position());
126  unsigned nhit = clus.hitsAndFractions().size();
127 
128  sigmaetaeta = 0.;
129  sigmaphiphi = 0.;
130  sigmaetaetal = 0.;
131  sigmaphiphil = 0.;
132 
133  double sumw = 0.;
134  double sumlogw = 0.;
135 
136  for (unsigned int ih = 0; ih < nhit; ++ih) {
137  const DetId& id = (clus.hitsAndFractions())[ih].first;
138  if ((clus.hitsAndFractions())[ih].second == 0.)
139  continue;
140 
141  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::Forward && id.subdetId() == HGCEE)) {
142  const HGCRecHit* theHit = &(*eerh_->find(id));
143 
144  GlobalPoint cellPos = rhtools_.getPosition(id);
145  double weight = theHit->energy();
146  // take w0=2 To be optimized
147  double logweight = 0;
148  if (clus.energy() != 0) {
149  logweight = std::max(0., 2 + log(theHit->energy() / clus.energy()));
150  }
151  double deltaetaeta2 = (cellPos.eta() - position.eta()) * (cellPos.eta() - position.eta());
152  double deltaphiphi2 = (cellPos.phi() - position.phi()) * (cellPos.phi() - position.phi());
153  sigmaetaeta += deltaetaeta2 * weight;
154  sigmaphiphi += deltaphiphi2 * weight;
155  sigmaetaetal += deltaetaeta2 * logweight;
156  sigmaphiphil += deltaphiphi2 * logweight;
157  sumw += weight;
158  sumlogw += logweight;
159  }
160  }
161 
162  if (sumw <= 0.)
163  return false;
164 
165  sigmaetaeta /= sumw;
166  sigmaetaeta = std::sqrt(sigmaetaeta);
167  sigmaphiphi /= sumw;
168  sigmaphiphi = std::sqrt(sigmaphiphi);
169 
170  if (sumlogw != 0) {
171  sigmaetaetal /= sumlogw;
172  sigmaetaetal = std::sqrt(sigmaetaetal);
173  sigmaphiphil /= sumlogw;
174  sigmaphiphil = std::sqrt(sigmaphiphil);
175  }
176 
177  return true;
178 }
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:154
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
T eta() const
Definition: PV3DBase.h:73
Definition: weight.py:1
constexpr float energy() const
Definition: CaloRecHit.h:29
int getLayer(const DetId) const
U second(std::pair< T, U > const &p)
RecHitTools rhtools_
Definition: ClusterTools.h:61
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:140
T sqrt(T t)
Definition: SSEVec.h:23
const HGCRecHitCollection * eerh_
Definition: ClusterTools.h:65
double energy() const
cluster energy
Definition: CaloCluster.h:149
Definition: DetId.h:17
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
iterator find(key_type k)
static int position[264][3]
Definition: ReadPGInfo.cc:289
unsigned int lastLayerEE(bool nose=false) const
Definition: RecHitTools.h:76

◆ sort_by_z()

std::vector<size_t> hgcal::ClusterTools::sort_by_z ( const reco::HGCalMultiCluster v) const
inlineprivate

Definition at line 53 of file ClusterTools.h.

References testProducerWithPsetDescEmpty_cfi::i1, testProducerWithPsetDescEmpty_cfi::i2, heavyIonCSV_trainingSettings::idx, jetUpdater_cfi::sort, and findQualityFiles::v.

53  {
54  std::vector<size_t> idx(v.size());
55  std::iota(std::begin(idx), std::end(idx), 0);
56  sort(
57  idx.begin(), idx.end(), [&v](size_t i1, size_t i2) { return v.clusters()[i1]->z() < v.clusters()[i2]->z(); });
58  return idx;
59  }

Member Data Documentation

◆ bhrh_

const HGCRecHitCollection * hgcal::ClusterTools::bhrh_
private

Definition at line 65 of file ClusterTools.h.

Referenced by getClusterHadronFraction(), and getEvent().

◆ bhtok

const edm::EDGetTokenT<HGCRecHitCollection> hgcal::ClusterTools::bhtok
private

Definition at line 62 of file ClusterTools.h.

Referenced by getEvent().

◆ caloGeometryToken_

const edm::ESGetToken<CaloGeometry, CaloGeometryRecord> hgcal::ClusterTools::caloGeometryToken_
private

Definition at line 63 of file ClusterTools.h.

Referenced by getEventSetup().

◆ eerh_

const HGCRecHitCollection* hgcal::ClusterTools::eerh_
private

Definition at line 65 of file ClusterTools.h.

Referenced by getClusterHadronFraction(), getEvent(), and getWidths().

◆ eetok

const edm::EDGetTokenT<HGCRecHitCollection> hgcal::ClusterTools::eetok
private

Definition at line 62 of file ClusterTools.h.

Referenced by getEvent().

◆ fhrh_

const HGCRecHitCollection * hgcal::ClusterTools::fhrh_
private

Definition at line 65 of file ClusterTools.h.

Referenced by getClusterHadronFraction(), and getEvent().

◆ fhtok

const edm::EDGetTokenT<HGCRecHitCollection> hgcal::ClusterTools::fhtok
private

Definition at line 62 of file ClusterTools.h.

Referenced by getEvent().

◆ rhtools_

RecHitTools hgcal::ClusterTools::rhtools_
private

Definition at line 61 of file ClusterTools.h.

Referenced by getEventSetup(), getLayer(), and getWidths().