CMS 3D CMS Logo

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

#include <HcalNumberingFromPS.h>

Public Member Functions

std::pair< int, int > getEta (const int &det, const math::XYZVectorD &pos) const
 
std::pair< int, int > getPhi (const int &det, const int &ieta, const double &phi) const
 
 HcalNumberingFromPS (const edm::ParameterSet &)
 
HcalNumberingFromDDD::HcalID unitID (int det, int layer, int depth, const math::XYZVectorD &pos) const
 
 ~HcalNumberingFromPS ()
 

Private Attributes

int depth29Mx_
 
std::vector< int > depthHBHE_
 
int etaHBHE_
 
std::vector< int > etaMax_
 
std::vector< int > etaMin_
 
std::vector< double > etaTable_
 
std::vector< double > phibin_
 
std::vector< double > phioff_
 
double rMinHO_
 
std::vector< std::vector< int > > segmentation_
 
std::vector< double > zHO_
 

Static Private Attributes

static const int nEtas_ = 29
 

Detailed Description

Definition at line 11 of file HcalNumberingFromPS.h.

Constructor & Destructor Documentation

◆ HcalNumberingFromPS()

HcalNumberingFromPS::HcalNumberingFromPS ( const edm::ParameterSet conf)

Definition at line 13 of file HcalNumberingFromPS.cc.

13  {
14  etaTable_ = conf.getParameter<std::vector<double> >("EtaTable");
15  phibin_ = conf.getParameter<std::vector<double> >("PhiBin");
16  phioff_ = conf.getParameter<std::vector<double> >("PhiOffset");
17  etaMin_ = conf.getParameter<std::vector<int> >("EtaMin");
18  etaMax_ = conf.getParameter<std::vector<int> >("EtaMax");
19  etaHBHE_ = conf.getParameter<int>("EtaHBHE");
20  depthHBHE_ = conf.getParameter<std::vector<int> >("DepthHBHE");
21  depth29Mx_ = conf.getParameter<int>("Depth29Max");
22  rMinHO_ = conf.getParameter<double>("RMinHO");
23  zHO_ = conf.getParameter<std::vector<double> >("ZHO");
24  const double deg = M_PI / 180.0;
25  for (auto& phi : phibin_)
26  phi *= deg;
27  for (auto& phi : phioff_)
28  phi *= deg;
29 #ifdef EDM_ML_DEBUG
30  edm::LogVerbatim("HcalSim") << "HcalNumberingFromPS:: EtaTable with " << etaTable_.size() << " elements";
31  for (unsigned k = 0; k < etaTable_.size(); ++k)
32  edm::LogVerbatim("HcalSim") << "EtaTable[" << k << "] = " << etaTable_[k];
33  edm::LogVerbatim("HcalSim") << "HcalNumberingFromPS:: PhiBin with " << phibin_.size() << " elements";
34  for (unsigned k = 0; k < phibin_.size(); ++k)
35  edm::LogVerbatim("HcalSim") << "PhiBin[" << k << "] = " << phibin_[k];
36  edm::LogVerbatim("HcalSim") << "HcalNumberingFromPS:: PhiOff with " << phioff_.size() << " elements";
37  for (unsigned k = 0; k < phioff_.size(); ++k)
38  edm::LogVerbatim("HcalSim") << "PhiOff[" << k << "] = " << phioff_[k];
39  edm::LogVerbatim("HcalSim") << "HcalNumberingFromPS:: EtaMin/EtaMax with " << etaMin_.size() << ":" << etaMax_.size()
40  << " elements";
41  for (unsigned k = 0; k < etaMin_.size(); ++k)
42  edm::LogVerbatim("HcalSim") << "EtaMin[" << k << "] = " << etaMin_[k] << " EtaMax[" << k << "] = " << etaMax_[k];
43  edm::LogVerbatim("HcalSim") << "HcalNumberingFromPS:: EtaHBHE " << etaHBHE_ << " DepthHBHE " << depthHBHE_[0] << ":"
44  << depthHBHE_[1] << " RMinHO " << rMinHO_ << " zHO with " << zHO_.size() << " elements";
45  for (unsigned k = 0; k < zHO_.size(); ++k)
46  edm::LogVerbatim("HcalSim") << "ZHO[" << k << "] = " << zHO_[k];
47 #endif
48 
49  segmentation_.resize(nEtas_);
50  for (int ring = 0; ring < nEtas_; ++ring) {
51  char name[10];
52  snprintf(name, 10, "Eta%d", ring + 1);
53  if (ring > 0) {
54  segmentation_[ring] = conf.getUntrackedParameter<std::vector<int> >(name, segmentation_[ring - 1]);
55  } else {
56  segmentation_[ring] = conf.getUntrackedParameter<std::vector<int> >(name);
57  }
58 #ifdef EDM_ML_DEBUG
59  edm::LogVerbatim("HcalSim") << "HcalNumberingFromPS:: Ring " << ring + 1 << " with " << segmentation_[ring].size()
60  << " layers";
61  for (unsigned int k = 0; k < segmentation_[ring].size(); ++k)
62  edm::LogVerbatim("HcalSim") << "Layer[" << k << "] = " << segmentation_[ring][k];
63 #endif
64  }
65 }

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), dqmdumpme::k, M_PI, Skims_PA_cff::name, and relativeConstraints::ring.

