CMS 3D CMS Logo

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

#include <HcalHitRelabeller.h>

Public Member Functions

 HcalHitRelabeller (bool nd=false)
 
void process (std::vector< PCaloHit > &hcalHits)
 
DetId relabel (const uint32_t testId) const
 
void setGeometry (const HcalDDDRecConstants *&)
 

Static Public Member Functions

static DetId relabel (const uint32_t testId, const HcalDDDRecConstants *theRecNumber)
 

Private Member Functions

double energyWt (const uint32_t testId) const
 

Private Attributes

bool neutralDensity_
 
const HcalDDDRecConstantstheRecNumber
 

Detailed Description

Definition at line 10 of file HcalHitRelabeller.h.

Constructor & Destructor Documentation

HcalHitRelabeller::HcalHitRelabeller ( bool  nd = false)

Definition at line 9 of file HcalHitRelabeller.cc.

References neutralDensity_.

9  :
10  theRecNumber(nullptr),
11  neutralDensity_(nd) {
12 #ifdef EDM_ML_DEBUG
13  edm::LogVerbatim("HcalSim") << "HcalHitRelabeller initialized with"
14  << " neutralDensity " << neutralDensity_;
15 #endif
16 }
const HcalDDDRecConstants * theRecNumber

Member Function Documentation

double HcalHitRelabeller::energyWt ( const uint32_t  testId) const
private

Definition at line 101 of file HcalHitRelabeller.cc.

References egammaForCoreTracking_cff::depth, PVValHelper::eta, HcalDDDRecConstants::getLayer0Wt(), phi, theRecNumber, HcalTestNumbering::unpackHcalIndex(), z, and ecaldqm::zside().

Referenced by process().

101  {
102 
103  int det, z, depth, eta, phi, layer;
104  HcalTestNumbering::unpackHcalIndex(testId,det,z,depth,eta,phi,layer);
105  int zside = (z==0) ? (-1) : (1);
106  double wt = (((det==1) || (det==2)) && (depth == 1)) ?
107  theRecNumber->getLayer0Wt(det,phi,zside) : 1.0;
108 #ifdef EDM_ML_DEBUG
109  edm::LogVerbatim("HcalSim") << "EnergyWT::det: " << det << " z: " << z
110  << ":" << zside << " depth: " << depth
111  << " ieta: " << eta << " iphi: " << phi
112  << " layer: " << layer << " wt " << wt;
113 #endif
114  return wt;
115 }
int zside(DetId const &)
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
double getLayer0Wt(const int &det, const int &phi, const int &zside) const
const HcalDDDRecConstants * theRecNumber
void HcalHitRelabeller::process ( std::vector< PCaloHit > &  hcalHits)

Definition at line 18 of file HcalHitRelabeller.cc.

References TauDecayModes::dec, randomXiThetaGunProducer_cfi::energy, energyWt(), cuy::ii, neutralDensity_, DetId::rawId(), relabel(), and theRecNumber.

18  {
19 
20  if (theRecNumber) {
21 #ifdef EDM_ML_DEBUG
22  int ii(0);
23 #endif
24  for (auto & hcalHit : hcalHits) {
25 #ifdef EDM_ML_DEBUG
26  edm::LogVerbatim("HcalSim") << "Hit[" << ii << "] " << std::hex
27  << hcalHit.id() << std::dec
28  << " Neutral density " << neutralDensity_;
29 #endif
30  double energy = (hcalHit.energy());
31  if (neutralDensity_) {
32  energy *= (energyWt(hcalHit.id()));
33  hcalHit.setEnergy(energy);
34  }
35  DetId newid = relabel(hcalHit.id());
36 #ifdef EDM_ML_DEBUG
37  edm::LogVerbatim("HcalSim") << "Hit " << ii << " out of "
38  << hcalHits.size() << " " << std::hex
39  << newid.rawId() << std::dec << " E "
40  << energy;
41 #endif
42  hcalHit.setID(newid.rawId());
43 #ifdef EDM_ML_DEBUG
44  edm::LogVerbatim("HcalSim") << "Modified Hit "
45  << HcalDetId(hcalHit.id());
46  ++ii;
47 #endif
48  }
49  } else {
50  edm::LogWarning("HcalSim") << "HcalHitRelabeller: no valid HcalDDDRecConstants";
51  }
52 
53 }
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
ii
Definition: cuy.py:590
Definition: DetId.h:18
const HcalDDDRecConstants * theRecNumber
DetId relabel(const uint32_t testId) const
double energyWt(const uint32_t testId) const
DetId HcalHitRelabeller::relabel ( const uint32_t  testId) const
DetId HcalHitRelabeller::relabel ( const uint32_t  testId,
const HcalDDDRecConstants theRecNumber 
)
static

