CMS 3D CMS Logo

HcalGeometry.cc
Go to the documentation of this file.
3 #include <algorithm>
4 #include <iostream>
5 
6 #include <Math/Transform3D.h>
7 #include <Math/EulerAngles.h>
8 
13 
14 //#define EDM_ML_DEBUG
15 
17  m_topology(topology), m_mergePosition(topology.getMergePositionFlag()) {
18  init();
19 }
20 
22 
25 #ifdef EDM_ML_DEBUG
26  std::cout << "HcalGeometry: " << "HcalGeometry::init() "
27  << " HBSize " << m_topology.getHBSize()
28  << " HESize " << m_topology.getHESize()
29  << " HOSize " << m_topology.getHOSize()
30  << " HFSize " << m_topology.getHFSize() << std::endl;
31 #endif
32 
37 }
38 
40  // this must test the last record filled to avoid a race condition
41  if(!m_emptyIds.isSet()) {
42  std::unique_ptr<std::vector<DetId>> p_hbIds{new std::vector<DetId>};
43  std::unique_ptr<std::vector<DetId>> p_heIds{new std::vector<DetId>};
44  std::unique_ptr<std::vector<DetId>> p_hoIds{new std::vector<DetId>};
45  std::unique_ptr<std::vector<DetId>> p_hfIds{new std::vector<DetId>};
46  std::unique_ptr<std::vector<DetId>> p_emptyIds{new std::vector<DetId>};
47 
48  const std::vector<DetId>& baseIds (CaloSubdetectorGeometry::getValidDetIds());
49  for (unsigned int i ( 0 ) ; i != baseIds.size() ; ++i) {
50  const DetId id ( baseIds[i] );
51  if (id.subdetId() == HcalBarrel) {
52  p_hbIds->emplace_back( id ) ;
53  } else if (id.subdetId() == HcalEndcap) {
54  p_heIds->emplace_back( id ) ;
55  } else if (id.subdetId() == HcalOuter) {
56  p_hoIds->emplace_back( id ) ;
57  } else if (id.subdetId() == HcalForward) {
58  p_hfIds->emplace_back( id ) ;
59  }
60  }
61  std::sort( p_hbIds->begin(), p_hbIds->end() ) ;
62  std::sort( p_heIds->begin(), p_heIds->end() ) ;
63  std::sort( p_hoIds->begin(), p_hoIds->end() ) ;
64  std::sort( p_hfIds->begin(), p_hfIds->end() ) ;
65  p_emptyIds->resize( 0 ) ;
66 
67  m_hbIds.set(std::move(p_hbIds));
68  m_heIds.set(std::move(p_heIds));
69  m_hoIds.set(std::move(p_hoIds));
70  m_hfIds.set(std::move(p_hfIds));
71  m_emptyIds.set(std::move(p_emptyIds));
72  }
73 }
74 
75 const std::vector<DetId>&
77  int subdet ) const {
78  if( 0 != subdet ) fillDetIds() ;
79  return ( 0 == subdet ? CaloSubdetectorGeometry::getValidDetIds() :
80  ( HcalBarrel == subdet ? *m_hbIds.load() :
81  ( HcalEndcap == subdet ? *m_heIds.load() :
82  ( HcalOuter == subdet ? *m_hoIds.load() :
83  ( HcalForward == subdet ? *m_hfIds.load() : *m_emptyIds.load() ) ) ) ) ) ;
84 }
85 
86 std::shared_ptr<const CaloCellGeometry> HcalGeometry::getGeometry(const DetId& id) const {
87 #ifdef EDM_ML_DEBUG
88  std::cout << "HcalGeometry::getGeometry for " << HcalDetId(id) << " "
89  << m_mergePosition << " ";
90 #endif
91  if (!m_mergePosition) {
92 #ifdef EDM_ML_DEBUG
94  << "\n";
95 #endif
96  return getGeometryBase(id);
97  } else {
98 #ifdef EDM_ML_DEBUG
100  << getGeometryBase(m_topology.idFront(id)) << "\n";
101 #endif
102  return getGeometryBase(m_topology.idFront(id));
103  }
104 }
105 
107 
108  // Now find the closest eta_bin, eta value of a bin i is average
109  // of eta[i] and eta[i-1]
110  static const double z_long=1100.0;
111  double abseta = fabs(r.eta());
112  double absz = fabs(r.z());
113 
114  // figure out subdetector, giving preference to HE in HE/HF overlap region
116  if (abseta <= m_topology.etaMax(HcalBarrel) ) {
117  bc = HcalBarrel;
118  } else if (absz >= z_long) {
119  bc = HcalForward;
120  } else if (abseta <= m_topology.etaMax(HcalEndcap) ) {
121  bc = HcalEndcap;
122  } else {
123  bc = HcalForward;
124  }
125 
126  // find eta bin
127  int etaring = etaRing(bc, abseta);
128 
129  int phibin = phiBin(bc, etaring, r.phi());
130 
131  // add a sign to the etaring
132  int etabin = (r.z() > 0) ? etaring : -etaring;
133 
134  if (bc == HcalForward) {
135  static const double z_short=1137.0;
136  // Next line is premature depth 1 and 2 can coexist for large z-extent
137  // HcalDetId bestId(bc,etabin,phibin,((fabs(r.z())>=z_short)?(2):(1)));
138  // above line is no good with finite precision
139  HcalDetId bestId(bc,etabin,phibin,((fabs(r.z()) - z_short >-0.1)?(2):(1)));
140  return correctId(bestId);
141  } else {
142 
143  //Now do depth if required
144  int zside = (r.z() > 0) ? 1 : -1;
145  int dbin = m_topology.dddConstants()->getMinDepth(((int)(bc)-1),etaring,phibin,zside);
146  double pointrz(0), drz(99999.);
147  HcalDetId currentId(bc, etabin, phibin, dbin);
148  if (bc == HcalBarrel) pointrz = r.mag();
149  else pointrz = std::abs(r.z());
150  HcalDetId bestId;
151  for ( ; currentId != HcalDetId(); m_topology.incrementDepth(currentId)) {
152  std::shared_ptr<const CaloCellGeometry> cell = getGeometry(currentId);
153  if (cell == nullptr) {
154  assert (bestId != HcalDetId());
155  break;
156  } else {
157  double rz;
158  if (bc == HcalEndcap) rz = std::abs(cell->getPosition().z());
159  else rz = cell->getPosition().mag();
160  if (std::abs(pointrz-rz)<drz) {
161  bestId = currentId;
162  drz = std::abs(pointrz-rz);
163  }
164  }
165  }
166 
167  return correctId(bestId);
168  }
169 }
170 
172  if (!m_mergePosition) {
173  return (getGeometryBase(id)->getPosition());
174  } else {
176  }
177 }
178 
180  if (!m_mergePosition) {
181  return (getGeometryBase(id)->getBackPoint());
182  } else {
183  std::vector<HcalDetId> ids;
185  return (getGeometryBase(ids.back())->getBackPoint());
186  }
187 }
188 
190  if (!m_mergePosition) {
191  return (getGeometryBase(id)->getCorners());
192  } else {
193  std::vector<HcalDetId> ids;
198  for (unsigned int k=0; k<4; ++k) {
199  mcorners[k] = mcf[k];
200  mcorners[k+4] = mcb[k+4];
201  }
202  return mcorners;
203  }
204 }
205 
206 int HcalGeometry::etaRing(HcalSubdetector bc, double abseta) const {
207  return m_topology.etaRing(bc, abseta);
208 }
209 
210 int HcalGeometry::phiBin(HcalSubdetector bc, int etaring, double phi) const {
211  return m_topology.phiBin(bc, etaring, phi);
212 }
213 
215  double dR ) const {
216  CaloSubdetectorGeometry::DetIdSet dis; // this is the return object
217 
218  if (0.000001 < dR) {
219  if (dR > M_PI/2.) {// this version needs "small" dR
220  dis = CaloSubdetectorGeometry::getCells(r, dR); // base class version
221  } else {
222  const double dR2 ( dR*dR ) ;
223  const double reta ( r.eta() ) ;
224  const double rphi ( r.phi() ) ;
225  const double lowEta ( reta - dR ) ;
226  const double highEta ( reta + dR ) ;
227  const double lowPhi ( rphi - dR ) ;
228  const double highPhi ( rphi + dR ) ;
229 
230  const double hfEtaHi (m_topology.etaMax(HcalForward));
231 
232  if (highEta > -hfEtaHi &&
233  lowEta < hfEtaHi ) { // in hcal
235 
236  for (unsigned int is ( 0 ) ; is != 4 ; ++is ) {
237  const int sign ( reta>0 ? 1 : -1 ) ;
238  const int ieta_center ( sign*etaRing( hs[is], fabs( reta ) ) ) ;
239  const int ieta_lo ( ( 0 < lowEta*sign ? sign : -sign )*etaRing( hs[is], fabs( lowEta ) ) ) ;
240  const int ieta_hi ( ( 0 < highEta*sign ? sign : -sign )*etaRing( hs[is], fabs( highEta ) ) ) ;
241  const int iphi_lo ( phiBin( hs[is], ieta_center, lowPhi ) ) ;
242  const int iphi_hi ( phiBin( hs[is], ieta_center, highPhi ) ) ;
243  const int jphi_lo ( iphi_lo>iphi_hi ? iphi_lo - 72 : iphi_lo ) ;
244  const int jphi_hi ( iphi_hi ) ;
245 
246  const int idep_lo ( 1 == is ? 4 : 1 ) ;
247  const int idep_hi ( m_topology.maxDepth(hs[is]) );
248  for (int ieta ( ieta_lo ) ; ieta <= ieta_hi ; ++ieta) {// over eta limits
249  if (ieta != 0) {
250  for (int jphi ( jphi_lo ) ; jphi <= jphi_hi ; ++jphi) { // over phi limits
251  const int iphi ( 1 > jphi ? jphi+72 : jphi ) ;
252  for (int idep ( idep_lo ) ; idep <= idep_hi ; ++idep ) {
253  const HcalDetId did ( hs[is], ieta, iphi, idep ) ;
254  if (m_topology.valid(did)) {
255  std::shared_ptr<const CaloCellGeometry> cell ( getGeometryBase( did ) );
256  if (nullptr != cell ) {
257  const GlobalPoint& p ( cell->getPosition() ) ;
258  const double eta ( p.eta() ) ;
259  const double phi ( p.phi() ) ;
260  if (reco::deltaR2(eta, phi, reta, rphi ) < dR2) dis.insert( did ) ;
261  }
262  }
263  }
264  }
265  }
266  }
267  }
268  }
269  }
270  }
271  return dis;
272 }
273 
274 unsigned int HcalGeometry::getHxSize(const int type) const {
275  unsigned int hxsize(0);
276  if (!m_hbIds.isSet()) fillDetIds();
277  if (type == 1) hxsize = (m_hbIds.isSet()) ? m_hbIds->size() : 0;
278  else if (type == 2) hxsize = (m_heIds.isSet()) ? m_heIds->size() : 0;
279  else if (type == 3) hxsize = (m_hoIds.isSet()) ? m_hoIds->size() : 0;
280  else if (type == 4) hxsize = (m_hfIds.isSet()) ? m_hfIds->size() : 0;
281  else if (type == 0) hxsize = (m_emptyIds.isSet()) ? m_emptyIds->size() : 0;
282  return hxsize;
283 }
284 
286  assert( i < numberOfBarrelAlignments() ) ;
287  const int ieta ( i < numberOfBarrelAlignments()/2 ? -1 : 1 ) ;
288  const int iphi ( 1 + (4*i)%72 ) ;
289  return HcalDetId( HcalBarrel, ieta, iphi, 1 ) ;
290 }
291 
293  assert( i < numberOfEndcapAlignments() ) ;
294  const int ieta ( i < numberOfEndcapAlignments()/2 ? -16 : 16 ) ;
295  const int iphi ( 1 + (4*i)%72 ) ;
296  return HcalDetId( HcalEndcap, ieta, iphi, 1 ) ;
297 }
298 
300  assert( i < numberOfForwardAlignments() ) ;
301  const int ieta ( i < numberOfForwardAlignments()/2 ? -29 : 29 ) ;
302  const int iphi ( 1 + (4*i)%72 ) ;
303  return HcalDetId( HcalForward, ieta, iphi, 1 ) ;
304 }
305 
307  assert( i < numberOfOuterAlignments() ) ;
308  const int ring ( i/12 ) ;
309  const int ieta ( 0 == ring ? -11 :
310  1 == ring ? -5 :
311  2 == ring ? 1 :
312  3 == ring ? 5 : 11 ) ;
313  const int iphi ( 1 + ( i - ring*12 )*6 ) ;
314  return HcalDetId( HcalOuter, ieta, iphi, 4 ) ;
315 }
316 
318  assert( i < numberOfAlignments() ) ;
319 
320  const unsigned int nB ( numberOfBarrelAlignments() ) ;
321  const unsigned int nE ( numberOfEndcapAlignments() ) ;
322  const unsigned int nF ( numberOfForwardAlignments() ) ;
323 // const unsigned int nO ( numberOfOuterAlignments() ) ;
324 
325  return ( i < nB ? detIdFromBarrelAlignmentIndex( i ) :
326  i < nB+nE ? detIdFromEndcapAlignmentIndex( i - nB ) :
327  i < nB+nE+nF ? detIdFromForwardAlignmentIndex( i - nB - nE ) :
328  detIdFromOuterAlignmentIndex( i - nB - nE - nF ) ) ;
329 }
330 
332  unsigned int nD) {
333  const HcalDetId hid ( id ) ;
334  const unsigned int iphi ( hid.iphi() ) ;
335  const int ieta ( hid.ieta() ) ;
336  const unsigned int index ( ( 0 < ieta ? nD/2 : 0 ) + ( iphi + 1 )%72/4 ) ;
337  assert( index < nD ) ;
338  return index ;
339 }
340 
343 }
344 
347 }
348 
351 }
352 
354  const HcalDetId hid ( id ) ;
355  const int ieta ( hid.ieta() ) ;
356  const int iphi ( hid.iphi() ) ;
357  const int ring ( ieta < -10 ? 0 :
358  ( ieta < -4 ? 1 :
359  ( ieta < 5 ? 2 :
360  ( ieta < 11 ? 3 : 4 ) ) ) ) ;
361 
362  const unsigned int index ( 12*ring + ( iphi - 1 )/6 ) ;
363  assert( index < numberOfOuterAlignments() ) ;
364  return index ;
365 }
366 
368  assert(id.det() == DetId::Hcal) ;
369 
370  const HcalDetId hid ( id ) ;
371  bool isHB = (hid.subdet() == HcalBarrel);
372  bool isHE = (hid.subdet() == HcalEndcap);
373  bool isHF = (hid.subdet() == HcalForward);
374  // bool isHO = (hid.subdet() == HcalOuter);
375 
376  const unsigned int nB ( numberOfBarrelAlignments() ) ;
377  const unsigned int nE ( numberOfEndcapAlignments() ) ;
378  const unsigned int nF ( numberOfForwardAlignments() ) ;
379  // const unsigned int nO ( numberOfOuterAlignments() ) ;
380 
381  const unsigned int index (isHB ? alignmentBarrelIndexLocal(id) :
382  isHE ? alignmentEndcapIndexLocal(id) + nB :
383  isHF ? alignmentForwardIndexLocal( id ) + nB + nE :
384  alignmentOuterIndexLocal(id) + nB + nE + nF
385  );
386 
387  assert( index < numberOfAlignments() ) ;
388  return index ;
389 }
390 
392  return (unsigned int)DetId::Hcal - 1 ;
393 }
394 
396  const CCGFloat* pv,
397  unsigned int i,
398  Pt3D& ref) {
400 
401  if (hid.subdet() == HcalForward ) {
402  IdealZPrism::localCorners( lc, pv, ref ) ;
403  } else {
404  IdealObliquePrism::localCorners( lc, pv, ref ) ;
405  }
406 }
407 
409  const GlobalPoint& f2 ,
410  const GlobalPoint& f3 ,
411  const CCGFloat* parm ,
412  const DetId& detId) {
413 
414  assert (detId.det()==DetId::Hcal);
415 
416  const HcalDetId hid ( detId ) ;
417  unsigned int din=m_topology.detId2denseId(detId);
418 
419 #ifdef EDM_ML_DEBUG
420  std::cout << "HcalGeometry: " << " newCell subdet "
421  << detId.subdetId() << ", raw ID "
422  << detId.rawId() << ", hid " << hid << ", din "
423  << din << ", index " << std::endl;
424 #endif
425 
426  if (hid.subdet()==HcalBarrel) {
427  m_hbCellVec.at( din ) = IdealObliquePrism( f1, cornersMgr(), parm ) ;
428  } else if (hid.subdet()==HcalEndcap) {
429  const unsigned int index ( din - m_hbCellVec.size() ) ;
430  m_heCellVec.at( index ) = IdealObliquePrism( f1, cornersMgr(), parm ) ;
431  } else if (hid.subdet()==HcalOuter) {
432  const unsigned int index ( din
433  - m_hbCellVec.size()
434  - m_heCellVec.size() ) ;
435  m_hoCellVec.at( index ) = IdealObliquePrism( f1, cornersMgr(), parm ) ;
436  } else {
437  const unsigned int index ( din
438  - m_hbCellVec.size()
439  - m_heCellVec.size()
440  - m_hoCellVec.size() ) ;
441  m_hfCellVec.at( index ) = IdealZPrism( f1, cornersMgr(), parm, hid.depth()==1 ? IdealZPrism::EM : IdealZPrism::HADR ) ;
442  }
443 
444  return din;
445 }
446 
448  const GlobalPoint& f2 ,
449  const GlobalPoint& f3 ,
450  const CCGFloat* parm ,
451  const DetId& detId) {
452 
453  unsigned int din = newCellImpl(f1,f2,f3,parm,detId);
454 
455  addValidID( detId ) ;
456  m_dins.emplace_back( din );
457 }
458 
460  const GlobalPoint& f2 ,
461  const GlobalPoint& f3 ,
462  const CCGFloat* parm ,
463  const DetId& detId) {
464 
465  unsigned int din = newCellImpl(f1,f2,f3,parm,detId);
466 
467  m_validIds.emplace_back( detId ) ;
468  m_dins.emplace_back( din );
469 }
470 
472  // Modify the RawPtr class
473  const CaloCellGeometry* cell(nullptr);
474  if (m_hbCellVec.size() > din) {
475  cell = (&m_hbCellVec[din]);
476  } else if (m_hbCellVec.size()+m_heCellVec.size() > din) {
477  const unsigned int ind (din - m_hbCellVec.size() ) ;
478  cell = (&m_heCellVec[ind]);
479  } else if (m_hbCellVec.size()+m_heCellVec.size()+m_hoCellVec.size() > din) {
480  const unsigned int ind (din - m_hbCellVec.size() - m_heCellVec.size());
481  cell = (&m_hoCellVec[ind]);
482  } else if (m_hbCellVec.size()+m_heCellVec.size()+m_hoCellVec.size()+
483  m_hfCellVec.size() > din) {
484  const unsigned int ind (din - m_hbCellVec.size() - m_heCellVec.size() -
485  m_hoCellVec.size() ) ;
486  cell = (&m_hfCellVec[ind]);
487  }
488 
489  return (( nullptr == cell || nullptr == cell->param()) ? nullptr : cell ) ;
490 }
491 
495  CaloSubdetectorGeometry::IVec& dinsVec ) const {
496 
497  tVec.reserve( m_topology.ncells()*numberOfTransformParms());
498  iVec.reserve( numberOfShapes()==1 ? 1 : m_topology.ncells());
499  dVec.reserve( numberOfShapes()*numberOfParametersPerShape());
500  dinsVec.reserve( m_topology.ncells());
501 
502  for (const auto & pv : parVecVec()) {
503  for (float iv : pv) {
504  dVec.emplace_back( iv ) ;
505  }
506  }
507 
508  for( auto i : m_dins ) {
509  Tr3D tr ;
510  auto ptr = cellGeomPtr( i );
511 
512  if (nullptr != ptr) {
513  dinsVec.emplace_back( i );
514 
515  const CCGFloat* par ( ptr->param() ) ;
516 
517  unsigned int ishape ( 9999 ) ;
518 
519  for( unsigned int ivv ( 0 ) ; ivv != parVecVec().size() ; ++ivv ) {
520  bool ok ( true ) ;
521  const CCGFloat* pv ( &(*parVecVec()[ivv].begin() ) ) ;
522  for( unsigned int k ( 0 ) ; k != numberOfParametersPerShape() ; ++k ) {
523  ok = ok && ( fabs( par[k] - pv[k] ) < 1.e-6 ) ;
524  }
525  if( ok ) {
526  ishape = ivv;
527  break ;
528  }
529  }
530  assert( 9999 != ishape ) ;
531 
532  const unsigned int nn (( numberOfShapes()==1) ? (unsigned int)1 : m_dins.size() ) ;
533  if( iVec.size() < nn ) iVec.emplace_back( ishape ) ;
534 
535  ptr->getTransform( tr, ( Pt3DVec* ) nullptr ) ;
536 
537  if( Tr3D() == tr ) { // for preshower there is no rotation
538  const GlobalPoint& gp ( ptr->getPosition() ) ;
539  tr = HepGeom::Translate3D( gp.x(), gp.y(), gp.z() ) ;
540  }
541 
542  const CLHEP::Hep3Vector tt ( tr.getTranslation() ) ;
543  tVec.emplace_back( tt.x() ) ;
544  tVec.emplace_back( tt.y() ) ;
545  tVec.emplace_back( tt.z() ) ;
546  if (6 == numberOfTransformParms()) {
547  const CLHEP::HepRotation rr ( tr.getRotation() ) ;
548  const ROOT::Math::Transform3D rtr (rr.xx(), rr.xy(), rr.xz(), tt.x(),
549  rr.yx(), rr.yy(), rr.yz(), tt.y(),
550  rr.zx(), rr.zy(), rr.zz(), tt.z());
552  rtr.GetRotation( ea ) ;
553  tVec.emplace_back( ea.Phi() ) ;
554  tVec.emplace_back( ea.Theta() ) ;
555  tVec.emplace_back( ea.Psi() ) ;
556  }
557  }
558  }
559 }
560 
562 
563  if (m_mergePosition) {
564  HcalDetId hid(id);
565  return ((DetId)(m_topology.mergedDepthDetId(hid)));
566  } else {
567  return id;
568  }
569 }
570 
571 void HcalGeometry::increaseReserve(unsigned int extra) {
572  m_validIds.reserve(m_validIds.size()+extra);
573 }
574 
576  std::sort(m_validIds.begin(),m_validIds.end());
577 }
578 
DetId denseId2detId(unsigned int) const override
return a linear packed id
static unsigned int numberOfBarrelAlignments()
Definition: HcalGeometry.h:68
type
Definition: HCALResponse.h:21
unsigned int getHFSize() const
Definition: HcalTopology.h:141
static unsigned int alignmentTransformIndexLocal(const DetId &id)
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:167
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const override
Get the cell geometry of a given detector id. Should return false if not found.
Definition: HcalGeometry.cc:86
std::vector< CCGFloat > DimVec
unsigned int getHOSize() const
Definition: HcalTopology.h:140
CaloCellGeometry::Pt3DVec Pt3DVec
Definition: HcalGeometry.cc:11
CaloSubdetectorGeometry::DetIdSet getCells(const GlobalPoint &r, double dR) const override
Get a list of all cells within a dR of the given cell.
edm::AtomicPtrCache< std::vector< DetId > > m_heIds
Definition: HcalGeometry.h:163
const CaloCellGeometry * getGeometryRawPtr(uint32_t index) const override
static unsigned int alignmentEndcapIndexLocal(const DetId &id)
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
int phiBin(HcalSubdetector subdet, int etaRing, double phi) const
bool valid(const DetId &id) const override
CaloTopology const * topology(0)
unsigned int getHxSize(const int type) const
static DetId detIdFromOuterAlignmentIndex(unsigned int i)
unsigned int detId2denseId(const DetId &id) const override
return a linear packed id
const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const override
Get a list of valid detector ids (for the given subdetector)
Definition: HcalGeometry.cc:76
unsigned int numberOfShapes() const override
Definition: HcalGeometry.h:45
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
const HcalTopology & m_topology
Definition: HcalGeometry.h:159
std::vector< unsigned int > IVec
edm::AtomicPtrCache< std::vector< DetId > > m_emptyIds
Definition: HcalGeometry.h:166
bool isHE(int etabin, int depth)
std::vector< CCGFloat > TrVec
HepGeom::Transform3D Tr3D
std::vector< IdealObliquePrism > HECellVec
Definition: HcalGeometry.h:27
int zside(DetId const &)
std::vector< Pt3D > Pt3DVec
bool isHB(int etabin, int depth)
CaloCellGeometry::CCGFloat CCGFloat
bool m_mergePosition
Definition: HcalGeometry.h:160
bool getMergePositionFlag() const
Definition: HcalTopology.h:170
bool isSet() const
virtual const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const
Get a list of valid detector ids (for the given subdetector)
HcalDetId mergedDepthDetId(const HcalDetId &id) const
Definition: HcalTopology.h:169
std::shared_ptr< const CaloCellGeometry > getGeometryBase(const DetId &id) const
Definition: HcalGeometry.h:129
CaloCellGeometry::Tr3D Tr3D
Definition: HcalGeometry.cc:12
uint32_t rawId() const
get the raw id
Definition: DetId.h:44
void newCellFast(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId)
CaloCellGeometry::Pt3DVec Pt3DVec
Definition: HcalGeometry.h:33
virtual std::shared_ptr< const CaloCellGeometry > cellGeomPtr(uint32_t index) const
virtual DetIdSet getCells(const GlobalPoint &r, double dR) const
Get a list of all cells within a dR of the given cell.
const CCGFloat * param() const
T mag() const
Definition: PV3DBase.h:67
HBCellVec m_hbCellVec
Definition: HcalGeometry.h:169
HOCellVec m_hoCellVec
Definition: HcalGeometry.h:171
static void localCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
HFCellVec m_hfCellVec
Definition: HcalGeometry.h:172
unsigned int numberOfParametersPerShape() const override
Definition: HcalGeometry.h:46
std::vector< DetId > m_validIds
std::vector< IdealObliquePrism > HOCellVec
Definition: HcalGeometry.h:28
static unsigned int numberOfOuterAlignments()
Definition: HcalGeometry.h:74
int maxDepth(HcalSubdetector subdet) const
int getMinDepth(const int &itype, const int &ieta, const int &iphi, const int &zside) const
void sortValidIds()
int etaRing(HcalSubdetector subdet, double eta) const
eta and phi index from eta, phi values
DetId correctId(const DetId &id) const
T z() const
Definition: PV3DBase.h:64
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
static DetId detIdFromForwardAlignmentIndex(unsigned int i)
static unsigned int alignmentTransformIndexGlobal(const DetId &id)
static unsigned int alignmentBarrelIndexLocal(const DetId &id)
HcalSubdetector
Definition: HcalAssistant.h:31
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static DetId detIdFromBarrelAlignmentIndex(unsigned int i)
void unmergeDepthDetId(const HcalDetId &id, std::vector< HcalDetId > &ids) const
Definition: HcalTopology.h:171
CaloCellGeometry::CCGFloat CCGFloat
Definition: HcalGeometry.cc:9
int phiBin(HcalSubdetector bc, int etaring, double phi) const
static DetId detIdFromLocalAlignmentIndex(unsigned int i)
int etaRing(HcalSubdetector bc, double abseta) const
helper methods for getClosestCell
CaloSubdetectorGeometry::IVec m_dins
Definition: HcalGeometry.h:167
HcalGeometry(const HcalTopology &topology)
Definition: HcalGeometry.cc:16
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:38
GlobalPoint getPosition(const DetId &id) const
void localCorners(Pt3DVec &lc, const CCGFloat *pv, unsigned int i, Pt3D &ref)
#define M_PI
int k[5][pyjets_maxn]
static unsigned int alignmentOuterIndexLocal(const DetId &id)
std::vector< IdealZPrism > HFCellVec
Definition: HcalGeometry.h:29
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:103
~HcalGeometry() override
The HcalGeometry will delete all its cell geometries at destruction time.
Definition: HcalGeometry.cc:21
Definition: DetId.h:18
bool incrementDepth(HcalDetId &id) const
static unsigned int numberOfEndcapAlignments()
Definition: HcalGeometry.h:70
CaloCellGeometry::Pt3D Pt3D
Definition: HcalGeometry.cc:10
AlgebraicVector EulerAngles
Definition: Definitions.h:36
void addValidID(const DetId &id)
bool set(std::unique_ptr< T > iNewValue) const
CaloCellGeometry::CornersMgr * cornersMgr()
Detector
Definition: DetId.h:24
static void localCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
Definition: IdealZPrism.cc:122
CaloCellGeometry::CornersVec getCorners(const DetId &id) const
HepGeom::Point3D< CCGFloat > Pt3D
bool isHF(int etabin, int depth)
std::vector< IdealObliquePrism > HBCellVec
Definition: HcalGeometry.h:26
unsigned int getHESize() const
Definition: HcalTopology.h:139
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
HECellVec m_heCellVec
Definition: HcalGeometry.h:170
CaloCellGeometry::Tr3D Tr3D
void increaseReserve(unsigned int extra)
T eta() const
Definition: PV3DBase.h:76
void getSummary(CaloSubdetectorGeometry::TrVec &trVector, CaloSubdetectorGeometry::IVec &iVector, CaloSubdetectorGeometry::DimVec &dimVector, CaloSubdetectorGeometry::IVec &dinsVector) const override
edm::AtomicPtrCache< std::vector< DetId > > m_hfIds
Definition: HcalGeometry.h:165
#define begin
Definition: vmac.h:32
void fillDetIds() const
Definition: HcalGeometry.cc:39
void newCell(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId) override
static unsigned int numberOfForwardAlignments()
Definition: HcalGeometry.h:72
unsigned int ncells() const override
return a count of valid cells (for dense indexing use)
static unsigned int alignmentForwardIndexLocal(const DetId &id)
GlobalPoint getBackPosition(const DetId &id) const
HcalDetId idFront(const HcalDetId &id) const
Definition: HcalTopology.h:176
static DetId detIdFromEndcapAlignmentIndex(unsigned int i)
static unsigned int alignmentBarEndForIndexLocal(const DetId &id, unsigned int nD)
virtual unsigned int numberOfTransformParms() const
edm::AtomicPtrCache< std::vector< DetId > > m_hoIds
Definition: HcalGeometry.h:164
CaloCellGeometry::Pt3D Pt3D
Definition: HcalGeometry.h:32
Detector det() const
get the detector field from this detid
Definition: DetId.h:36
unsigned int getHBSize() const
Definition: HcalTopology.h:138
unsigned int newCellImpl(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId)
T const * load() const
edm::AtomicPtrCache< std::vector< DetId > > m_hbIds
Definition: HcalGeometry.h:162
static unsigned int numberOfAlignments()
Definition: HcalGeometry.h:78
def move(src, dest)
Definition: eostools.py:510
DetId getClosestCell(const GlobalPoint &r) const override
double etaMax(HcalSubdetector subdet) const