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