00001 #include "Geometry/CaloGeometry/interface/CaloGenericDetId.h"
00002 #include "Geometry/CaloGeometry/interface/TruncatedPyramid.h"
00003 #include "Geometry/EcalAlgo/interface/EcalBarrelGeometry.h"
00004 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00005 #include "FWCore/Utilities/interface/Exception.h"
00006 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00007
00008 #include <CLHEP/Geometry/Point3D.h>
00009 #include <CLHEP/Geometry/Plane3D.h>
00010 #include <CLHEP/Geometry/Vector3D.h>
00011
00012 #include <iomanip>
00013 #include <iostream>
00014
00015 EcalBarrelGeometry::EcalBarrelGeometry() :
00016 _nnxtalEta ( 85 ) ,
00017 _nnxtalPhi ( 360 ) ,
00018 _PhiBaskets ( 18 ) ,
00019 m_borderMgr ( 0 ),
00020 m_borderPtrVec ( 0 ) ,
00021 m_radius ( -1. )
00022 {
00023 const int neba[] = {25,45,65,85} ;
00024 _EtaBaskets = std::vector<int>( &neba[0], &neba[3] ) ;
00025 }
00026
00027
00028 EcalBarrelGeometry::~EcalBarrelGeometry()
00029 {
00030 delete m_borderPtrVec ;
00031 delete m_borderMgr ;
00032 }
00033
00034
00035 unsigned int
00036 EcalBarrelGeometry::alignmentTransformIndexLocal( const DetId& id )
00037 {
00038 const CaloGenericDetId gid ( id ) ;
00039
00040 assert( gid.isEB() ) ;
00041
00042 unsigned int index ( EBDetId(id).ism() - 1 ) ;
00043
00044 return index ;
00045 }
00046
00047 DetId
00048 EcalBarrelGeometry::detIdFromLocalAlignmentIndex( unsigned int iLoc )
00049 {
00050 return EBDetId( iLoc + 1, 1, EBDetId::SMCRYSTALMODE ) ;
00051 }
00052
00053 unsigned int
00054 EcalBarrelGeometry::alignmentTransformIndexGlobal( const DetId& id )
00055 {
00056 return (unsigned int)DetId::Ecal - 1 ;
00057 }
00058
00059 DetId
00060 EcalBarrelGeometry::getClosestCell(const GlobalPoint& r) const
00061 {
00062
00063
00064 int leverx = 1;
00065 int levery = 1;
00066 double pointz = r.z();
00067 int zbin=1;
00068 if(pointz<0)
00069 zbin=-1;
00070
00071
00072 double pointeta = r.eta();
00073
00074 double deta=999.;
00075 int etabin=1;
00076
00077 int guessed_eta = (int)( fabs(pointeta) / 0.0174)+1;
00078 int guessed_eta_begin = guessed_eta-1;
00079 int guessed_eta_end = guessed_eta+1;
00080 if (guessed_eta_begin < 1) guessed_eta_begin = 1;
00081 if (guessed_eta_end > 85) guessed_eta_end = 85;
00082
00083 for(int bin=guessed_eta_begin; bin<= guessed_eta_end; bin++)
00084 {
00085 try
00086 {
00087 if (!present(EBDetId(zbin*bin,1,EBDetId::ETAPHIMODE)))
00088 continue;
00089
00090 double eta = getGeometry(EBDetId(zbin*bin,1,EBDetId::ETAPHIMODE))->getPosition().eta();
00091
00092 if(fabs(pointeta-eta)<deta)
00093 {
00094 deta=fabs(pointeta-eta);
00095 etabin=bin;
00096 }
00097 else break;
00098 }
00099 catch ( cms::Exception &e )
00100 {
00101 }
00102 }
00103
00104
00105
00106 const double twopi = M_PI+M_PI;
00107
00108
00109 const double tilt=twopi/36.;
00110 double pointphi = r.phi()+tilt;
00111
00112
00113 if(pointphi > twopi)
00114 pointphi -= twopi;
00115 if(pointphi < 0)
00116 pointphi += twopi;
00117
00118
00119 int phibin = static_cast<int>(pointphi / (twopi/_nnxtalPhi)) + 1;
00120
00121
00122
00123
00124
00125
00126 try
00127 {
00128 EBDetId myCell(zbin*etabin,phibin,EBDetId::ETAPHIMODE);
00129
00130 if (!present(myCell))
00131 return DetId(0);
00132
00133 HepGeom::Point3D<double> A;
00134 HepGeom::Point3D<double> B;
00135 HepGeom::Point3D<double> C;
00136 HepGeom::Point3D<double> point(r.x(),r.y(),r.z());
00137
00138
00139
00140
00141
00142
00143 std::vector<double> history;
00144 history.resize(4,0.);
00145
00146
00147 int start = 1;
00148 int counter = 0;
00149
00150 while (leverx==1 || levery == 1)
00151 {
00152 leverx = 0;
00153 levery = 0;
00154 const CaloCellGeometry::CornersVec& corners
00155 ( getGeometry(myCell)->getCorners() ) ;
00156 std::vector<double> SS;
00157
00158
00159 for (short i=0; i < 4 ; ++i)
00160 {
00161 A = HepGeom::Point3D<double> (corners[i%4].x(),corners[i%4].y(),corners[i%4].z());
00162 B = HepGeom::Point3D<double> (corners[(i+1)%4].x(),corners[(i+1)%4].y(),corners[(i+1)%4].z());
00163 C = HepGeom::Point3D<double> (corners[4+(i+1)%4].x(),corners[4+(i+1)%4].y(),corners[4+(i+1)%4].z());
00164 HepGeom::Plane3D<double> plane(A,B,C);
00165 plane.normalize();
00166 double distance = plane.distance(point);
00167 if(plane.d()>0.) distance=-distance;
00168 if (corners[0].z()<0.) distance=-distance;
00169 SS.push_back(distance);
00170 }
00171
00172
00173
00174
00175
00176 if ( ( SS[0]>0.&&SS[2]>0. )||( SS[0]<0.&&SS[2]<0. ) )
00177 {
00178 levery = 1;
00179 if ( history[0]>0. && history[2]>0. && SS[0]<0 && SS[2]<0 &&
00180 (fabs(SS[0])+fabs(SS[2]))> (fabs(history[0])+fabs(history[2]))) levery = 0 ;
00181 if ( history[0]<0. && history[2]<0. && SS[0]>0 && SS[2]>0 &&
00182 (fabs(SS[0])+fabs(SS[2]))> (fabs(history[0])+fabs(history[2]))) levery = 0 ;
00183
00184
00185 if (SS[0]>0. )
00186 {
00187 EBDetId nextPoint;
00188 if (myCell.iphi()==EBDetId::MIN_IPHI)
00189 nextPoint=EBDetId(myCell.ieta(),EBDetId::MAX_IPHI);
00190 else
00191 nextPoint=EBDetId(myCell.ieta(),myCell.iphi()-1);
00192 if (present(nextPoint))
00193 myCell=nextPoint;
00194 else
00195 levery=0;
00196 }
00197 else
00198 {
00199 EBDetId nextPoint;
00200 if (myCell.iphi()==EBDetId::MAX_IPHI)
00201 nextPoint=EBDetId(myCell.ieta(),EBDetId::MIN_IPHI);
00202 else
00203 nextPoint=EBDetId(myCell.ieta(),myCell.iphi()+1);
00204 if (present(nextPoint))
00205 myCell=nextPoint;
00206 else
00207 levery=0;
00208 }
00209 }
00210
00211
00212 if ( ( ( SS[1]>0.&&SS[3]>0. )||( SS[1]<0.&&SS[3]<0. )) && start==1 )
00213 {
00214 leverx = 1;
00215
00216 if ( history[1]>0. && history[3]>0. && SS[1]<0 && SS[3]<0 &&
00217 (fabs(SS[1])+fabs(SS[3]))> (fabs(history[1])+fabs(history[3])) )
00218 {
00219 leverx = 0;
00220 start = 0;
00221 }
00222
00223 if ( history[1]<0. && history[3]<0. && SS[1]>0 && SS[3]>0 &&
00224 (fabs(SS[1])+fabs(SS[3]))> (fabs(history[1])+fabs(history[3])) )
00225 {
00226 leverx = 0;
00227 start = 0;
00228 }
00229
00230
00231 if (SS[1]>0.)
00232 {
00233 EBDetId nextPoint;
00234 if (myCell.ieta()==-1)
00235 nextPoint=EBDetId (1,myCell.iphi());
00236 else
00237 {
00238 int nieta= myCell.ieta()+1;
00239 if(nieta==86) nieta=85;
00240 nextPoint=EBDetId(nieta,myCell.iphi());
00241 }
00242 if (present(nextPoint))
00243 myCell = nextPoint;
00244 else
00245 leverx = 0;
00246 }
00247 else
00248 {
00249 EBDetId nextPoint;
00250 if (myCell.ieta()==1)
00251 nextPoint=EBDetId(-1,myCell.iphi());
00252 else
00253 {
00254 int nieta=myCell.ieta()-1;
00255 if(nieta==-86) nieta=-85;
00256 nextPoint=EBDetId(nieta,myCell.iphi());
00257 }
00258 if (present(nextPoint))
00259 myCell = nextPoint;
00260 else
00261 leverx = 0;
00262 }
00263 }
00264
00265
00266
00267 history =SS;
00268
00269 counter++;
00270 if (counter == 10)
00271 {
00272 leverx=0;
00273 levery=0;
00274 }
00275 }
00276
00277 return DetId(myCell);
00278 }
00279 catch ( cms::Exception &e )
00280 {
00281 return DetId(0);
00282 }
00283
00284 }
00285
00286 CaloSubdetectorGeometry::DetIdSet
00287 EcalBarrelGeometry::getCells( const GlobalPoint& r,
00288 double dR ) const
00289 {
00290 static const int maxphi ( EBDetId::MAX_IPHI ) ;
00291 static const int maxeta ( EBDetId::MAX_IETA ) ;
00292 CaloSubdetectorGeometry::DetIdSet dis;
00293
00294 if( 0.000001 < dR )
00295 {
00296 if( dR > M_PI/2. )
00297 {
00298 dis = CaloSubdetectorGeometry::getCells( r, dR ) ;
00299 }
00300 else
00301 {
00302 const double dR2 ( dR*dR ) ;
00303 const double reta ( r.eta() ) ;
00304 const double rz ( r.z() ) ;
00305 const double rphi ( r.phi() ) ;
00306 const double lowEta ( reta - dR ) ;
00307 const double highEta ( reta + dR ) ;
00308
00309 if( highEta > -1.5 &&
00310 lowEta < 1.5 )
00311 {
00312 const double scale ( maxphi/(2*M_PI) ) ;
00313 const int ieta_center ( int( reta*scale + ((rz<0)?(-1):(1))) ) ;
00314 const double phi ( rphi<0 ? rphi + 2*M_PI : rphi ) ;
00315 const int iphi_center ( int( phi*scale + 11. ) ) ;
00316
00317 const double fr ( dR*scale ) ;
00318 const double frp ( 1.08*fr + 1. ) ;
00319 const double frm ( 0.92*fr - 1. ) ;
00320 const int idr ( (int)frp ) ;
00321 const int idr2p ( (int)(frp*frp) ) ;
00322 const int idr2m ( frm > 0 ? int(frm*frm) : 0 ) ;
00323
00324 for( int de ( -idr ) ; de <= idr ; ++de )
00325 {
00326 int ieta ( de + ieta_center ) ;
00327
00328 if( std::abs(ieta) <= maxeta &&
00329 ieta != 0 )
00330 {
00331 const int de2 ( de*de ) ;
00332 for( int dp ( -idr ) ; dp <= idr ; ++dp )
00333 {
00334 const int irange2 ( dp*dp + de2 ) ;
00335
00336 if( irange2 <= idr2p )
00337 {
00338 const int iphi ( ( iphi_center + dp + maxphi - 1 )%maxphi + 1 ) ;
00339
00340 if( iphi != 0 )
00341 {
00342 const EBDetId id ( ieta, iphi ) ;
00343
00344 bool ok ( irange2 < idr2m ) ;
00345
00346 if( !ok )
00347 {
00348 const CaloCellGeometry* cell ( getGeometry( id ) );
00349 if( 0 != cell )
00350 {
00351 const GlobalPoint& p ( cell->getPosition() ) ;
00352 const double eta ( p.eta() ) ;
00353 const double phi ( p.phi() ) ;
00354 ok = ( reco::deltaR2( eta, phi, reta, rphi ) < dR2 ) ;
00355 }
00356 }
00357 if( ok ) dis.insert( id ) ;
00358 }
00359 }
00360 }
00361 }
00362 }
00363 }
00364 }
00365 }
00366 return dis;
00367 }
00368
00369 const EcalBarrelGeometry::OrderedListOfEEDetId*
00370 EcalBarrelGeometry::getClosestEndcapCells( EBDetId id ) const
00371 {
00372 OrderedListOfEEDetId* ptr ( 0 ) ;
00373 if( 0 != id.rawId() )
00374 {
00375 const int iPhi ( id.iphi() ) ;
00376
00377 const int iz ( id.ieta()>0 ? 1 : -1 ) ;
00378 const EEDetId eeid ( EEDetId::idOuterRing( iPhi, iz ) ) ;
00379
00380
00381
00382
00383 const int iq ( eeid.iquadrant() ) ;
00384 const int xout ( 1==iq || 4==iq ? 1 : -1 ) ;
00385 const int yout ( 1==iq || 2==iq ? 1 : -1 ) ;
00386 if( 0 == m_borderMgr )
00387 {
00388 m_borderMgr = new EZMgrFL<EEDetId>( 720*9, 9 ) ;
00389 }
00390 if( 0 == m_borderPtrVec )
00391 {
00392 m_borderPtrVec = new VecOrdListEEDetIdPtr() ;
00393 m_borderPtrVec->reserve( 720 ) ;
00394 for( unsigned int i ( 0 ) ; i != 720 ; ++i )
00395 {
00396 const int kz ( 360>i ? -1 : 1 ) ;
00397 const EEDetId eeid ( EEDetId::idOuterRing( i%360+1, kz ) ) ;
00398
00399 const int jx ( eeid.ix() ) ;
00400 const int jy ( eeid.iy() ) ;
00401
00402 OrderedListOfEEDetId& olist ( *new OrderedListOfEEDetId( m_borderMgr ) );
00403 int il ( 0 ) ;
00404
00405 for( unsigned int k ( 1 ) ; k <= 25 ; ++k )
00406 {
00407 const int kx ( 1==k || 2==k || 3==k || 12==k || 13==k ? 0 :
00408 ( 4==k || 6==k || 8==k || 15==k || 20==k ? 1 :
00409 ( 5==k || 7==k || 9==k || 16==k || 19==k ? -1 :
00410 ( 10==k || 14==k || 21==k || 22==k || 25==k ? 2 : -2 )))) ;
00411 const int ky ( 1==k || 4==k || 5==k || 10==k || 11==k ? 0 :
00412 ( 2==k || 6==k || 7==k || 14==k || 17==k ? 1 :
00413 ( 3==k || 8==k || 9==k || 18==k || 21==k ? -1 :
00414 ( 12==k || 15==k || 16==k || 22==k || 23==k ? 2 : -2 )))) ;
00415
00416 if( 8>=il && EEDetId::validDetId( jx + kx*xout ,
00417 jy + ky*yout , kz ) )
00418 {
00419 olist[il++]=EEDetId( jx + kx*xout ,
00420 jy + ky*yout , kz ) ;
00421 }
00422 }
00423 m_borderPtrVec->push_back( &olist ) ;
00424 }
00425 }
00426 ptr = (*m_borderPtrVec)[ iPhi - 1 + ( 0>iz ? 0 : 360 ) ] ;
00427 }
00428 return ptr ;
00429 }
00430
00431 std::vector<HepGeom::Point3D<double> >
00432 EcalBarrelGeometry::localCorners( const double* pv,
00433 unsigned int i,
00434 HepGeom::Point3D<double> & ref )
00435 {
00436 const bool negz ( EBDetId::kSizeForDenseIndexing/2 > i ) ;
00437 const bool odd ( 1 == i%2 ) ;
00438
00439 return ( ( ( negz && !odd ) ||
00440 ( !negz && odd ) ) ? TruncatedPyramid::localCornersReflection( pv, ref ) :
00441 TruncatedPyramid::localCornersSwap( pv, ref ) ) ;
00442 }
00443
00444 CaloCellGeometry*
00445 EcalBarrelGeometry::newCell( const GlobalPoint& f1 ,
00446 const GlobalPoint& f2 ,
00447 const GlobalPoint& f3 ,
00448 CaloCellGeometry::CornersMgr* mgr,
00449 const double* parm ,
00450 const DetId& detId )
00451 {
00452 return ( new TruncatedPyramid( mgr, f1, f2, f3, parm ) ) ;
00453 }
00454
00455 double
00456 EcalBarrelGeometry::avgRadiusXYFrontFaceCenter() const
00457 {
00458 if( 0 > m_radius )
00459 {
00460 double sum ( 0 ) ;
00461 const CaloSubdetectorGeometry::CellCont& cells ( cellGeometries() ) ;
00462 for( unsigned int i ( 0 ) ; i != cells.size() ; ++i )
00463 {
00464 sum += cells[i]->getPosition().perp() ;
00465 }
00466 m_radius = sum/cells.size() ;
00467 }
00468 return m_radius ;
00469 }
00470