CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoEgamma/EgammaTools/src/ECALPositionCalculator.cc

Go to the documentation of this file.
00001 
00002 #include "RecoEgamma/EgammaTools/interface/ECALPositionCalculator.h"
00003 
00004 // Framework includes
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "DataFormats/GeometryVector/interface/Pi.h"
00007 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00008 #include "MagneticField/Engine/interface/MagneticField.h"
00009 
00010 static const float R_ECAL           = 136.5;
00011 static const float Z_Endcap         = 328.0;
00012 static const float etaBarrelEndcap  = 1.479;
00013 
00014 double ECALPositionCalculator::ecalPhi(const MagneticField *magField,const math::XYZVector &momentum, const math::XYZPoint &vertex, int charge)
00015 {
00016 
00017    // Get kinematic variables
00018 
00019    float ptParticle = momentum.Rho();
00020    float etaParticle = momentum.eta();
00021    float phiParticle = momentum.phi();
00022    float vRho = vertex.Rho();
00023 
00024    // Magnetic field
00025 
00026    const float RBARM = 1.357 ;          // was 1.31 , updated on 16122003
00027    const float ZENDM = 3.186 ;          // was 3.15 , updated on 16122003
00028 
00029    float rbend = RBARM-(vRho/100.0);    //Assumed vRho in cm
00030    float bend  = 0.3 * magField->inTesla(GlobalPoint(0.,0.,0.)).z() * rbend / 2.0; 
00031    float phi = 0.0;
00032 
00033    if( fabs(etaParticle) <=  etaBarrelEndcap)
00034    {
00035       if (fabs(bend/ptParticle) <= 1.) 
00036       {
00037          phi = phiParticle - asin(bend/ptParticle)*charge;
00038          if(phi >  Geom::pi()) phi = phi - Geom::twoPi();
00039          if(phi < -Geom::pi()) phi = phi + Geom::twoPi();
00040       } else {
00041          edm::LogWarning("") << "[EcalPositionFromTrack::phiTransformation] Warning: "
00042                              << "Too low Pt, giving up";
00043          return phiParticle;
00044       }
00045 
00046    } // end if in the barrel
00047   
00048    if(fabs(etaParticle) > etaBarrelEndcap)
00049    {
00050       float rHit = 0.0;
00051       rHit = ZENDM / sinh(fabs(etaParticle));
00052       if (fabs(((rHit-(vRho/100.0))/rbend)*bend/ptParticle) <= 1.0)
00053       {
00054          phi = phiParticle - asin(((rHit-(vRho/100.0)) / rbend)*bend/ptParticle)*charge;
00055          if(phi >  Geom::pi()) phi = phi - Geom::twoPi();
00056          if(phi < -Geom::pi()) phi = phi + Geom::twoPi();
00057       } else {
00058          edm::LogWarning("") << "[EcalPositionFromTrack::phiTransformation] Warning: "
00059                              << "Too low Pt, giving up";
00060          return phiParticle;
00061       }  
00062 
00063     } // end if in the endcap
00064   
00065    return phi;
00066 
00067 }
00068 
00069 double ECALPositionCalculator::ecalEta(const math::XYZVector &momentum, const math::XYZPoint &vertex)
00070 {
00071 
00072    // Get kinematic variables
00073 
00074    float etaParticle = momentum.eta();
00075    float vZ = vertex.z();
00076    float vRho = vertex.Rho();
00077 
00078    if (etaParticle != 0.0)
00079    {
00080       float theta = 0.0;
00081       float zEcal = (R_ECAL-vRho)*sinh(etaParticle)+vZ;
00082       
00083       if(zEcal != 0.0) theta = atan(R_ECAL/zEcal);
00084       if(theta < 0.0) theta = theta + Geom::pi() ;
00085 
00086       float ETA = - log(tan(0.5*theta));
00087       
00088       if( fabs(ETA) > etaBarrelEndcap )
00089       {
00090          float Zend = Z_Endcap ;
00091          if(etaParticle < 0.0 )  Zend = - Zend;
00092          float Zlen = Zend - vZ ;
00093          float RR = Zlen/sinh(etaParticle);
00094          theta = atan((RR+vRho)/Zend);
00095          if(theta < 0.0) theta = theta+Geom::pi();
00096          ETA = -log(tan(0.5*theta));
00097       }
00098       return ETA;
00099 
00100     } else {
00101        edm::LogWarning("")  << "[EcalPositionFromTrack::etaTransformation] Warning: "
00102                             << "Eta equals to zero, not correcting";
00103        return etaParticle;
00104     }
00105 
00106 }
00107