CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
hgcal::ClusterTools Class Reference

#include <ClusterTools.h>

Public Member Functions

 ClusterTools (const edm::ParameterSet &, edm::ConsumesCollector &)
 
float getClusterHadronFraction (const reco::CaloCluster &) const
 
void getEvent (const edm::Event &)
 
void getEventSetup (const edm::EventSetup &)
 
double getMultiClusterEnergy (const reco::HGCalMultiCluster &) const
 
math::XYZPoint getMultiClusterPosition (const reco::HGCalMultiCluster &, double vz=0.) const
 
 ~ClusterTools ()
 

Private Attributes

const HGCRecHitCollectionbhrh_
 
const edm::EDGetTokenT
< HGCRecHitCollection
bhtok
 
const HGCRecHitCollectioneerh_
 
const edm::EDGetTokenT
< HGCRecHitCollection
eetok
 
const HGCRecHitCollectionfhrh_
 
const edm::EDGetTokenT
< HGCRecHitCollection
fhtok
 
RecHitTools rhtools_
 

Detailed Description

Definition at line 26 of file ClusterTools.h.

Constructor & Destructor Documentation

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

Definition at line 18 of file ClusterTools.cc.

19  :
20  eetok( sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCEEInput")) ),
21  fhtok( sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCFHInput")) ),
22  bhtok( sumes.consumes<HGCRecHitCollection>(conf.getParameter<edm::InputTag>("HGCBHInput")) ) {
23 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
const edm::EDGetTokenT< HGCRecHitCollection > fhtok
Definition: ClusterTools.h:43
const edm::EDGetTokenT< HGCRecHitCollection > eetok
Definition: ClusterTools.h:43
const edm::EDGetTokenT< HGCRecHitCollection > bhtok
Definition: ClusterTools.h:43
hgcal::ClusterTools::~ClusterTools ( )
inline

Definition at line 29 of file ClusterTools.h.

29 {}

Member Function Documentation

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

Definition at line 40 of file ClusterTools.cc.

References bhrh_, eerh_, relval_parameters_module::energy, Exception, f, fhrh_, edm::SortedCollection< T, SORT >::find(), DetId::Forward, HLT_FULL_cff::fraction, DetId::Hcal, HcalEndcap, HGCEE, HGCHEF, reco::CaloCluster::hitsAndFractions(), and groupFilesInBlocks::temp.

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

Definition at line 25 of file ClusterTools.cc.

References bhrh_, bhtok, eerh_, eetok, fhrh_, fhtok, edm::Event::getByToken(), hgcal::RecHitTools::getEvent(), edm::Handle< T >::product(), rhtools_, and groupFilesInBlocks::temp.

25  {
26  rhtools_.getEvent(ev);
28  ev.getByToken(eetok, temp);
29  eerh_ = temp.product();
30  ev.getByToken(fhtok, temp);
31  fhrh_ = temp.product();
32  ev.getByToken(bhtok, temp);
33  bhrh_ = temp.product();
34 }
void getEvent(const edm::Event &)
Definition: RecHitTools.cc:63
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
RecHitTools rhtools_
Definition: ClusterTools.h:42
const HGCRecHitCollection * eerh_
Definition: ClusterTools.h:44
T const * product() const
Definition: Handle.h:81
const edm::EDGetTokenT< HGCRecHitCollection > fhtok
Definition: ClusterTools.h:43
const edm::EDGetTokenT< HGCRecHitCollection > eetok
Definition: ClusterTools.h:43
const HGCRecHitCollection * bhrh_
Definition: ClusterTools.h:44
const edm::EDGetTokenT< HGCRecHitCollection > bhtok
Definition: ClusterTools.h:43
const HGCRecHitCollection * fhrh_
Definition: ClusterTools.h:44
void ClusterTools::getEventSetup ( const edm::EventSetup es)

Definition at line 36 of file ClusterTools.cc.

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

36  {
38 }
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:66
RecHitTools rhtools_
Definition: ClusterTools.h:42
double ClusterTools::getMultiClusterEnergy ( const reco::HGCalMultiCluster clu) const

Definition at line 105 of file ClusterTools.cc.

References reco::HGCalMultiCluster::clusters().

105  {
106  double acc = 0.0;
107  for(const auto& ptr : clu.clusters() ) {
108  acc += ptr->energy();
109  }
110  return acc;
111 }
const edm::PtrVector< reco::BasicCluster > & clusters() const
math::XYZPoint ClusterTools::getMultiClusterPosition ( const reco::HGCalMultiCluster clu,
double  vz = 0. 
) const

Definition at line 79 of file ClusterTools.cc.

References assert(), reco::HGCalMultiCluster::clusters(), f, mathSSE::sqrt(), groupFilesInBlocks::temp, histoStyle::weight, and x().

79  {
80  if( clu.clusters().size() == 0 ) return math::XYZPoint();
81  double acc_rho = 0.0;
82  double acc_eta = 0.0;
83  double acc_phi = 0.0;
84  double totweight = 0.;
85  for( const auto& ptr : clu.clusters() ) {
86  const double x = ptr->x();
87  const double y = ptr->y();
88  const double point_r2 = (x*x + y*y);
89  const double point_z = ptr->z()-vz;
90  const double point_h = std::sqrt(point_r2 + point_z*point_z);
91  const double weight = ptr->energy() * ptr->size();
92  assert((y != 0. || x != 0.) && "Cluster position somehow in beampipe.");
93  assert(point_z != 0.f && "Layer-cluster position given as reference point.");
94  const double point_r = std::sqrt(point_r2);
95  acc_rho += point_r * weight;
96  acc_phi += vdt::fast_atan2(y,x) * weight;
97  acc_eta += -1. * vdt::fast_log(point_r/(point_z + point_h)) * weight;
98  totweight += weight;
99  }
100  const double invweight = 1.0/totweight;
101  reco::PFCluster::REPPoint temp(acc_rho*invweight,acc_eta*invweight,acc_phi*invweight);
102  return math::XYZPoint(temp.x(),temp.y(),temp.z());
103 }
assert(m_qm.get())
const edm::PtrVector< reco::BasicCluster > & clusters() const
T sqrt(T t)
Definition: SSEVec.h:18
double f[11][100]
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:54
int weight
Definition: histoStyle.py:50

Member Data Documentation

const HGCRecHitCollection * hgcal::ClusterTools::bhrh_
private

Definition at line 44 of file ClusterTools.h.

Referenced by getClusterHadronFraction(), and getEvent().

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

Definition at line 43 of file ClusterTools.h.

Referenced by getEvent().

const HGCRecHitCollection* hgcal::ClusterTools::eerh_
private

Definition at line 44 of file ClusterTools.h.

Referenced by getClusterHadronFraction(), and getEvent().

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

Definition at line 43 of file ClusterTools.h.

Referenced by getEvent().

const HGCRecHitCollection * hgcal::ClusterTools::fhrh_
private

Definition at line 44 of file ClusterTools.h.

Referenced by getClusterHadronFraction(), and getEvent().

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

Definition at line 43 of file ClusterTools.h.

Referenced by getEvent().

RecHitTools hgcal::ClusterTools::rhtools_
private

Definition at line 42 of file ClusterTools.h.

Referenced by getEvent(), and getEventSetup().