Definition at line 63 of file HcalHitRelabeller.cc.

References TauDecayModes::dec, egammaForCoreTracking_cff::depth, PVValHelper::eta, HcalDDDRecConstants::getHCID(), HcalBarrel, HcalEndcap, HcalForward, HcalOuter, phi, DetId::rawId(), Validation_hcalonly_cfi::sign, HcalTestNumbering::unpackHcalIndex(), and z.

63  {
64 
65 #ifdef EDM_ML_DEBUG
66  edm::LogVerbatim("HcalSim") << "Enter HcalHitRelabeller::relabel";
67 #endif
68  HcalDetId hid;
69  int det, z, depth, eta, phi, layer, sign;
70  HcalTestNumbering::unpackHcalIndex(testId,det,z,depth,eta,phi,layer);
71 #ifdef EDM_ML_DEBUG
72  edm::LogVerbatim("HcalSim") << "det: " << det << " "
73  << "z: " << z << " "
74  << "depth: " << depth << " "
75  << "ieta: " << eta << " "
76  << "iphi: " << phi << " "
77  << "layer: " << layer;
78 #endif
79  sign=(z==0)?(-1):(1);
80  HcalDDDRecConstants::HcalID id = theRecNumber->getHCID(det,sign*eta,phi,layer,depth);
81 
82  if (id.subdet==int(HcalBarrel)) {
83  hid=HcalDetId(HcalBarrel,sign*id.eta,id.phi,id.depth);
84  } else if (id.subdet==int(HcalEndcap)) {
85  hid=HcalDetId(HcalEndcap,sign*id.eta,id.phi,id.depth);
86  } else if (id.subdet==int(HcalOuter)) {
87  hid=HcalDetId(HcalOuter,sign*id.eta,id.phi,id.depth);
88  } else if (id.subdet==int(HcalForward)) {
89  hid=HcalDetId(HcalForward,sign*id.eta,id.phi,id.depth);
90  }
91 #ifdef EDM_ML_DEBUG
92  edm::LogVerbatim("HcalSim") << " new HcalDetId -> hex.RawID = "
93  << std::hex << hid.rawId() << std::dec
94  << " det, z, depth, eta, phi = " << det << " "
95  << z << " "<< id.depth << " " << id.eta << " "
96  << id.phi << " ---> " << hid;
97 #endif
98  return hid;
99 }
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
HcalID getHCID(int subdet, int ieta, int iphi, int lay, int idepth) const
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
void HcalHitRelabeller::setGeometry ( const HcalDDDRecConstants *&  recNum)

Definition at line 55 of file HcalHitRelabeller.cc.

References theRecNumber.

55  {
56  theRecNumber = recNum;
57 }
const HcalDDDRecConstants * theRecNumber

Member Data Documentation

bool HcalHitRelabeller::neutralDensity_
private

Definition at line 22 of file HcalHitRelabeller.h.

Referenced by HcalHitRelabeller(), and process().

const HcalDDDRecConstants* HcalHitRelabeller::theRecNumber
private

Definition at line 21 of file HcalHitRelabeller.h.

Referenced by energyWt(), process(), relabel(), and setGeometry().