◆ ~HcalNumberingFromPS()

HcalNumberingFromPS::~HcalNumberingFromPS ( )
inline

Definition at line 14 of file HcalNumberingFromPS.h.

14 {}

Member Function Documentation

◆ getEta()

std::pair< int, int > HcalNumberingFromPS::getEta ( const int &  det,
const math::XYZVectorD pos 
) const

Definition at line 98 of file HcalNumberingFromPS.cc.

98  {
99  int ieta(1);
100  int subdet(det);
101  double eta = std::abs(pos.Eta());
102  if (pos.Rho() > rMinHO_) {
103  subdet = static_cast<int>(HcalOuter);
104  double z = std::abs(pos.z());
105  if (z > zHO_[3]) {
106  if (eta <= etaTable_[10])
107  eta = etaTable_[10] + 0.001;
108  } else if (z > zHO_[1]) {
109  if (eta <= etaTable_[4])
110  eta = etaTable_[4] + 0.001;
111  }
112  }
113  for (unsigned int i = 1; i < etaTable_.size(); i++) {
114  if (eta < etaTable_[i]) {
115  ieta = i;
116  break;
117  }
118  }
119 
120  if ((subdet == static_cast<int>(HcalBarrel)) || (subdet == static_cast<int>(HcalOuter))) {
121  if (ieta > etaMax_[0])
122  ieta = etaMax_[0];
123  } else if (det == static_cast<int>(HcalEndcap)) {
124  if (ieta <= etaMin_[1])
125  ieta = etaMin_[1];
126  }
127  return std::make_pair(subdet, ieta);
128 }

References funct::abs(), PVValHelper::eta, HcalBarrel, HcalEndcap, HcalOuter, mps_fire::i, and LEDCalibrationChannels::ieta.

◆ getPhi()

std::pair< int, int > HcalNumberingFromPS::getPhi ( const int &  det,
const int &  ieta,
const double &  phi 
) const

Definition at line 130 of file HcalNumberingFromPS.cc.

130  {
131  double fioff = ((det == static_cast<int>(HcalEndcap)) ? phioff_[1] : phioff_[0]);
132  double fibin = phibin_[ieta - 1];
133  int nphi = int((2._pi + 0.1 * fibin) / fibin);
134  double hphi = phi + fioff;
135  if (hphi < 0)
136  hphi += (2._pi);
137  int iphi = int(hphi / fibin) + 1;
138  if (iphi > nphi)
139  iphi = 1;
140  const double fiveDegInRad = 5._deg;
141  int units = int(fibin / fiveDegInRad + 0.5);
142  if (units < 1)
143  units = 1;
144  int iphi_skip = iphi;
145  if (units == 2)
146  iphi_skip = (iphi - 1) * 2 + 1;
147  else if (units == 4)
148  iphi_skip = (iphi - 1) * 4 - 1;
149  if (iphi_skip < 0)
150  iphi_skip += 72;
151  return std::make_pair(iphi, iphi_skip);
152 }

