CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalBarrelGeometry.cc
Go to the documentation of this file.
7 
8 #include <CLHEP/Geometry/Point3D.h>
9 #include <CLHEP/Geometry/Plane3D.h>
10 #include <CLHEP/Geometry/Vector3D.h>
11 
12 #include <iomanip>
13 #include <iostream>
14 
18 typedef HepGeom::Plane3D<CCGFloat> Pl3D ;
19 
21  _nnxtalEta ( 85 ) ,
22  _nnxtalPhi ( 360 ) ,
23  _PhiBaskets ( 18 ) ,
24  m_borderMgr ( nullptr ),
25  m_borderPtrVec ( nullptr ) ,
26  m_radius ( -1. ),
27  m_check ( false ),
28  m_cellVec ( k_NumberOfCellsForCorners )
29 {
30  const int neba[] = {25,45,65,85} ;
31  _EtaBaskets = std::vector<int>( neba, neba+4 ) ;
32 }
33 
34 
36 {
37  if(m_borderPtrVec)
38  {
39  auto ptr = m_borderPtrVec.load(std::memory_order_acquire);
40  for(auto& v: (*ptr)) {
41  if(v) delete v;
42  v = nullptr;
43  }
44  delete m_borderPtrVec.load() ;
45  }
46  delete m_borderMgr.load() ;
47 }
48 
49 
50 unsigned int
52 {
53  const CaloGenericDetId gid ( id ) ;
54 
55  assert( gid.isEB() ) ;
56 
57  unsigned int index ( EBDetId(id).ism() - 1 ) ;
58 
59  return index ;
60 }
61 
62 DetId
64 {
65  return EBDetId( iLoc + 1, 1, EBDetId::SMCRYSTALMODE ) ;
66 }
67 
68 unsigned int
70 {
71  return (unsigned int)DetId::Ecal - 1 ;
72 }
73 // Get closest cell, etc...
74 DetId
76 {
77 
78  // z is the easy one
79  int leverx = 1;
80  int levery = 1;
81  CCGFloat pointz = r.z();
82  int zbin=1;
83  if(pointz<0)
84  zbin=-1;
85 
86  // Now find the closest eta
87  CCGFloat pointeta = r.eta();
88  // double eta;
89  CCGFloat deta=999.;
90  int etabin=1;
91 
92  int guessed_eta = (int)( fabs(pointeta) / 0.0174)+1;
93  int guessed_eta_begin = guessed_eta-1;
94  int guessed_eta_end = guessed_eta+1;
95  if (guessed_eta_begin < 1) guessed_eta_begin = 1;
96  if (guessed_eta_end > 85) guessed_eta_end = 85;
97 
98  for(int bin=guessed_eta_begin; bin<= guessed_eta_end; bin++)
99  {
100  try
101  {
102  if (!present(EBDetId(zbin*bin,1,EBDetId::ETAPHIMODE)))
103  continue;
104 
106 
107  if(fabs(pointeta-eta)<deta)
108  {
109  deta=fabs(pointeta-eta);
110  etabin=bin;
111  }
112  else break;
113  }
114  catch ( cms::Exception &e )
115  {
116  }
117  }
118 
119 
120  // Now the closest phi. always same number of phi bins(!?)
121  constexpr CCGFloat twopi = M_PI+M_PI;
122  // 10 degree tilt
123  constexpr CCGFloat tilt=twopi/36.;
124 
125  CCGFloat pointphi = r.phi()+tilt;
126 
127  // put phi in correct range (0->2pi)
128  if(pointphi > twopi)
129  pointphi -= twopi;
130  if(pointphi < 0)
131  pointphi += twopi;
132 
133  //calculate phi bin, distinguish + and - eta
134  int phibin = static_cast<int>(pointphi / (twopi/_nnxtalPhi)) + 1;
135  // if(point.z()<0.0)
136  // {
137  // phibin = nxtalPhi/2 - 1 - phibin;
138  // if(phibin<0)
139  // phibin += nxtalPhi;
140  // }
141  try
142  {
143  EBDetId myCell(zbin*etabin,phibin,EBDetId::ETAPHIMODE);
144 
145  if (!present(myCell))
146  return DetId(0);
147 
148  Pt3D A;
149  Pt3D B;
150  Pt3D C;
151  Pt3D point(r.x(),r.y(),r.z());
152 
153  // D.K. : equation of plane : AA*x+BB*y+CC*z+DD=0;
154  // finding equation for each edge
155 
156  // Since the point can lie between crystals, it is necessary to keep track of the movements
157  // to avoid infinite loops
158  CCGFloat history[4]{0.f};
159 
160  //
161  // stop movement in eta direction when closest cell was found (point between crystals)
162  int start = 1;
163  int counter = 0;
164  // Moving until find closest crystal in eta and phi directions (leverx and levery)
165  while (leverx==1 || levery == 1)
166  {
167  leverx = 0;
168  levery = 0;
169  const CaloCellGeometry::CornersVec& corners
170  ( getGeometry(myCell)->getCorners() ) ;
171  CCGFloat SS[4];
172 
173  // compute the distance of the point with respect of the 4 crystal lateral planes
174  for (short i=0; i < 4 ; ++i)
175  {
176  A = Pt3D(corners[i%4].x(),corners[i%4].y(),corners[i%4].z());
177  B = Pt3D(corners[(i+1)%4].x(),corners[(i+1)%4].y(),corners[(i+1)%4].z());
178  C = Pt3D(corners[4+(i+1)%4].x(),corners[4+(i+1)%4].y(),corners[4+(i+1)%4].z());
179  Pl3D plane(A,B,C);
180  plane.normalize();
181  CCGFloat distance = plane.distance(point);
182  if(plane.d()>0.) distance=-distance;
183  if (corners[0].z()<0.) distance=-distance;
184  SS[i] = distance;
185  }
186 
187  // SS's - normals
188  // check position of the point with respect to opposite side of crystal
189  // if SS's have opposite sign, the point lies inside that crystal
190 
191  if ( ( SS[0]>0.&&SS[2]>0. )||( SS[0]<0.&&SS[2]<0. ) )
192  {
193  levery = 1;
194  if ( history[0]>0. && history[2]>0. && SS[0]<0 && SS[2]<0 &&
195  (fabs(SS[0])+fabs(SS[2]))> (fabs(history[0])+fabs(history[2]))) levery = 0 ;
196  if ( history[0]<0. && history[2]<0. && SS[0]>0 && SS[2]>0 &&
197  (fabs(SS[0])+fabs(SS[2]))> (fabs(history[0])+fabs(history[2]))) levery = 0 ;
198 
199 
200  if (SS[0]>0. )
201  {
202  EBDetId nextPoint;
203  if (myCell.iphi()==EBDetId::MIN_IPHI)
204  nextPoint=EBDetId(myCell.ieta(),EBDetId::MAX_IPHI);
205  else
206  nextPoint=EBDetId(myCell.ieta(),myCell.iphi()-1);
207  if (present(nextPoint))
208  myCell=nextPoint;
209  else
210  levery=0;
211  }
212  else
213  {
214  EBDetId nextPoint;
215  if (myCell.iphi()==EBDetId::MAX_IPHI)
216  nextPoint=EBDetId(myCell.ieta(),EBDetId::MIN_IPHI);
217  else
218  nextPoint=EBDetId(myCell.ieta(),myCell.iphi()+1);
219  if (present(nextPoint))
220  myCell=nextPoint;
221  else
222  levery=0;
223  }
224  }
225 
226 
227  if ( ( ( SS[1]>0.&&SS[3]>0. )||( SS[1]<0.&&SS[3]<0. )) && start==1 )
228  {
229  leverx = 1;
230 
231  if ( history[1]>0. && history[3]>0. && SS[1]<0 && SS[3]<0 &&
232  (fabs(SS[1])+fabs(SS[3]))> (fabs(history[1])+fabs(history[3])) )
233  {
234  leverx = 0;
235  start = 0;
236  }
237 
238  if ( history[1]<0. && history[3]<0. && SS[1]>0 && SS[3]>0 &&
239  (fabs(SS[1])+fabs(SS[3]))> (fabs(history[1])+fabs(history[3])) )
240  {
241  leverx = 0;
242  start = 0;
243  }
244 
245 
246  if (SS[1]>0.)
247  {
248  EBDetId nextPoint;
249  if (myCell.ieta()==-1)
250  nextPoint=EBDetId (1,myCell.iphi());
251  else
252  {
253  int nieta= myCell.ieta()+1;
254  if(nieta==86) nieta=85;
255  nextPoint=EBDetId(nieta,myCell.iphi());
256  }
257  if (present(nextPoint))
258  myCell = nextPoint;
259  else
260  leverx = 0;
261  }
262  else
263  {
264  EBDetId nextPoint;
265  if (myCell.ieta()==1)
266  nextPoint=EBDetId(-1,myCell.iphi());
267  else
268  {
269  int nieta=myCell.ieta()-1;
270  if(nieta==-86) nieta=-85;
271  nextPoint=EBDetId(nieta,myCell.iphi());
272  }
273  if (present(nextPoint))
274  myCell = nextPoint;
275  else
276  leverx = 0;
277  }
278  }
279 
280  // Update the history. If the point lies between crystals, the closest one
281  // is returned
282  std::copy(SS,SS+4,history);
283 
284  counter++;
285  if (counter == 10)
286  {
287  leverx=0;
288  levery=0;
289  }
290  }
291  // D.K. if point lies netween cells, take a closest cell.
292  return DetId(myCell);
293  }
294  catch ( cms::Exception &e )
295  {
296  return DetId(0);
297  }
298 
299 }
300 
303  double dR ) const
304 {
305  constexpr int maxphi ( EBDetId::MAX_IPHI ) ;
306  constexpr int maxeta ( EBDetId::MAX_IETA ) ;
307  constexpr float scale( maxphi/(2*M_PI) ) ; // angle to index
308 
309  CaloSubdetectorGeometry::DetIdSet dis; // this is the return object
310 
311  if( 0.000001 < dR )
312  {
313  if( dR > M_PI/2. ) // this version needs "small" dR
314  {
315  dis = CaloSubdetectorGeometry::getCells( r, dR ) ; // base class version
316  }
317  else
318  {
319  const float dR2 ( dR*dR ) ;
320  const float reta ( r.eta() ) ;
321  const float rz ( r.z() ) ;
322  const float rphi ( r.phi() ) ;
323  const float lowEta ( reta - dR ) ;
324  const float highEta ( reta + dR ) ;
325 
326  if( highEta > -1.5 &&
327  lowEta < 1.5 ) // in barrel
328  {
329  const int ieta_center ( int( reta*scale + ((rz<0)?(-1):(1))) ) ;
330  const float phi ( rphi<0 ? rphi + float(2*M_PI) : rphi ) ;
331  const int iphi_center ( int( phi*scale + 11.f ) ) ; // phi=-9.4deg is iphi=1
332 
333  const float fr ( dR*scale ) ; // # crystal widths in dR
334  const float frp ( 1.08f*fr + 1.f ) ; // conservatively above fr
335  const float frm ( 0.92f*fr - 1.f ) ; // conservatively below fr
336  const int idr ( (int)frp ) ; // integerize
337  const int idr2p ( (int)(frp*frp) ) ;
338  const int idr2m ( frm > 0 ? int(frm*frm) : 0 ) ;
339 
340  for( int de ( -idr ) ; de <= idr ; ++de ) // over eta limits
341  {
342  int ieta ( de + ieta_center ) ;
343 
344  if( std::abs(ieta) <= maxeta &&
345  ieta != 0 ) // eta is in EB
346  {
347  const int de2 ( de*de ) ;
348  for( int dp ( -idr ) ; dp <= idr ; ++dp ) // over phi limits
349  {
350  const int irange2 ( dp*dp + de2 ) ;
351 
352  if( irange2 <= idr2p ) // cut off corners that must be too far away
353  {
354  const int iphi ( ( iphi_center + dp + maxphi - 1 )%maxphi + 1 ) ;
355 
356  if( iphi != 0 )
357  {
358  const EBDetId id ( ieta, iphi ) ;
359 
360  bool ok ( irange2 < idr2m ) ; // no more calculation necessary if inside this radius
361 
362  if( !ok ) // if not ok, then we have to test this cell for being inside cone
363  {
364  const CaloCellGeometry* cell = &m_cellVec[ id.denseIndex()];
365  const float eta ( cell->etaPos() ) ;
366  const float phi ( cell->phiPos() ) ;
367  ok = ( reco::deltaR2( eta, phi, reta, rphi ) < dR2 ) ;
368  }
369  if( ok ) dis.insert( id ) ;
370  }
371  }
372  }
373  }
374  }
375  }
376  }
377  }
378  return dis;
379 }
380 
383 {
384  OrderedListOfEEDetId* ptr ( nullptr ) ;
385  auto ptrVec = m_borderPtrVec.load(std::memory_order_acquire);
386  if(!ptrVec) {
387  if( 0 != id.rawId() ) {
388  const int iPhi ( id.iphi() ) ;
389  const int iz ( id.ieta()>0 ? 1 : -1 ) ;
390  const EEDetId eeid ( EEDetId::idOuterRing( iPhi, iz ) ) ;
391  const int iq ( eeid.iquadrant() ) ;
392  const int xout ( 1==iq || 4==iq ? 1 : -1 ) ;
393  const int yout ( 1==iq || 2==iq ? 1 : -1 ) ;
394  if (!m_borderMgr.load(std::memory_order_acquire)) {
395  EZMgrFL<EEDetId>* expect = nullptr;
396  auto ptrMgr = new EZMgrFL<EEDetId>( 720*9, 9 ) ;
397  bool exchanged = m_borderMgr.compare_exchange_strong(expect, ptrMgr, std::memory_order_acq_rel);
398  if(!exchanged) delete ptrMgr;
399  }
400  VecOrdListEEDetIdPtr* expect = nullptr;
401  auto ptrVec = new VecOrdListEEDetIdPtr();
402  ptrVec->reserve(720);
403  for( unsigned int i ( 0 ) ; i != 720 ; ++i )
404  {
405  const int kz ( 360>i ? -1 : 1 ) ;
406  const EEDetId eeid ( EEDetId::idOuterRing( i%360+1, kz ) ) ;
407 
408  const int jx ( eeid.ix() ) ;
409  const int jy ( eeid.iy() ) ;
410 
411  OrderedListOfEEDetId& olist ( *new OrderedListOfEEDetId( m_borderMgr.load(std::memory_order_acquire) ) );
412  int il ( 0 ) ;
413 
414  for( unsigned int k ( 1 ) ; k <= 25 ; ++k )
415  {
416  const int kx ( 1==k || 2==k || 3==k || 12==k || 13==k ? 0 :
417  ( 4==k || 6==k || 8==k || 15==k || 20==k ? 1 :
418  ( 5==k || 7==k || 9==k || 16==k || 19==k ? -1 :
419  ( 10==k || 14==k || 21==k || 22==k || 25==k ? 2 : -2 )))) ;
420  const int ky ( 1==k || 4==k || 5==k || 10==k || 11==k ? 0 :
421  ( 2==k || 6==k || 7==k || 14==k || 17==k ? 1 :
422  ( 3==k || 8==k || 9==k || 18==k || 21==k ? -1 :
423  ( 12==k || 15==k || 16==k || 22==k || 23==k ? 2 : -2 )))) ;
424 
425  if( 8>=il && EEDetId::validDetId( jx + kx*xout ,
426  jy + ky*yout , kz ) )
427  {
428  olist[il++]=EEDetId( jx + kx*xout, jy + ky*yout, kz ) ;
429  }
430  }
431  ptrVec->push_back( &olist ) ;
432  }
433  bool exchanged = m_borderPtrVec.compare_exchange_strong(expect, ptrVec, std::memory_order_acq_rel);
434  if(!exchanged) delete ptrVec;
435  ptrVec = m_borderPtrVec.load(std::memory_order_acquire);
436  ptr = (*ptrVec)[iPhi - 1 + ( 0>iz ? 0 : 360 )];
437  }
438  }
439  return ptr;
440 }
441 
442 void
444  const CCGFloat* pv ,
445  unsigned int i ,
446  Pt3D& ref )
447 {
448  const bool negz ( EBDetId::kSizeForDenseIndexing/2 > i ) ;
449  const bool odd ( 1 == i%2 ) ;
450 
451  if( ( ( negz && !odd ) ||
452  ( !negz && odd ) ) )
453  {
455  }
456  else
457  {
458  TruncatedPyramid::localCornersSwap( lc, pv, ref ) ;
459  }
460 }
461 
462 void
464  const GlobalPoint& f2 ,
465  const GlobalPoint& f3 ,
466  const CCGFloat* parm ,
467  const DetId& detId )
468 {
469  const unsigned int cellIndex ( EBDetId( detId ).denseIndex() ) ;
470  m_cellVec[ cellIndex ] =
471  TruncatedPyramid( cornersMgr(), f1, f2, f3, parm ) ;
472  addValidID( detId ) ;
473 }
474 
475 CCGFloat
477 {
478  if(!m_check.load(std::memory_order_acquire))
479  {
480  CCGFloat sum ( 0 ) ;
481  for( uint32_t i ( 0 ) ; i != m_cellVec.size() ; ++i )
482  {
483  const CaloCellGeometry* cell ( cellGeomPtr(i) ) ;
484  if( 0 != cell )
485  {
486  const GlobalPoint& pos ( cell->getPosition() ) ;
487  sum += pos.perp() ;
488  }
489  }
490  m_radius = sum/m_cellVec.size() ;
491  m_check.store(true, std::memory_order_release);
492  }
493  return m_radius ;
494 }
495 
496 const CaloCellGeometry*
498 {
499  const CaloCellGeometry* cell ( &m_cellVec[ index ] ) ;
500  return ( m_cellVec.size() < index ||
501  0 == cell->param() ? 0 : cell ) ;
502 }
static unsigned int alignmentTransformIndexGlobal(const DetId &id)
int i
Definition: DBlmapReader.cc:9
const OrderedListOfEEDetId * getClosestEndcapCells(EBDetId id) const
static const int MIN_IPHI
Definition: EBDetId.h:142
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
virtual CaloSubdetectorGeometry::DetIdSet getCells(const GlobalPoint &r, double dR) const
Get a list of all cells within a dR of the given cell.
int ix() const
Definition: EEDetId.h:76
T perp() const
Definition: PV3DBase.h:72
virtual DetId getClosestCell(const GlobalPoint &r) const
CaloCellGeometry::Pt3DVec Pt3DVec
double_binary B
Definition: DDStreamer.cc:234
float phiPos() const
std::atomic< bool > m_check
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
std::vector< int > _EtaBaskets
int iquadrant() const
Definition: EEDetId.cc:264
std::vector< Pt3D > Pt3DVec
T eta() const
#define nullptr
EZArrayFL< EEDetId > OrderedListOfEEDetId
CaloCellGeometry::CCGFloat CCGFloat
#define constexpr
float float float z
virtual const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
int iphi() const
get the crystal iphi
Definition: EBDetId.h:53
HepGeom::Plane3D< CCGFloat > Pl3D
virtual bool present(const DetId &id) const
is this detid present in the geometry?
static void localCornersSwap(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
bool isEB() 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
static DetId detIdFromLocalAlignmentIndex(unsigned int iLoc)
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CaloCellGeometry::CCGFloat CCGFloat
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
double f[11][100]
int iy() const
Definition: EEDetId.h:82
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
std::atomic< VecOrdListEEDetIdPtr * > m_borderPtrVec
static const int ETAPHIMODE
Definition: EBDetId.h:166
#define M_PI
CCGFloat avgRadiusXYFrontFaceCenter() const
CaloCellGeometry::Pt3D Pt3D
Definition: DetId.h:18
void addValidID(const DetId &id)
static const int MAX_IPHI
Definition: EBDetId.h:144
virtual const CaloCellGeometry * cellGeomPtr(uint32_t index) const
static bool validDetId(int crystal_ix, int crystal_iy, int iz)
Definition: EEDetId.h:248
CaloCellGeometry::CornersMgr * cornersMgr()
static const int MAX_IETA
Definition: EBDetId.h:143
CaloCellGeometry::Pt3DVec Pt3DVec
HepGeom::Point3D< CCGFloat > Pt3D
Definition: EZMgrFL.h:8
std::atomic< EZMgrFL< EEDetId > * > m_borderMgr
virtual void newCell(const GlobalPoint &f1, const GlobalPoint &f2, const GlobalPoint &f3, const CCGFloat *parm, const DetId &detId)
A base class to handle the particular shape of Ecal Xtals. Taken from ORCA Calorimetry Code...
T eta() const
Definition: PV3DBase.h:76
std::vector< OrderedListOfEEDetId * > VecOrdListEEDetIdPtr
static std::atomic< unsigned int > counter
float etaPos() const
static EEDetId idOuterRing(int iPhi, int zEnd)
Definition: EEDetId.cc:435
volatile std::atomic< bool > shutdown_flag false
Definition: DDAxes.h:10
static const int SMCRYSTALMODE
Definition: EBDetId.h:167
const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
int ism(int ieta, int iphi)
Definition: EcalPyUtils.cc:56
T x() const
Definition: PV3DBase.h:62
static void localCornersReflection(Pt3DVec &vec, const CCGFloat *pv, Pt3D &ref)
static unsigned int alignmentTransformIndexLocal(const DetId &id)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
CaloCellGeometry::Pt3D Pt3D
static void localCorners(Pt3DVec &lc, const CCGFloat *pv, unsigned int i, Pt3D &ref)
Definition: DDAxes.h:10