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 
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  return getClosestCell(r,false);
108 }
109 
111  bool ignoreCorrect) const {
112 
113  // Now find the closest eta_bin, eta value of a bin i is average
114  // of eta[i] and eta[i-1]
115  static const double z_long=1100.0;
116  double abseta = fabs(r.eta());
117  double absz = fabs(r.z());
118 
119  // figure out subdetector, giving preference to HE in HE/HF overlap region
121  if (abseta <= m_topology.etaMax(HcalBarrel) ) {
122  bc = HcalBarrel;
123  } else if (absz >= z_long) {
124  bc = HcalForward;
125  } else if (abseta <= m_topology.etaMax(HcalEndcap) ) {
126  bc = HcalEndcap;
127  } else {
128  bc = HcalForward;
129  }
130 
131  // find eta bin
132  int etaring = etaRing(bc, abseta);
133 
134  int phibin = phiBin(bc, etaring, r.phi());
135 
136  // add a sign to the etaring
137  int etabin = (r.z() > 0) ? etaring : -etaring;
138 
139  if (bc == HcalForward) {
140  static const double z_short=1137.0;
141  // Next line is premature depth 1 and 2 can coexist for large z-extent
142  // HcalDetId bestId(bc,etabin,phibin,((fabs(r.z())>=z_short)?(2):(1)));
143  // above line is no good with finite precision
144  HcalDetId bestId(bc,etabin,phibin,((fabs(r.z()) - z_short >-0.1)?(2):(1)));
145  return correctId(bestId);
146  } else {
147 
148  //Now do depth if required
149  int zside = (r.z() > 0) ? 1 : -1;
150  int dbin = m_topology.dddConstants()->getMinDepth(((int)(bc)-1),etaring,phibin,zside);
151  double pointrz(0), drz(99999.);
152  HcalDetId currentId(bc, etabin, phibin, dbin);
153  if (bc == HcalBarrel) pointrz = r.mag();
154  else pointrz = std::abs(r.z());
155  HcalDetId bestId;
156  for ( ; currentId != HcalDetId(); m_topology.incrementDepth(currentId)) {
157  const CaloCellGeometry * cell = getGeometry(currentId);
158  if (cell == nullptr) {
159  assert (bestId != HcalDetId());
160  break;
161  } else {
162  double rz;
163  if (bc == HcalEndcap) rz = std::abs(cell->getPosition().z());
164  else rz = cell->getPosition().mag();
165  if (std::abs(pointrz-rz)<drz) {
166  bestId = currentId;
167  drz = std::abs(pointrz-rz);
168  }
169  }
170  }
171 #ifdef EDM_ML_DEBUG
172  std::cout << bestId << " Corrected to " << HcalDetId(correctId(bestId))
173  << std::endl;
174 #endif
175 
176  return (ignoreCorrect ? bestId : correctId(bestId));
177  }
178 }
179 
181  if (!m_mergePosition) {
182  return (getGeometryBase(id)->getPosition());
183  } else {
185  }
186 }
187 
189  if (!m_mergePosition) {
190  return (getGeometryBase(id)->getBackPoint());
191  } else {
192  std::vector<HcalDetId> ids;
194  return (getGeometryBase(ids.back())->getBackPoint());
195  }
196 }
197 
199  if (!m_mergePosition) {
200  return (getGeometryBase(id)->getCorners());
201  } else {
202  std::vector<HcalDetId> ids;
207  for (unsigned int k=0; k<4; ++k) {
208  mcorners[k] = mcf[k];
209  mcorners[k+4] = mcb[k+4];
210  }
211  return mcorners;
212  }
213 }
214 
215 int HcalGeometry::etaRing(HcalSubdetector bc, double abseta) const {
216  return m_topology.etaRing(bc, abseta);
217 }
218 
219 int HcalGeometry::phiBin(HcalSubdetector bc, int etaring, double phi) const {
220  return m_topology.phiBin(bc, etaring, phi);
221 }
222 
224  double dR ) const {
225  CaloSubdetectorGeometry::DetIdSet dis; // this is the return object
226 
227  if (0.000001 < dR) {
228  if (dR > M_PI/2.) {// this version needs "small" dR
229  dis = CaloSubdetectorGeometry::getCells(r, dR); // base class version
230  } else {
231  const double dR2 ( dR*dR ) ;
232  const double reta ( r.eta() ) ;
233  const double rphi ( r.phi() ) ;
234  const double lowEta ( reta - dR ) ;
235  const double highEta ( reta + dR ) ;
236  const double lowPhi ( rphi - dR ) ;
237  const double highPhi ( rphi + dR ) ;
238 
239  const double hfEtaHi (m_topology.etaMax(HcalForward));
240 
241  if (highEta > -hfEtaHi &&
242  lowEta < hfEtaHi ) { // in hcal
244 
245  for (unsigned int is ( 0 ) ; is != 4 ; ++is ) {
246  const int sign ( reta>0 ? 1 : -1 ) ;
247  const int ieta_center ( sign*etaRing( hs[is], fabs( reta ) ) ) ;
248  const int ieta_lo ( ( 0 < lowEta*sign ? sign : -sign )*etaRing( hs[is], fabs( lowEta ) ) ) ;
249  const int ieta_hi ( ( 0 < highEta*sign ? sign : -sign )*etaRing( hs[is], fabs( highEta ) ) ) ;
250  const int iphi_lo ( phiBin( hs[is], ieta_center, lowPhi ) ) ;
251  const int iphi_hi ( phiBin( hs[is], ieta_center, highPhi ) ) ;
252  const int jphi_lo ( iphi_lo>iphi_hi ? iphi_lo - 72 : iphi_lo ) ;
253  const int jphi_hi ( iphi_hi ) ;
254 
255  const int idep_lo ( 1 == is ? 4 : 1 ) ;
256  const int idep_hi ( m_topology.maxDepth(hs[is]) );
257  for (int ieta ( ieta_lo ) ; ieta <= ieta_hi ; ++ieta) {// over eta limits
258  if (ieta != 0) {
259  for (int jphi ( jphi_lo ) ; jphi <= jphi_hi ; ++jphi) { // over phi limits
260  const int iphi ( 1 > jphi ? jphi+72 : jphi ) ;
261  for (int idep ( idep_lo ) ; idep <= idep_hi ; ++idep ) {
262  const HcalDetId did ( hs[is], ieta, iphi, idep ) ;
263  if (m_topology.valid(did)) {
264  const CaloCellGeometry* cell ( getGeometryBase( did ) );
265  if (nullptr != cell ) {
266  const GlobalPoint& p ( cell->getPosition() ) ;
267  const double eta ( p.eta() ) ;
268  const double phi ( p.phi() ) ;
269  if (reco::deltaR2(eta, phi, reta, rphi ) < dR2) dis.insert( did ) ;
270  }
271  }
272  }
273  }
274  }
275  }
276  }
277  }
278  }
279  }
280  return dis;
281 }
282 
283 unsigned int HcalGeometry::getHxSize(const int type) const {
284  unsigned int hxsize(0);
285  if (!m_hbIds.isSet()) fillDetIds();
286  if (type == 1) hxsize = (m_hbIds.isSet()) ? m_hbIds->size() : 0;
287  else if (type == 2) hxsize = (m_heIds.isSet()) ? m_heIds->size() : 0;
288  else if (type == 3) hxsize = (m_hoIds.isSet()) ? m_hoIds->size() : 0;
289  else if (type == 4) hxsize = (m_hfIds.isSet()) ? m_hfIds->size() : 0;
290  else if (type == 0) hxsize = (m_emptyIds.isSet()) ? m_emptyIds->size() : 0;
291  return hxsize;
292 }
293 
295  assert( i < numberOfBarrelAlignments() ) ;
296  const int ieta ( i < numberOfBarrelAlignments()/2 ? -1 : 1 ) ;
297  const int iphi ( 1 + (4*i)%72 ) ;
298  return HcalDetId( HcalBarrel, ieta, iphi, 1 ) ;
299 }
300 
302  assert( i < numberOfEndcapAlignments() ) ;
303  const int ieta ( i < numberOfEndcapAlignments()/2 ? -16 : 16 ) ;
304  const int iphi ( 1 + (4*i)%72 ) ;
305  return HcalDetId( HcalEndcap, ieta, iphi, 1 ) ;
306 }
307 
309  assert( i < numberOfForwardAlignments() ) ;
310  const int ieta ( i < numberOfForwardAlignments()/2 ? -29 : 29 ) ;
311  const int iphi ( 1 + (4*i)%72 ) ;
312  return HcalDetId( HcalForward, ieta, iphi, 1 ) ;
313 }
314 
316  assert( i < numberOfOuterAlignments() ) ;
317  const int ring ( i/12 ) ;
318  const int ieta ( 0 == ring ? -11 :
319  1 == ring ? -5 :
320  2 == ring ? 1 :
321  3 == ring ? 5 : 11 ) ;
322  const int iphi ( 1 + ( i - ring*12 )*6 ) ;
323  return HcalDetId( HcalOuter, ieta, iphi, 4 ) ;
324 }
325 
327  assert( i < numberOfAlignments() ) ;
328 
329  const unsigned int nB ( numberOfBarrelAlignments() ) ;
330  const unsigned int nE ( numberOfEndcapAlignments() ) ;
331  const unsigned int nF ( numberOfForwardAlignments() ) ;
332 // const unsigned int nO ( numberOfOuterAlignments() ) ;
333 
334  return ( i < nB ? detIdFromBarrelAlignmentIndex( i ) :
335  i < nB+nE ? detIdFromEndcapAlignmentIndex( i - nB ) :
336  i < nB+nE+nF ? detIdFromForwardAlignmentIndex( i - nB - nE ) :
337  detIdFromOuterAlignmentIndex( i - nB - nE - nF ) ) ;
338 }
339 
341  unsigned int nD) {
342  const HcalDetId hid ( id ) ;
343  const unsigned int iphi ( hid.iphi() ) ;
344  const int ieta ( hid.ieta() ) ;
345  const unsigned int index ( ( 0 < ieta ? nD/2 : 0 ) + ( iphi + 1 )%72/4 ) ;
346  assert( index < nD ) ;
347  return index ;
348 }
349 
352 }
353 
356 }
357 
360 }
361 
363  const HcalDetId hid ( id ) ;
364  const int ieta ( hid.ieta() ) ;
365  const int iphi ( hid.iphi() ) ;
366  const int ring ( ieta < -10 ? 0 :
367  ( ieta < -4 ? 1 :
368  ( ieta < 5 ? 2 :
369  ( ieta < 11 ? 3 : 4 ) ) ) ) ;
370 
371  const unsigned int index ( 12*ring + ( iphi - 1 )/6 ) ;
372  assert( index < numberOfOuterAlignments() ) ;
373  return index ;
374 }
375 
377  assert(id.det() == DetId::Hcal) ;
378 
379  const HcalDetId hid ( id ) ;
380  bool isHB = (hid.subdet() == HcalBarrel);
381  bool isHE = (hid.subdet() == HcalEndcap);
382  bool isHF = (hid.subdet() == HcalForward);
383  // bool isHO = (hid.subdet() == HcalOuter);
384 
385  const unsigned int nB ( numberOfBarrelAlignments() ) ;
386  const unsigned int nE ( numberOfEndcapAlignments() ) ;
387  const unsigned int nF ( numberOfForwardAlignments() ) ;
388  // const unsigned int nO ( numberOfOuterAlignments() ) ;
389 
390  const unsigned int index (isHB ? alignmentBarrelIndexLocal(id) :
391  isHE ? alignmentEndcapIndexLocal(id) + nB :
392  isHF ? alignmentForwardIndexLocal( id ) + nB + nE :
393  alignmentOuterIndexLocal(id) + nB + nE + nF
394  );
395 
396  assert( index < numberOfAlignments() ) ;
397  return index ;
398 }
399 
401  return (unsigned int)DetId::Hcal - 1 ;
402 }
403 
405  const CCGFloat* pv,
406  unsigned int i,
407  Pt3D& ref) {
409 
410  if (hid.subdet() == HcalForward ) {
411  IdealZPrism::localCorners( lc, pv, ref ) ;
412  } else {
413  IdealObliquePrism::localCorners( lc, pv, ref ) ;
414  }
415 }
416 
418  const GlobalPoint& f2 ,
419  const GlobalPoint& f3 ,
420  const CCGFloat* parm ,
421  const DetId& detId) {
422 
423  assert (detId.det()==DetId::Hcal);
424 
425  const HcalDetId hid ( detId ) ;
426  unsigned int din=m_topology.detId2denseId(detId);
427 
428 #ifdef EDM_ML_DEBUG
429  std::cout << "HcalGeometry: " << " newCell subdet "
430  << detId.subdetId() << ", raw ID "
431  << detId.rawId() << ", hid " << hid << ", din "
432  << din << ", index " << std::endl;
433 #endif
434 
435  if (hid.subdet()==HcalBarrel) {
436  m_hbCellVec.at( din ) = IdealObliquePrism( f1, cornersMgr(), parm ) ;
437  } else if (hid.subdet()==HcalEndcap) {
438  const unsigned int index ( din - m_hbCellVec.size() ) ;
439  m_heCellVec.at( index ) = IdealObliquePrism( f1, cornersMgr(), parm ) ;
440  } else if (hid.subdet()==HcalOuter) {
441  const unsigned int index ( din
442  - m_hbCellVec.size()
443  - m_heCellVec.size() ) ;
444  m_hoCellVec.at( index ) = IdealObliquePrism( f1, cornersMgr(), parm ) ;
445  } else {
446  const unsigned int index ( din
447  - m_hbCellVec.size()
448  - m_heCellVec.size()
449  - m_hoCellVec.size() ) ;
450  m_hfCellVec.at( index ) = IdealZPrism( f1, cornersMgr(), parm, hid.depth()==1 ? IdealZPrism::EM : IdealZPrism::HADR ) ;
451  }
452 
453  return din;
454 }
455 
457  const GlobalPoint& f2 ,
458  const GlobalPoint& f3 ,
459  const CCGFloat* parm ,
460  const DetId& detId) {
461 
462  unsigned int din = newCellImpl(f1,f2,f3,parm,detId);
463 
464  addValidID( detId ) ;
465  m_dins.emplace_back( din );
466 }
467 
469  const GlobalPoint& f2 ,
470  const GlobalPoint& f3 ,
471  const CCGFloat* parm ,
472  const DetId& detId) {
473 
474  unsigned int din = newCellImpl(f1,f2,f3,parm,detId);
475 
476  m_validIds.emplace_back( detId ) ;
477  m_dins.emplace_back( din );
478 }
479 
480 const CaloCellGeometry* HcalGeometry::cellGeomPtr( unsigned int din ) const {
481  const CaloCellGeometry* cell ( nullptr ) ;
482  if( m_hbCellVec.size() > din ) {
483  cell = &m_hbCellVec[ din ] ;
484  } else {
485  if (m_hbCellVec.size() + m_heCellVec.size() > din) {
486  const unsigned int index (din - m_hbCellVec.size() ) ;
487  cell = &m_heCellVec[ index ] ;
488  } else if (m_hbCellVec.size()+m_heCellVec.size()+m_hoCellVec.size() > din) {
489  const unsigned int index (din - m_hbCellVec.size() - m_heCellVec.size());
490  cell = &m_hoCellVec[ index ] ;
491  } else if (m_hbCellVec.size() + m_heCellVec.size() + m_hoCellVec.size() +
492  m_hfCellVec.size() > din) {
493  const unsigned int index (din - m_hbCellVec.size() - m_heCellVec.size() -
494  m_hoCellVec.size() ) ;
495  cell = &m_hfCellVec[ index ] ;
496  }
497  }
498 
499  return (( nullptr == cell || nullptr == cell->param()) ? nullptr : cell ) ;
500 }
501 
505  CaloSubdetectorGeometry::IVec& dinsVec ) const {
506 
507  tVec.reserve( m_topology.ncells()*numberOfTransformParms());
508  iVec.reserve( numberOfShapes()==1 ? 1 : m_topology.ncells());
509  dVec.reserve( numberOfShapes()*numberOfParametersPerShape());
510  dinsVec.reserve( m_topology.ncells());
511 
512  for (const auto & pv : parVecVec()) {
513  for (float iv : pv) {
514  dVec.emplace_back( iv ) ;
515  }
516  }
517 
518  for( auto i : m_dins ) {
519  Tr3D tr ;
520  const CaloCellGeometry* ptr ( cellGeomPtr( i ) ) ;
521 
522  if (nullptr != ptr) {
523  dinsVec.emplace_back( i );
524 
525  const CCGFloat* par ( ptr->param() ) ;
526 
527  unsigned int ishape ( 9999 ) ;
528 
529  for( unsigned int ivv ( 0 ) ; ivv != parVecVec().size() ; ++ivv ) {
530  bool ok ( true ) ;
531  const CCGFloat* pv ( &(*parVecVec()[ivv].begin() ) ) ;
532  for( unsigned int k ( 0 ) ; k != numberOfParametersPerShape() ; ++k ) {
533  ok = ok && ( fabs( par[k] - pv[k] ) < 1.e-6 ) ;
534  }
535  if( ok ) {
536  ishape = ivv;
537  break ;
538  }
539  }
540  assert( 9999 != ishape ) ;
541 
542  const unsigned int nn (( numberOfShapes()==1) ? (unsigned int)1 : m_dins.size() ) ;
543  if( iVec.size() < nn ) iVec.emplace_back( ishape ) ;
544 
545  ptr->getTransform( tr, ( Pt3DVec* ) nullptr ) ;
546 
547  if( Tr3D() == tr ) { // for preshower there is no rotation
548  const GlobalPoint& gp ( ptr->getPosition() ) ;
549  tr = HepGeom::Translate3D( gp.x(), gp.y(), gp.z() ) ;
550  }
551 
552  const CLHEP::Hep3Vector tt ( tr.getTranslation() ) ;
553  tVec.emplace_back( tt.x() ) ;
554  tVec.emplace_back( tt.y() ) ;
555  tVec.emplace_back( tt.z() ) ;
556  if (6 == numberOfTransformParms()) {
557  const CLHEP::HepRotation rr ( tr.getRotation() ) ;
558  const ROOT::Math::Transform3D rtr (rr.xx(), rr.xy(), rr.xz(), tt.x(),
559  rr.yx(), rr.yy(), rr.yz(), tt.y(),
560  rr.zx(), rr.zy(), rr.zz(), tt.z());
562  rtr.GetRotation( ea ) ;
563  tVec.emplace_back( ea.Phi() ) ;
564  tVec.emplace_back( ea.Theta() ) ;
565  tVec.emplace_back( ea.Psi() ) ;
566  }
567  }
568  }
569 }
570 
572 
573  if (m_mergePosition) {
574  HcalDetId hid(id);
575  return ((DetId)(m_topology.mergedDepthDetId(hid)));
576  } else {
577  return id;
578  }
579 }
580 
581 void HcalGeometry::increaseReserve(unsigned int extra) {
582  m_validIds.reserve(m_validIds.size()+extra);
583 }
584 
586  std::sort(m_validIds.begin(),m_validIds.end());
587 }
588 
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:135
static unsigned int alignmentTransformIndexLocal(const DetId &id)
const HcalDDDRecConstants * dddConstants() const
Definition: HcalTopology.h:161
std::vector< CCGFloat > DimVec
unsigned int getHOSize() const
Definition: HcalTopology.h:134
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:161
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
virtual const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
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:43
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
const HcalTopology & m_topology
Definition: HcalGeometry.h:157
std::vector< unsigned int > IVec
edm::AtomicPtrCache< std::vector< DetId > > m_emptyIds
Definition: HcalGeometry.h:164
const CaloCellGeometry * getGeometryBase(const DetId &id) const
Definition: HcalGeometry.h:111
bool isHE(int etabin, int depth)
std::vector< CCGFloat > TrVec
HepGeom::Transform3D Tr3D
std::vector< IdealObliquePrism > HECellVec
Definition: HcalGeometry.h:25
std::vector< Pt3D > Pt3DVec
bool isHB(int etabin, int depth)
CaloCellGeometry::CCGFloat CCGFloat
bool m_mergePosition
Definition: HcalGeometry.h:158
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:163
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
CaloCellGeometry::Tr3D Tr3D
Definition: HcalGeometry.cc:12
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
void newCellFast(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId)
CaloCellGeometry::Pt3DVec Pt3DVec
Definition: HcalGeometry.h:31
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:167
HOCellVec m_hoCellVec
Definition: HcalGeometry.h:169
static void localCorners(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
HFCellVec m_hfCellVec
Definition: HcalGeometry.h:170
virtual void getTransform(Tr3D &tr, Pt3DVec *lptr) const
--------— only needed by specific utility; overloaded when needed -—
unsigned int numberOfParametersPerShape() const override
Definition: HcalGeometry.h:44
std::vector< DetId > m_validIds
std::vector< IdealObliquePrism > HOCellVec
Definition: HcalGeometry.h:26
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:165
CaloCellGeometry::CCGFloat CCGFloat
Definition: HcalGeometry.cc:9
int phiBin(HcalSubdetector bc, int etaring, double phi) const
const CaloCellGeometry * cellGeomPtr(unsigned int index) const override
static DetId detIdFromLocalAlignmentIndex(unsigned int i)
int etaRing(HcalSubdetector bc, double abseta) const
helper methods for getClosestCell
CaloSubdetectorGeometry::IVec m_dins
Definition: HcalGeometry.h:165
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:37
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:27
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:124
~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:24
unsigned int getHESize() const
Definition: HcalTopology.h:133
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
HECellVec m_heCellVec
Definition: HcalGeometry.h:168
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:163
#define begin
Definition: vmac.h:30
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:170
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:162
CaloCellGeometry::Pt3D Pt3D
Definition: HcalGeometry.h:30
bool withSpecialRBXHBHE() const
Definition: HcalTopology.h:162
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
unsigned int getHBSize() const
Definition: HcalTopology.h:132
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:160
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