References HcalEndcap, LEDCalibrationChannels::ieta, createfilelist::int, LEDCalibrationChannels::iphi, nphi, and units().

◆ unitID()

HcalNumberingFromDDD::HcalID HcalNumberingFromPS::unitID ( int  det,
int  layer,
int  depth,
const math::XYZVectorD pos 
) const

Definition at line 67 of file HcalNumberingFromPS.cc.

70  {
71  int subdet = ((det == 3) ? static_cast<int>(HcalBarrel) : static_cast<int>(HcalEndcap));
72  std::pair<int, int> deteta = getEta(subdet, pos);
73  std::pair<int, int> iphi = getPhi(deteta.first, deteta.second, pos.Phi());
74  int newDepth(depth);
75  int zside = ((pos.z() > 0) ? 1 : 0);
76  if (deteta.first == static_cast<int>(HcalBarrel)) {
77  newDepth = segmentation_[deteta.second - 1][layer - 1];
78  if ((deteta.second == etaHBHE_) && (newDepth > depthHBHE_[0]))
79  newDepth = depthHBHE_[0];
80  } else if (deteta.first == static_cast<int>(HcalEndcap)) {
81  newDepth = segmentation_[deteta.second - 1][layer - 1];
82  if ((deteta.second == etaHBHE_) && (newDepth < depthHBHE_[1]))
83  newDepth = depthHBHE_[1];
84  if ((deteta.second == etaMax_[1]) && (newDepth > depth29Mx_))
85  --(deteta.second);
86  } else if (deteta.first == static_cast<int>(HcalOuter)) {
87  newDepth = 4;
88  }
89 #ifdef EDM_ML_DEBUG
90  edm::LogVerbatim("HcalSim") << "HcalNumberingFromPS:: det " << det << ":" << subdet << ":" << deteta.first
91  << "\t Eta " << pos.Eta() << ":" << deteta.second << "\t Phi " << pos.Phi() << ":"
92  << iphi.first << ":" << iphi.second << "\t Layer|Depth " << layer << ":" << depth << ":"
93  << newDepth;
94 #endif
95  return HcalNumberingFromDDD::HcalID(deteta.first, zside, newDepth, deteta.second, iphi.first, iphi.second, layer);
96 }

References LEDCalibrationChannels::depth, hitfit::deteta(), HcalBarrel, HcalEndcap, HcalOuter, LEDCalibrationChannels::iphi, and ecaldqm::zside().

Member Data Documentation

◆ depth29Mx_

int HcalNumberingFromPS::depth29Mx_
private

Definition at line 24 of file HcalNumberingFromPS.h.

◆ depthHBHE_

std::vector<int> HcalNumberingFromPS::depthHBHE_
private

Definition at line 23 of file HcalNumberingFromPS.h.

◆ etaHBHE_

int HcalNumberingFromPS::etaHBHE_
private

Definition at line 24 of file HcalNumberingFromPS.h.

◆ etaMax_

std::vector<int> HcalNumberingFromPS::etaMax_
private

Definition at line 23 of file HcalNumberingFromPS.h.

◆ etaMin_

std::vector<int> HcalNumberingFromPS::etaMin_
private

Definition at line 23 of file HcalNumberingFromPS.h.

◆ etaTable_

std::vector<double> HcalNumberingFromPS::etaTable_
private

Definition at line 22 of file HcalNumberingFromPS.h.

◆ nEtas_

const int HcalNumberingFromPS::nEtas_ = 29
staticprivate

Definition at line 21 of file HcalNumberingFromPS.h.

◆ phibin_

std::vector<double> HcalNumberingFromPS::phibin_
private

Definition at line 22 of file HcalNumberingFromPS.h.

