CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/Geometry/EcalAlgo/src/EcalBarrelGeometry.cc

Go to the documentation of this file.
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 typedef CaloCellGeometry::CCGFloat CCGFloat ;
00016 typedef CaloCellGeometry::Pt3D     Pt3D     ;
00017 typedef CaloCellGeometry::Pt3DVec  Pt3DVec  ;
00018 typedef HepGeom::Plane3D<CCGFloat> Pl3D     ;
00019 
00020 EcalBarrelGeometry::EcalBarrelGeometry() :
00021    _nnxtalEta     ( 85 ) ,
00022    _nnxtalPhi     ( 360 ) ,
00023    _PhiBaskets    ( 18 ) ,
00024    m_borderMgr    ( 0 ),
00025    m_borderPtrVec ( 0 ) ,
00026    m_radius       ( -1. ),
00027    m_cellVec      ( k_NumberOfCellsForCorners )
00028 {
00029    const int neba[] = {25,45,65,85} ;
00030    _EtaBaskets = std::vector<int>( &neba[0], &neba[3] ) ;
00031 }
00032 
00033 
00034 EcalBarrelGeometry::~EcalBarrelGeometry() 
00035 {
00036    delete m_borderPtrVec ;
00037    delete m_borderMgr ;
00038 }
00039 
00040 
00041 unsigned int
00042 EcalBarrelGeometry::alignmentTransformIndexLocal( const DetId& id )
00043 {
00044    const CaloGenericDetId gid ( id ) ;
00045 
00046    assert( gid.isEB() ) ;
00047 
00048    unsigned int index ( EBDetId(id).ism() - 1 ) ;
00049 
00050    return index ;
00051 }
00052 
00053 DetId 
00054 EcalBarrelGeometry::detIdFromLocalAlignmentIndex( unsigned int iLoc )
00055 {
00056    return EBDetId( iLoc + 1, 1, EBDetId::SMCRYSTALMODE ) ;
00057 }
00058 
00059 unsigned int
00060 EcalBarrelGeometry::alignmentTransformIndexGlobal( const DetId& /*id*/ )
00061 {
00062    return (unsigned int)DetId::Ecal - 1 ;
00063 }
00064 // Get closest cell, etc...
00065 DetId 
00066 EcalBarrelGeometry::getClosestCell(const GlobalPoint& r) const 
00067 {
00068 
00069   // z is the easy one
00070   int leverx = 1;
00071   int levery = 1;
00072   CCGFloat pointz = r.z();
00073   int zbin=1;
00074   if(pointz<0)
00075     zbin=-1;
00076 
00077   // Now find the closest eta
00078   CCGFloat pointeta = r.eta();
00079   //  double eta;
00080   CCGFloat deta=999.;
00081   int etabin=1;
00082   
00083   int guessed_eta = (int)( fabs(pointeta) / 0.0174)+1;
00084   int guessed_eta_begin = guessed_eta-1;
00085   int guessed_eta_end   = guessed_eta+1;
00086   if (guessed_eta_begin < 1) guessed_eta_begin = 1;
00087   if (guessed_eta_end > 85) guessed_eta_end = 85;
00088 
00089     for(int bin=guessed_eta_begin; bin<= guessed_eta_end; bin++)
00090       {
00091         try
00092           {
00093             if (!present(EBDetId(zbin*bin,1,EBDetId::ETAPHIMODE)))
00094               continue;
00095 
00096             CCGFloat eta = getGeometry(EBDetId(zbin*bin,1,EBDetId::ETAPHIMODE))->getPosition().eta();
00097 
00098             if(fabs(pointeta-eta)<deta)
00099               {
00100                 deta=fabs(pointeta-eta);
00101                 etabin=bin;
00102               }
00103             else break;
00104           }
00105         catch ( cms::Exception &e ) 
00106           {
00107           }
00108       }
00109     
00110 
00111   // Now the closest phi. always same number of phi bins(!?)
00112   const CCGFloat twopi = M_PI+M_PI;
00113 
00114   // 10 degree tilt
00115   const CCGFloat tilt=twopi/36.;
00116   CCGFloat pointphi = r.phi()+tilt;
00117 
00118   // put phi in correct range (0->2pi)
00119   if(pointphi > twopi)
00120     pointphi -= twopi;
00121   if(pointphi < 0)
00122     pointphi += twopi;
00123 
00124   //calculate phi bin, distinguish + and - eta
00125   int phibin = static_cast<int>(pointphi / (twopi/_nnxtalPhi)) + 1;
00126   //   if(point.z()<0.0)
00127   //     {
00128   //       phibin = nxtalPhi/2 - 1 - phibin;
00129   //       if(phibin<0)
00130   //         phibin += nxtalPhi;
00131   //     }
00132   try
00133     {
00134       EBDetId myCell(zbin*etabin,phibin,EBDetId::ETAPHIMODE);
00135 
00136       if (!present(myCell))
00137         return DetId(0);
00138       
00139       Pt3D A;
00140       Pt3D B;
00141       Pt3D C;
00142       Pt3D point(r.x(),r.y(),r.z());
00143 
00144       // D.K. : equation of plane : AA*x+BB*y+CC*z+DD=0;
00145       // finding equation for each edge
00146 
00147       // Since the point can lie between crystals, it is necessary to keep track of the movements
00148       // to avoid infinite loops
00149       std::vector<CCGFloat> history;
00150       history.resize(4,0.);
00151       //
00152       // stop movement in eta direction when closest cell was found (point between crystals)
00153       int start = 1;
00154       int counter = 0;
00155       // Moving until find closest crystal in eta and phi directions (leverx and levery)
00156       while  (leverx==1 || levery == 1)
00157         {
00158           leverx = 0;
00159           levery = 0;
00160           const CaloCellGeometry::CornersVec& corners 
00161              ( getGeometry(myCell)->getCorners() ) ;
00162           std::vector<CCGFloat> SS;
00163 
00164           // compute the distance of the point with respect of the 4 crystal lateral planes
00165           for (short i=0; i < 4 ; ++i)
00166             {
00167               A = Pt3D(corners[i%4].x(),corners[i%4].y(),corners[i%4].z());
00168               B = Pt3D(corners[(i+1)%4].x(),corners[(i+1)%4].y(),corners[(i+1)%4].z());
00169               C = Pt3D(corners[4+(i+1)%4].x(),corners[4+(i+1)%4].y(),corners[4+(i+1)%4].z());
00170               Pl3D plane(A,B,C);
00171               plane.normalize();
00172               CCGFloat distance = plane.distance(point);
00173               if(plane.d()>0.) distance=-distance;
00174               if (corners[0].z()<0.) distance=-distance;
00175               SS.push_back(distance);
00176             }
00177 
00178           // SS's - normals
00179           // check position of the point with respect to opposite side of crystal
00180           // if SS's have opposite sign, the  point lies inside that crystal
00181 
00182           if ( ( SS[0]>0.&&SS[2]>0. )||( SS[0]<0.&&SS[2]<0. ) )
00183             {
00184               levery = 1;
00185               if ( history[0]>0. && history[2]>0. && SS[0]<0 && SS[2]<0 &&
00186                    (fabs(SS[0])+fabs(SS[2]))> (fabs(history[0])+fabs(history[2]))) levery = 0  ;
00187               if ( history[0]<0. && history[2]<0. && SS[0]>0 && SS[2]>0 &&
00188                    (fabs(SS[0])+fabs(SS[2]))> (fabs(history[0])+fabs(history[2]))) levery = 0  ;
00189 
00190 
00191               if (SS[0]>0. )
00192                 {
00193                   EBDetId nextPoint;
00194                   if (myCell.iphi()==EBDetId::MIN_IPHI) 
00195                     nextPoint=EBDetId(myCell.ieta(),EBDetId::MAX_IPHI);
00196                   else 
00197                     nextPoint=EBDetId(myCell.ieta(),myCell.iphi()-1);
00198                   if (present(nextPoint))
00199                     myCell=nextPoint;
00200                   else
00201                     levery=0;             
00202                 }
00203               else
00204                 {
00205                   EBDetId nextPoint;
00206                   if (myCell.iphi()==EBDetId::MAX_IPHI)
00207                     nextPoint=EBDetId(myCell.ieta(),EBDetId::MIN_IPHI);
00208                   else
00209                     nextPoint=EBDetId(myCell.ieta(),myCell.iphi()+1);
00210                   if (present(nextPoint))
00211                     myCell=nextPoint;
00212                   else
00213                     levery=0;
00214                 }
00215             }
00216 
00217 
00218           if ( ( ( SS[1]>0.&&SS[3]>0. )||( SS[1]<0.&&SS[3]<0. )) && start==1  )
00219             {
00220               leverx = 1;
00221 
00222               if ( history[1]>0. && history[3]>0. && SS[1]<0 && SS[3]<0 &&
00223                    (fabs(SS[1])+fabs(SS[3]))> (fabs(history[1])+fabs(history[3])) )
00224                 {
00225                   leverx = 0;
00226                   start = 0;
00227                 }
00228 
00229               if ( history[1]<0. && history[3]<0. && SS[1]>0 && SS[3]>0 &&
00230                    (fabs(SS[1])+fabs(SS[3]))> (fabs(history[1])+fabs(history[3])) )
00231                 {
00232                   leverx = 0;
00233                   start = 0;
00234                 }
00235 
00236 
00237               if (SS[1]>0.)
00238                 {
00239                   EBDetId nextPoint;
00240                   if (myCell.ieta()==-1) 
00241                     nextPoint=EBDetId (1,myCell.iphi());
00242                   else 
00243                     {
00244                       int nieta= myCell.ieta()+1;
00245                       if(nieta==86) nieta=85;
00246                       nextPoint=EBDetId(nieta,myCell.iphi());
00247                     }
00248                   if (present(nextPoint))
00249                     myCell = nextPoint;
00250                   else
00251                     leverx = 0;
00252                 }
00253               else
00254                 {
00255                   EBDetId nextPoint;
00256                   if (myCell.ieta()==1) 
00257                     nextPoint=EBDetId(-1,myCell.iphi());
00258                   else 
00259                     {
00260                       int nieta=myCell.ieta()-1;
00261                       if(nieta==-86) nieta=-85;
00262                       nextPoint=EBDetId(nieta,myCell.iphi());
00263                     }
00264                   if (present(nextPoint))
00265                     myCell = nextPoint;
00266                   else
00267                     leverx = 0;
00268                 }
00269             }
00270           
00271           // Update the history. If the point lies between crystals, the closest one
00272           // is returned
00273           history =SS;
00274           
00275           counter++;
00276           if (counter == 10)
00277             {
00278               leverx=0;
00279               levery=0;
00280             }
00281         }
00282       // D.K. if point lies netween cells, take a closest cell.
00283       return DetId(myCell);
00284     }
00285   catch ( cms::Exception &e ) 
00286     { 
00287       return DetId(0);
00288     }
00289 
00290 }
00291 
00292 CaloSubdetectorGeometry::DetIdSet 
00293 EcalBarrelGeometry::getCells( const GlobalPoint& r, 
00294                               double             dR ) const 
00295 {
00296    static const int maxphi ( EBDetId::MAX_IPHI ) ;
00297    static const int maxeta ( EBDetId::MAX_IETA ) ;
00298    CaloSubdetectorGeometry::DetIdSet dis;  // this is the return object
00299 
00300    if( 0.000001 < dR )
00301    {
00302       if( dR > M_PI/2. ) // this version needs "small" dR
00303       {
00304          dis = CaloSubdetectorGeometry::getCells( r, dR ) ; // base class version
00305       }
00306       else
00307       {
00308          const double dR2     ( dR*dR ) ;
00309          const double reta    ( r.eta() ) ;
00310          const double rz      ( r.z()   ) ;
00311          const double rphi    ( r.phi() ) ;
00312          const double lowEta  ( reta - dR ) ;
00313          const double highEta ( reta + dR ) ;
00314          
00315          if( highEta > -1.5 &&
00316              lowEta  <  1.5    ) // in barrel
00317          {
00318             const double scale       ( maxphi/(2*M_PI) ) ; // angle to index
00319             const int    ieta_center ( int( reta*scale + ((rz<0)?(-1):(1))) ) ;
00320             const double phi         ( rphi<0 ? rphi + 2*M_PI : rphi ) ;
00321             const int    iphi_center ( int( phi*scale + 11. ) ) ; // phi=-9.4deg is iphi=1
00322 
00323             const double fr    ( dR*scale    ) ; // # crystal widths in dR
00324             const double frp   ( 1.08*fr + 1. ) ; // conservatively above fr 
00325             const double frm   ( 0.92*fr - 1. ) ; // conservatively below fr
00326             const int    idr   ( (int)frp        ) ; // integerize
00327             const int    idr2p ( (int)(frp*frp)     ) ;
00328             const int    idr2m ( frm > 0 ? int(frm*frm) : 0 ) ;
00329 
00330             for( int de ( -idr ) ; de <= idr ; ++de ) // over eta limits
00331             {
00332                int ieta ( de + ieta_center ) ;
00333                
00334                if( std::abs(ieta) <= maxeta &&
00335                    ieta      != 0         ) // eta is in EB
00336                {
00337                   const int de2 ( de*de ) ;
00338                   for( int dp ( -idr ) ; dp <= idr ; ++dp )  // over phi limits
00339                   {
00340                      const int irange2 ( dp*dp + de2 ) ;
00341                      
00342                      if( irange2 <= idr2p ) // cut off corners that must be too far away
00343                      {
00344                         const int iphi ( ( iphi_center + dp + maxphi - 1 )%maxphi + 1 ) ;
00345                         
00346                         if( iphi != 0 )
00347                         {
00348                            const EBDetId id ( ieta, iphi ) ;
00349                            
00350                            bool ok ( irange2 < idr2m ) ;  // no more calculation necessary if inside this radius
00351                            
00352                            if( !ok ) // if not ok, then we have to test this cell for being inside cone
00353                            {
00354                               const CaloCellGeometry* cell ( getGeometry( id ) );
00355                               if( 0 != cell )
00356                               {
00357                                  const GlobalPoint& p   ( cell->getPosition() ) ;
00358                                  const double       eta ( p.eta() ) ;
00359                                  const double       phi ( p.phi() ) ;
00360                                  ok = ( reco::deltaR2( eta, phi, reta, rphi ) < dR2 ) ;
00361                               }
00362                            }
00363                            if( ok ) dis.insert( id ) ;
00364                         }
00365                      }
00366                   }
00367                }
00368             }
00369          }
00370       }
00371    }
00372    return dis;
00373 }
00374 
00375 const EcalBarrelGeometry::OrderedListOfEEDetId* 
00376 EcalBarrelGeometry::getClosestEndcapCells( EBDetId id ) const
00377 {
00378    OrderedListOfEEDetId* ptr ( 0 ) ;
00379    if( 0 != id.rawId() )
00380    {
00381       const int iPhi     ( id.iphi() ) ;
00382 
00383       const int iz       ( id.ieta()>0 ? 1 : -1 ) ;
00384       const EEDetId eeid ( EEDetId::idOuterRing( iPhi, iz ) ) ;
00385 
00386 //      const int ix ( eeid.ix() ) ;
00387 //      const int iy ( eeid.iy() ) ;
00388 
00389       const int iq ( eeid.iquadrant() ) ;
00390       const int xout ( 1==iq || 4==iq ? 1 : -1 ) ;
00391       const int yout ( 1==iq || 2==iq ? 1 : -1 ) ;
00392       if( 0 == m_borderMgr )
00393       {
00394          m_borderMgr = new EZMgrFL<EEDetId>( 720*9, 9 ) ;
00395       }
00396       if( 0 == m_borderPtrVec )
00397       {
00398          m_borderPtrVec = new VecOrdListEEDetIdPtr() ;
00399          m_borderPtrVec->reserve( 720 ) ;
00400          for( unsigned int i ( 0 ) ; i != 720 ; ++i )
00401          {
00402             const int kz ( 360>i ? -1 : 1 ) ;
00403             const EEDetId eeid ( EEDetId::idOuterRing( i%360+1, kz ) ) ;
00404 
00405             const int jx ( eeid.ix() ) ;
00406             const int jy ( eeid.iy() ) ;
00407 
00408             OrderedListOfEEDetId& olist ( *new OrderedListOfEEDetId( m_borderMgr ) );
00409             int il ( 0 ) ;
00410 
00411             for( unsigned int k ( 1 ) ; k <= 25 ; ++k )
00412             {
00413                const int kx ( 1==k || 2==k || 3==k || 12==k || 13==k ? 0 :
00414                               ( 4==k || 6==k || 8==k || 15==k || 20==k ? 1 :
00415                                 ( 5==k || 7==k || 9==k || 16==k || 19==k ? -1 :
00416                                   ( 10==k || 14==k || 21==k || 22==k || 25==k ? 2 : -2 )))) ;
00417                const int ky ( 1==k || 4==k || 5==k || 10==k || 11==k ? 0 :
00418                               ( 2==k || 6==k || 7==k || 14==k || 17==k ? 1 :
00419                                 ( 3==k || 8==k || 9==k || 18==k || 21==k ? -1 :
00420                                   ( 12==k || 15==k || 16==k || 22==k || 23==k ? 2 : -2 )))) ;
00421 
00422                if( 8>=il && EEDetId::validDetId( jx + kx*xout ,
00423                                                  jy + ky*yout , kz ) ) 
00424                {
00425                   olist[il++]=EEDetId( jx + kx*xout ,
00426                                        jy + ky*yout , kz ) ;
00427                }
00428             }
00429             m_borderPtrVec->push_back( &olist ) ;
00430          }
00431       }
00432       ptr = (*m_borderPtrVec)[ iPhi - 1 + ( 0>iz ? 0 : 360 ) ] ;
00433    }
00434    return ptr ;
00435 }
00436 
00437 void
00438 EcalBarrelGeometry::localCorners( Pt3DVec&        lc  ,
00439                                   const CCGFloat* pv  , 
00440                                   unsigned int    i   ,
00441                                   Pt3D&           ref   )
00442 {
00443    const bool negz ( EBDetId::kSizeForDenseIndexing/2 >  i ) ;
00444    const bool odd  ( 1 == i%2 ) ;
00445 
00446    if( ( ( negz  && !odd ) ||
00447          ( !negz && odd  )    ) )
00448    {
00449       TruncatedPyramid::localCornersReflection( lc, pv, ref ) ;
00450    }
00451    else
00452    {
00453       TruncatedPyramid::localCornersSwap( lc, pv, ref ) ;
00454    }
00455 }
00456 
00457 void
00458 EcalBarrelGeometry::newCell( const GlobalPoint& f1 ,
00459                              const GlobalPoint& f2 ,
00460                              const GlobalPoint& f3 ,
00461                              const CCGFloat*    parm ,
00462                              const DetId&       detId   ) 
00463 {
00464    const unsigned int cellIndex ( EBDetId( detId ).denseIndex() ) ;
00465    m_cellVec[ cellIndex ] =
00466       TruncatedPyramid( cornersMgr(), f1, f2, f3, parm ) ;
00467    m_validIds.push_back( detId ) ;
00468 }
00469 
00470 CCGFloat 
00471 EcalBarrelGeometry::avgRadiusXYFrontFaceCenter() const 
00472 {
00473    if( 0 > m_radius )
00474    {
00475       CCGFloat sum ( 0 ) ;
00476       for( uint32_t i ( 0 ) ; i != m_cellVec.size() ; ++i )
00477       {
00478          const CaloCellGeometry* cell ( cellGeomPtr(i) ) ;
00479          if( 0 != cell )
00480          {
00481             const GlobalPoint& pos ( cell->getPosition() ) ;
00482             sum += pos.perp() ;
00483          }
00484       }
00485       m_radius = sum/m_cellVec.size() ;
00486    }
00487    return m_radius ;
00488 }
00489 
00490 const CaloCellGeometry* 
00491 EcalBarrelGeometry::cellGeomPtr( uint32_t index ) const
00492 {
00493    const CaloCellGeometry* cell ( &m_cellVec[ index ] ) ;
00494    return ( m_cellVec.size() < index ||
00495             0 == cell->param() ? 0 : cell ) ;
00496 }