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