◆ phioff_

std::vector<double> HcalNumberingFromPS::phioff_
private

Definition at line 22 of file HcalNumberingFromPS.h.

◆ rMinHO_

double HcalNumberingFromPS::rMinHO_
private

Definition at line 25 of file HcalNumberingFromPS.h.

◆ segmentation_

std::vector<std::vector<int> > HcalNumberingFromPS::segmentation_
private

Definition at line 27 of file HcalNumberingFromPS.h.

◆ zHO_

std::vector<double> HcalNumberingFromPS::zHO_
private

Definition at line 26 of file HcalNumberingFromPS.h.

HcalNumberingFromPS::etaMax_
std::vector< int > etaMax_
Definition: HcalNumberingFromPS.h:23
mps_fire.i
i
Definition: mps_fire.py:428
HcalNumberingFromPS::etaHBHE_
int etaHBHE_
Definition: HcalNumberingFromPS.h:24
HcalNumberingFromDDD::HcalID
Definition: HcalNumberingFromDDD.h:21
hitfit::deteta
double deteta(const Fourvec &v, double zvert)
NOT USED ANYMORE: Get the detector (D0-specific), requires z-vertex.
Definition: fourvec.cc:205
ecaldqm::zside
int zside(DetId const &)
Definition: EcalDQMCommonUtils.cc:189
pos
Definition: PixelAliasList.h:18
HcalNumberingFromPS::getPhi
std::pair< int, int > getPhi(const int &det, const int &ieta, const double &phi) const
Definition: HcalNumberingFromPS.cc:130
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
HcalBarrel
Definition: HcalAssistant.h:33
HcalNumberingFromPS::depthHBHE_
std::vector< int > depthHBHE_
Definition: HcalNumberingFromPS.h:23
units
TString units(TString variable, Char_t axis)
LEDCalibrationChannels.iphi
iphi
Definition: LEDCalibrationChannels.py:64
HcalNumberingFromPS::etaTable_
std::vector< double > etaTable_
Definition: HcalNumberingFromPS.h:22
nphi
const int nphi
Definition: CMTRawAnalyzer.h:424
HcalNumberingFromPS::rMinHO_
double rMinHO_
Definition: HcalNumberingFromPS.h:25
HcalNumberingFromPS::phioff_
std::vector< double > phioff_
Definition: HcalNumberingFromPS.h:22
PVValHelper::eta
Definition: PVValidationHelpers.h:69
DDAxes::z
HcalOuter
Definition: HcalAssistant.h:35
dqmdumpme.k
k
Definition: dqmdumpme.py:60
LEDCalibrationChannels.depth
depth
Definition: LEDCalibrationChannels.py:65
LEDCalibrationChannels.ieta
ieta
Definition: LEDCalibrationChannels.py:63
HcalNumberingFromPS::segmentation_
std::vector< std::vector< int > > segmentation_
Definition: HcalNumberingFromPS.h:27
createfilelist.int
int
Definition: createfilelist.py:10
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
HcalNumberingFromPS::depth29Mx_
int depth29Mx_
Definition: HcalNumberingFromPS.h:24
HcalNumberingFromPS::nEtas_
static const int nEtas_
Definition: HcalNumberingFromPS.h:21
DDAxes::phi
HcalNumberingFromPS::getEta
std::pair< int, int > getEta(const int &det, const math::XYZVectorD &pos) const
Definition: HcalNumberingFromPS.cc:98
HcalEndcap
Definition: HcalAssistant.h:34
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
relativeConstraints.ring
ring
Definition: relativeConstraints.py:68
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HcalNumberingFromPS::zHO_
std::vector< double > zHO_
Definition: HcalNumberingFromPS.h:26
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Log
Definition: MessageLogger.h:70
HcalNumberingFromPS::etaMin_
std::vector< int > etaMin_
Definition: HcalNumberingFromPS.h:23
HcalNumberingFromPS::phibin_
std::vector< double > phibin_
Definition: HcalNumberingFromPS.h:22