CMS 3D CMS Logo

HcalFlexiHardcodeGeometryLoader.cc
Go to the documentation of this file.
9 
10 #include <vector>
11 #include <algorithm>
12 
14 
15 //#define EDM_ML_DEBUG
16 
17 // ==============> Loader Itself <==========================
18 
20  MAX_HCAL_PHI = 72;
21  DEGREE2RAD = M_PI / 180.;
22 }
23 
25  const HcalDDDRecConstants& hcons) {
26  HcalGeometry* hcalGeometry = new HcalGeometry(fTopology);
27  if (nullptr == hcalGeometry->cornersMgr())
28  hcalGeometry->allocateCorners(fTopology.ncells() + fTopology.getHFSize());
29  if (nullptr == hcalGeometry->parMgr())
31  isBH_ = hcons.isBH();
32 #ifdef EDM_ML_DEBUG
33  std::cout << "FlexiGeometryLoader initialize with ncells " << fTopology.ncells() << " and shapes "
34  << hcalGeometry->numberOfShapes() << ":" << HcalGeometry::k_NumberOfParametersPerShape << " with BH Flag "
35  << isBH_ << std::endl;
36 #endif
37  if (fTopology.mode() == HcalTopologyMode::H2) { // TB geometry
38  fillHBHO(hcalGeometry, makeHBCells(hcons), true);
39  fillHBHO(hcalGeometry, makeHOCells(), false);
40  fillHE(hcalGeometry, makeHECells_H2());
41  } else { // regular geometry
42  fillHBHO(hcalGeometry, makeHBCells(hcons), true);
43  fillHBHO(hcalGeometry, makeHOCells(), false);
44  fillHF(hcalGeometry, makeHFCells(hcons));
45  fillHE(hcalGeometry, makeHECells(hcons));
46  }
47  //fast insertion of valid ids requires sort at end
48  hcalGeometry->sortValidIds();
49  return hcalGeometry;
50 }
51 
52 // ----------> HB <-----------
53 std::vector<HcalFlexiHardcodeGeometryLoader::HBHOCellParameters> HcalFlexiHardcodeGeometryLoader::makeHBCells(
54  const HcalDDDRecConstants& hcons) {
55  std::vector<HcalFlexiHardcodeGeometryLoader::HBHOCellParameters> result;
56  std::vector<std::pair<double, double> > gconsHB = hcons.getConstHBHE(0);
57  std::vector<HcalDDDRecConstants::HcalEtaBin> etabins = hcons.getEtaBins(0);
58 
59 #ifdef EDM_ML_DEBUG
60  std::cout << "FlexiGeometryLoader called for " << etabins.size() << " Eta Bins" << std::endl;
61  for (unsigned int k = 0; k < gconsHB.size(); ++k) {
62  std::cout << "gconsHB[" << k << "] = " << gconsHB[k].first << " +- " << gconsHB[k].second << std::endl;
63  }
64 #endif
65  for (auto& etabin : etabins) {
66  int iring = (etabin.zside >= 0) ? etabin.ieta : -etabin.ieta;
67  int depth = etabin.depthStart;
68  double dphi =
69  (etabin.phis.size() > 1) ? (etabin.phis[1].second - etabin.phis[0].second) : ((2.0 * M_PI) / MAX_HCAL_PHI);
70  for (unsigned int k = 0; k < etabin.layer.size(); ++k) {
71  int layf = etabin.layer[k].first - 1;
72  int layl = etabin.layer[k].second - 1;
73  double rmin = gconsHB[layf].first - gconsHB[layf].second;
74  double rmax = gconsHB[layl].first + gconsHB[layl].second;
75  for (unsigned int j = 0; j < etabin.phis.size(); ++j) {
76 #ifdef EDM_ML_DEBUG
77  std::cout << "HBRing " << iring << " eta " << etabin.etaMin << ":" << etabin.etaMax << " depth " << depth
78  << " R " << rmin << ":" << rmax << " Phi " << etabin.phis[j].first << ":" << etabin.phis[j].second
79  << ":" << dphi << " layer[" << k << "]: " << etabin.layer[k].first - 1 << ":"
80  << etabin.layer[k].second << std::endl;
81 #endif
83  iring, depth, etabin.phis[j].first, etabin.phis[j].second, dphi, rmin, rmax, etabin.etaMin, etabin.etaMax));
84  }
85  depth++;
86  }
87  }
88  return result;
89 }
90 
91 // ----------> HO <-----------
92 std::vector<HcalFlexiHardcodeGeometryLoader::HBHOCellParameters> HcalFlexiHardcodeGeometryLoader::makeHOCells() {
93  const double HORMIN0 = 390.0;
94  const double HORMIN1 = 412.6;
95  const double HORMAX = 413.6;
96  const int nCells = 15;
97  const int nPhi = 72;
98  const double etamin[nCells] = {
99  0.000, 0.087, 0.174, 0.261, 0.3395, 0.435, 0.522, 0.609, 0.696, 0.783, 0.873, 0.957, 1.044, 1.131, 1.218};
100  const double etamax[nCells] = {
101  0.087, 0.174, 0.261, 0.3075, 0.435, 0.522, 0.609, 0.696, 0.783, 0.8494, 0.957, 1.044, 1.131, 1.218, 1.305};
102  std::vector<HcalFlexiHardcodeGeometryLoader::HBHOCellParameters> result;
103  result.reserve(nCells * nPhi);
104  double dphi = ((2.0 * M_PI) / nPhi);
105  for (int i = 0; i < nCells; ++i) {
106  double rmin = ((i < 4) ? HORMIN0 : HORMIN1);
107  for (int iside = -1; iside <= 1; iside += 2) {
108  for (int j = 0; j < nPhi; ++j) {
109  double phi = (j + 0.5) * dphi;
110  // eta, depth, phi, phi0, deltaPhi, rMin, rMax, etaMin, etaMax
112  iside * (i + 1), 4, j + 1, phi, dphi, rmin, HORMAX, etamin[i], etamax[i]));
113  }
114  }
115  }
116  return result;
117 }
118 
119 //
120 // Convert constants to appropriate cells
121 //
123  HcalGeometry* fGeometry, const std::vector<HcalFlexiHardcodeGeometryLoader::HBHOCellParameters>& fCells, bool fHB) {
124  fGeometry->increaseReserve(fCells.size());
125  for (const auto& param : fCells) {
126  HcalDetId hid(fHB ? HcalBarrel : HcalOuter, param.ieta, param.iphi, param.depth);
127  float phiCenter = param.phi; // middle of the cell
128  float etaCenter = 0.5 * (param.etaMin + param.etaMax);
129  float x = param.rMin * cos(phiCenter);
130  float y = param.rMin * sin(phiCenter);
131  float z = (param.ieta < 0) ? -(param.rMin * sinh(etaCenter)) : (param.rMin * sinh(etaCenter));
132  // make cell geometry
133  GlobalPoint refPoint(x, y, z); // center of the cell's face
134  std::vector<CCGFloat> cellParams;
135  cellParams.reserve(5);
136  cellParams.emplace_back(0.5 * (param.etaMax - param.etaMin)); // deta_half
137  cellParams.emplace_back(0.5 * param.dphi); // dphi_half
138  cellParams.emplace_back(0.5 * (param.rMax - param.rMin) * cosh(etaCenter)); // dr_half
139  cellParams.emplace_back(fabs(refPoint.eta()));
140  cellParams.emplace_back(fabs(refPoint.z()));
141 #ifdef EDM_ML_DEBUG
142  std::cout << "HcalFlexiHardcodeGeometryLoader::fillHBHO-> " << hid << " " << hid.rawId() << " " << std::hex
143  << hid.rawId() << std::dec << " " << hid << " " << refPoint << '/' << cellParams[0] << '/'
144  << cellParams[1] << '/' << cellParams[2] << std::endl;
145 #endif
146  fGeometry->newCellFast(refPoint,
147  refPoint,
148  refPoint,
149  CaloCellGeometry::getParmPtr(cellParams, fGeometry->parMgr(), fGeometry->parVecVec()),
150  hid);
151  }
152 }
153 
154 // ----------> HE <-----------
155 std::vector<HcalFlexiHardcodeGeometryLoader::HECellParameters> HcalFlexiHardcodeGeometryLoader::makeHECells(
156  const HcalDDDRecConstants& hcons) {
157  std::vector<HcalFlexiHardcodeGeometryLoader::HECellParameters> result;
158  std::vector<std::pair<double, double> > gconsHE = hcons.getConstHBHE(1);
159 #ifdef EDM_ML_DEBUG
160  std::cout << "HcalFlexiHardcodeGeometryLoader:HE with " << gconsHE.size() << " cells" << std::endl;
161 #endif
162  if (!gconsHE.empty()) {
163  std::vector<HcalDDDRecConstants::HcalEtaBin> etabins = hcons.getEtaBins(1);
164 
165 #ifdef EDM_ML_DEBUG
166  std::cout << "FlexiGeometryLoader called for HE with " << etabins.size() << " Eta Bins and " << gconsHE.size()
167  << " depths" << std::endl;
168  for (unsigned int i = 0; i < gconsHE.size(); ++i)
169  std::cout << " Depth[" << i << "] = " << gconsHE[i].first << " +- " << gconsHE[i].second;
170  std::cout << std::endl;
171 #endif
172  for (auto& etabin : etabins) {
173  int iring = (etabin.zside >= 0) ? etabin.ieta : -etabin.ieta;
174  int depth = etabin.depthStart;
175  double dphi =
176  (etabin.phis.size() > 1) ? (etabin.phis[1].second - etabin.phis[0].second) : ((2.0 * M_PI) / MAX_HCAL_PHI);
177 #ifdef EDM_ML_DEBUG
178  std::cout << "FlexiGeometryLoader::Ring " << iring << " nphi " << etabin.phis.size() << " dstart " << depth
179  << " dphi " << dphi << " layers " << etabin.layer.size() << std::endl;
180  for (unsigned int j = 0; j < etabin.phis.size(); ++j)
181  std::cout << " [" << j << "] " << etabin.phis[j].first << ":" << etabin.phis[j].second;
182  std::cout << std::endl;
183 #endif
184  for (unsigned int k = 0; k < etabin.layer.size(); ++k) {
185  int layf = etabin.layer[k].first - 1;
186  int layl = etabin.layer[k].second - 1;
187  double zmin = gconsHE[layf].first - gconsHE[layf].second;
188  double zmax = gconsHE[layl].first + gconsHE[layl].second;
189  if (zmin < 1.0) {
190  for (int k2 = layf; k2 <= layl; ++k2) {
191  if (gconsHE[k2].first > 10) {
192  zmin = gconsHE[k2].first - gconsHE[k2].second;
193  break;
194  }
195  }
196  }
197  if (zmin >= zmax)
198  zmax = zmin + 10.;
199  for (unsigned int j = 0; j < etabin.phis.size(); ++j) {
200 #ifdef EDM_ML_DEBUG
201  std::cout << "HERing " << iring << " eta " << etabin.etaMin << ":" << etabin.etaMax << " depth " << depth
202  << " Z " << zmin << ":" << zmax << " Phi :" << etabin.phis[j].first << ":" << etabin.phis[j].second
203  << ":" << dphi << " layer[" << k << "]: " << etabin.layer[k].first - 1 << ":"
204  << etabin.layer[k].second - 1 << std::endl;
205 #endif
207  iring, depth, etabin.phis[j].first, etabin.phis[j].second, dphi, zmin, zmax, etabin.etaMin, etabin.etaMax));
208  }
209  depth++;
210  }
211  }
212  }
213  return result;
214 }
215 
216 // ----------> HE @ H2 <-----------
217 std::vector<HcalFlexiHardcodeGeometryLoader::HECellParameters> HcalFlexiHardcodeGeometryLoader::makeHECells_H2() {
218  const double HEZMIN_H2 = 400.715;
219  const double HEZMID_H2 = 436.285;
220  const double HEZMAX_H2 = 541.885;
221  const int nEtas = 10;
222  const int nDepth[nEtas] = {1, 2, 2, 2, 2, 2, 2, 2, 3, 3};
223  const int dStart[nEtas] = {3, 1, 1, 1, 1, 1, 1, 1, 1, 1};
224  const int nPhis[nEtas] = {8, 8, 8, 8, 8, 8, 4, 4, 4, 4};
225  const double etas[nEtas + 1] = {1.305, 1.373, 1.444, 1.521, 1.603, 1.693, 1.790, 1.880, 1.980, 2.090, 2.210};
226  const double zval[4 * nEtas] = {
227  409.885, 462.685, 0., 0., HEZMIN_H2, 427.485, 506.685, 0.0, HEZMIN_H2, HEZMID_H2,
228  524.285, 0., HEZMIN_H2, HEZMID_H2, HEZMAX_H2, 0., HEZMIN_H2, HEZMID_H2, HEZMAX_H2, 0.,
229  HEZMIN_H2, HEZMID_H2, HEZMAX_H2, 0., HEZMIN_H2, HEZMID_H2, HEZMAX_H2, 0., HEZMIN_H2, HEZMID_H2,
230  HEZMAX_H2, 0., HEZMIN_H2, 418.685, HEZMID_H2, HEZMAX_H2, HEZMIN_H2, 418.685, HEZMID_H2, HEZMAX_H2};
231  std::vector<HcalFlexiHardcodeGeometryLoader::HECellParameters> result;
232 
233  for (int i = 0; i < nEtas; ++i) {
234  int ieta = 16 + i;
235  for (int k = 0; k < nDepth[i]; ++k) {
236  int depth = dStart[i] + k;
237  for (int j = 0; j < nPhis[i]; ++j) {
238  int iphi = (nPhis[i] == 8) ? (j + 1) : (2 * j + 1);
239  double dphi = (40.0 * DEGREE2RAD) / nPhis[i];
240  double phi0 = (j + 0.5) * dphi;
241  // ieta, depth, iphi, phi0, deltaPhi, zMin, zMax, etaMin, etaMax
243  ieta, depth, iphi, phi0, dphi, zval[4 * i + k + 1], zval[4 * i + k + 2], etas[i], etas[i + 1]));
244  }
245  }
246  }
247  return result;
248 }
249 
250 // ----------> HF <-----------
251 std::vector<HcalFlexiHardcodeGeometryLoader::HFCellParameters> HcalFlexiHardcodeGeometryLoader::makeHFCells(
252  const HcalDDDRecConstants& hcons) {
253  const float HFZMIN1 = 1115.;
254  const float HFZMIN2 = 1137.;
255  const float HFZMAX = 1280.1;
256  std::vector<HcalDDDRecConstants::HFCellParameters> cells = hcons.getHFCellParameters();
257  unsigned int nCells = cells.size();
258  std::vector<HcalFlexiHardcodeGeometryLoader::HFCellParameters> result;
259  result.reserve(nCells);
260  for (unsigned int i = 0; i < nCells; ++i) {
262  cells[i].depth,
263  cells[i].firstPhi,
264  cells[i].stepPhi,
265  cells[i].nPhi,
266  5 * cells[i].stepPhi,
267  HFZMIN1,
268  HFZMAX,
269  cells[i].rMin,
270  cells[i].rMax);
271  result.emplace_back(cell1);
273  1 + cells[i].depth,
274  cells[i].firstPhi,
275  cells[i].stepPhi,
276  cells[i].nPhi,
277  5 * cells[i].stepPhi,
278  HFZMIN2,
279  HFZMAX,
280  cells[i].rMin,
281  cells[i].rMax);
282  result.emplace_back(cell2);
283  }
284  return result;
285 }
286 
288  HcalGeometry* fGeometry, const std::vector<HcalFlexiHardcodeGeometryLoader::HECellParameters>& fCells) {
289  fGeometry->increaseReserve(fCells.size());
290  for (const auto& param : fCells) {
291  HcalDetId hid(HcalEndcap, param.ieta, param.iphi, param.depth);
292  float phiCenter = param.phi; // middle of the cell
293  float etaCenter = 0.5 * (param.etaMin + param.etaMax);
294  int iside = (param.ieta >= 0) ? 1 : -1;
295  float z = (isBH_) ? (iside * 0.5 * (param.zMin + param.zMax)) : (iside * param.zMin);
296  float perp = std::abs(z) / sinh(etaCenter);
297  float x = perp * cos(phiCenter);
298  float y = perp * sin(phiCenter);
299  // make cell geometry
300  GlobalPoint refPoint(x, y, z); // center of the cell's face
301  std::vector<CCGFloat> cellParams;
302  cellParams.reserve(5);
303  cellParams.emplace_back(0.5 * (param.etaMax - param.etaMin)); //deta_half
304  cellParams.emplace_back(0.5 * param.dphi); // dphi_half
305  cellParams.emplace_back(-0.5 * (param.zMax - param.zMin) / tanh(etaCenter)); // dz_half, "-" means edges in Z
306  cellParams.emplace_back(fabs(refPoint.eta()));
307  cellParams.emplace_back(fabs(refPoint.z()));
308 #ifdef EDM_ML_DEBUG
309  std::cout << "HcalFlexiHardcodeGeometryLoader::fillHE-> " << hid << " " << hid.rawId() << " " << std::hex
310  << hid.rawId() << std::dec << " " << hid << refPoint << '/' << cellParams[0] << '/' << cellParams[1]
311  << '/' << cellParams[2] << std::endl;
312 #endif
313  fGeometry->newCellFast(refPoint,
314  refPoint,
315  refPoint,
316  CaloCellGeometry::getParmPtr(cellParams, fGeometry->parMgr(), fGeometry->parVecVec()),
317  hid);
318  }
319 }
320 
322  HcalGeometry* fGeometry, const std::vector<HcalFlexiHardcodeGeometryLoader::HFCellParameters>& fCells) {
323  fGeometry->increaseReserve(fCells.size());
324  for (const auto& param : fCells) {
325  for (int kPhi = 0; kPhi < param.nPhi; ++kPhi) {
326  int iPhi = param.phiFirst + kPhi * param.phiStep;
327  HcalDetId hid(HcalForward, param.eta, iPhi, param.depth);
328  // middle of the cell
329  float phiCenter = ((iPhi - 1) * 360. / MAX_HCAL_PHI + 0.5 * param.dphi) * DEGREE2RAD;
330  GlobalPoint inner(param.rMin, 0, param.zMin);
331  GlobalPoint outer(param.rMax, 0, param.zMin);
332  float iEta = inner.eta();
333  float oEta = outer.eta();
334  float etaCenter = 0.5 * (iEta + oEta);
335 
336  float perp = param.zMin / sinh(etaCenter);
337  float x = perp * cos(phiCenter);
338  float y = perp * sin(phiCenter);
339  float z = (param.eta > 0) ? param.zMin : -param.zMin;
340  // make cell geometry
341  GlobalPoint refPoint(x, y, z); // center of the cell's face
342  std::vector<CCGFloat> cellParams;
343  cellParams.reserve(5);
344  cellParams.emplace_back(0.5 * (iEta - oEta)); // deta_half
345  cellParams.emplace_back(0.5 * param.dphi * DEGREE2RAD); // dphi_half
346  cellParams.emplace_back(0.5 * (param.zMax - param.zMin)); // dz_half
347  cellParams.emplace_back(fabs(refPoint.eta()));
348  cellParams.emplace_back(fabs(refPoint.z()));
349 #ifdef EDM_ML_DEBUG
350  std::cout << "HcalFlexiHardcodeGeometryLoader::fillHF-> " << hid << " " << hid.rawId() << " " << std::hex
351  << hid.rawId() << std::dec << " " << hid << " " << refPoint << '/' << cellParams[0] << '/'
352  << cellParams[1] << '/' << cellParams[2] << std::endl;
353 #endif
354  fGeometry->newCellFast(refPoint,
355  refPoint,
356  refPoint,
357  CaloCellGeometry::getParmPtr(cellParams, fGeometry->parMgr(), fGeometry->parVecVec()),
358  hid);
359  }
360  }
361 }
IdealObliquePrism.h
HcalGeometry::increaseReserve
void increaseReserve(unsigned int extra)
Definition: HcalGeometry.cc:589
DDAxes::y
photonAnalyzer_cfi.rMax
rMax
Definition: photonAnalyzer_cfi.py:91
mps_fire.i
i
Definition: mps_fire.py:428
HcalFlexiHardcodeGeometryLoader::fillHE
void fillHE(HcalGeometry *fGeometry, const std::vector< HECellParameters > &fCells)
Definition: HcalFlexiHardcodeGeometryLoader.cc:287
CaloCellGeometry::CCGFloat
float CCGFloat
Definition: CaloCellGeometry.h:52
MessageLogger.h
gather_cfg.cout
cout
Definition: gather_cfg.py:144
CCGFloat
CaloCellGeometry::CCGFloat CCGFloat
Definition: CaloCellGeometry.cc:4
HcalTopology
Definition: HcalTopology.h:26
HcalFlexiHardcodeGeometryLoader::makeHECells_H2
std::vector< HECellParameters > makeHECells_H2()
Definition: HcalFlexiHardcodeGeometryLoader.cc:217
cells
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ cells
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:222
HcalFlexiHardcodeGeometryLoader::HcalFlexiHardcodeGeometryLoader
HcalFlexiHardcodeGeometryLoader()
Definition: HcalFlexiHardcodeGeometryLoader.cc:19
HcalFlexiHardcodeGeometryLoader::fillHF
void fillHF(HcalGeometry *fGeometry, const std::vector< HFCellParameters > &fCells)
Definition: HcalFlexiHardcodeGeometryLoader.cc:321
HcalGeometry::newCellFast
void newCellFast(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId)
Definition: HcalGeometry.cc:485
HcalFlexiHardcodeGeometryLoader::isBH_
bool isBH_
Definition: HcalFlexiHardcodeGeometryLoader.h:133
DDAxes::x
HcalBarrel
Definition: HcalAssistant.h:33
CaloSubdetectorGeometry::allocatePar
void allocatePar(ParVec::size_type n, unsigned int m)
Definition: CaloSubdetectorGeometry.cc:113
perp
T perp() const
Magnitude of transverse component.
Definition: Basic3DVectorLD.h:133
nCells
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ const uint32_t *__restrict__ nCells
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
SiStripMonitorCluster_cfi.zmin
zmin
Definition: SiStripMonitorCluster_cfi.py:200
HcalDDDRecConstants::getHFCellParameters
std::vector< HFCellParameters > getHFCellParameters() const
Definition: HcalDDDRecConstants.cc:210
HcalFlexiHardcodeGeometryLoader::makeHFCells
std::vector< HFCellParameters > makeHFCells(const HcalDDDRecConstants &hcons)
Definition: HcalFlexiHardcodeGeometryLoader.cc:251
HcalGeometry.h
muonTiming_cfi.etamin
etamin
Definition: muonTiming_cfi.py:30
HcalFlexiHardcodeGeometryLoader::load
CaloSubdetectorGeometry * load(const HcalTopology &fTopology, const HcalDDDRecConstants &hcons)
Definition: HcalFlexiHardcodeGeometryLoader.cc:24
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
HcalGeometry::sortValidIds
void sortValidIds()
Definition: HcalGeometry.cc:591
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
LEDCalibrationChannels.iphi
iphi
Definition: LEDCalibrationChannels.py:64
HcalFlexiHardcodeGeometryLoader::makeHOCells
std::vector< HBHOCellParameters > makeHOCells()
Definition: HcalFlexiHardcodeGeometryLoader.cc:92
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
SiStripMonitorCluster_cfi.zmax
zmax
Definition: SiStripMonitorCluster_cfi.py:201
HcalFlexiHardcodeGeometryLoader::makeHBCells
std::vector< HBHOCellParameters > makeHBCells(const HcalDDDRecConstants &hcons)
Definition: HcalFlexiHardcodeGeometryLoader.cc:53
EnergyCorrector.etas
etas
Definition: EnergyCorrector.py:45
SurfaceOrientation::inner
Definition: Surface.h:19
IdealZPrism.h
DDAxes::z
HcalOuter
Definition: HcalAssistant.h:35
dqmdumpme.k
k
Definition: dqmdumpme.py:60
Point3DBase< float, GlobalTag >
LEDCalibrationChannels.depth
depth
Definition: LEDCalibrationChannels.py:65
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:112
HcalTopology::mode
HcalTopologyMode::Mode mode() const
Definition: HcalTopology.h:34
HcalTopologyMode::H2
Definition: HcalTopologyMode.h:26
HcalTopology::getHFSize
unsigned int getHFSize() const
Definition: HcalTopology.h:135
LEDCalibrationChannels.ieta
ieta
Definition: LEDCalibrationChannels.py:63
JetComb::kPhi
Definition: TtSemiLepJetComb.h:17
HcalDDDRecConstants::getConstHBHE
std::vector< std::pair< double, double > > getConstHBHE(const int &type) const
Definition: HcalDDDRecConstants.h:54
PV3DBase::eta
T eta() const
Definition: PV3DBase.h:73
HLT_FULL_cff.nPhi
nPhi
Definition: HLT_FULL_cff.py:6587
HcalDetId
Definition: HcalDetId.h:12
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
hgcalTopologyTester_cfi.cell2
cell2
Definition: hgcalTopologyTester_cfi.py:12
HcalFlexiHardcodeGeometryLoader::DEGREE2RAD
double DEGREE2RAD
Definition: HcalFlexiHardcodeGeometryLoader.h:132
HcalGeometry::k_NumberOfParametersPerShape
Definition: HcalGeometry.h:38
HcalFlexiHardcodeGeometryLoader::makeHECells
std::vector< HECellParameters > makeHECells(const HcalDDDRecConstants &hcons)
Definition: HcalFlexiHardcodeGeometryLoader.cc:155
cms::cuda::for
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
Definition: HistoContainer.h:27
hgcalTopologyTester_cfi.cell1
cell1
Definition: hgcalTopologyTester_cfi.py:11
HcalFlexiHardcodeGeometryLoader::HFCellParameters
Definition: HcalFlexiHardcodeGeometryLoader.h:87
HcalFlexiHardcodeGeometryLoader::HECellParameters
Definition: HcalFlexiHardcodeGeometryLoader.h:56
HcalFlexiHardcodeGeometryLoader::HBHOCellParameters
Definition: HcalFlexiHardcodeGeometryLoader.h:25
HcalForward
Definition: HcalAssistant.h:36
DDAxes::phi
HcalTopology.h
CaloCellGeometry.h
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
CaloSubdetectorGeometry::cornersMgr
CaloCellGeometry::CornersMgr * cornersMgr()
Definition: CaloSubdetectorGeometry.h:82
HcalEndcap
Definition: HcalAssistant.h:34
HcalDDDRecConstants::getEtaBins
std::vector< HcalEtaBin > getEtaBins(const int &itype) const
Definition: HcalDDDRecConstants.cc:68
HcalFlexiHardcodeGeometryLoader::MAX_HCAL_PHI
int MAX_HCAL_PHI
Definition: HcalFlexiHardcodeGeometryLoader.h:131
muonTiming_cfi.etamax
etamax
Definition: muonTiming_cfi.py:23
HcalDDDRecConstants.h
HcalTopology::ncells
unsigned int ncells() const override
return a count of valid cells (for dense indexing use)
Definition: HcalTopology.cc:1687
CaloCellGeometry::getParmPtr
static const CCGFloat * getParmPtr(const std::vector< CCGFloat > &vd, ParMgr *mgr, ParVecVec &pvv)
Definition: CaloCellGeometry.cc:117
CaloSubdetectorGeometry
Definition: CaloSubdetectorGeometry.h:22
HcalFlexiHardcodeGeometryLoader.h
CaloSubdetectorGeometry::parMgr
ParMgr * parMgr()
Definition: CaloSubdetectorGeometry.h:86
HcalDDDRecConstants
Definition: HcalDDDRecConstants.h:23
HcalFlexiHardcodeGeometryLoader::fillHBHO
void fillHBHO(HcalGeometry *fGeometry, const std::vector< HBHOCellParameters > &fCells, bool fHB)
Definition: HcalFlexiHardcodeGeometryLoader.cc:122
photonAnalyzer_cfi.rMin
rMin
Definition: photonAnalyzer_cfi.py:90
L1TowerCalibrationProducer_cfi.iEta
iEta
Definition: L1TowerCalibrationProducer_cfi.py:60
mps_fire.result
result
Definition: mps_fire.py:311
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CaloSubdetectorGeometry::allocateCorners
void allocateCorners(CaloCellGeometry::CornersVec::size_type n)
Definition: CaloSubdetectorGeometry.cc:106
SurfaceOrientation::outer
Definition: Surface.h:19
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
CaloSubdetectorGeometry::parVecVec
ParVecVec & parVecVec()
Definition: CaloSubdetectorGeometry.h:89
HcalGeometry
Definition: HcalGeometry.h:17
TauDecayModes.dec
dec
Definition: TauDecayModes.py:143
defaults_cfi.etabins
etabins
Definition: defaults_cfi.py:7
HcalGeometry::numberOfShapes
unsigned int numberOfShapes() const override
Definition: HcalGeometry.h:42
HcalDDDRecConstants::isBH
bool isBH() const
Definition: HcalDDDRecConstants.h:109