CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Static Public Member Functions | Static Public Attributes
HICaloUtil Class Reference

#include <HICaloUtil.h>

Public Member Functions

 HICaloUtil ()
 

Static Public Member Functions

static double EcalEta (const reco::Candidate &p)
 
static double EcalEta (double EtaParticle, double Zvertex, double plane_Radius)
 
static double EcalPhi (const reco::Candidate &p)
 
static double EcalPhi (double PtParticle, double EtaParticle, double PhiParticle, int ChargeParticle, double Rstart)
 

Static Public Attributes

static const double kEEtaBarrelEndcap = 1.479
 
static const double kER_ECAL = 136.5
 
static const double kERBARM = 1.357
 
static const double kEZ_Endcap = 328.0
 
static const double kEZENDM = 3.186
 

Detailed Description

Definition at line 15 of file HICaloUtil.h.

Constructor & Destructor Documentation

HICaloUtil::HICaloUtil ( )
inline

Definition at line 17 of file HICaloUtil.h.

17 {}

Member Function Documentation

double HICaloUtil::EcalEta ( const reco::Candidate p)
static

Definition at line 24 of file HICaloUtil.cc.

References reco::Candidate::eta(), mathSSE::sqrt(), and reco::Candidate::vertex().

25 {
26  // Calculate the eta of the particle wrt the center of the detector.
27 
28  double rhovert = sqrt(p.vertex().X()*p.vertex().X()+p.vertex().Y()*p.vertex().Y());
29  double zvert = p.vertex().Z();
30  return EcalEta(p.eta(), zvert, rhovert);
31 }
static double EcalEta(const reco::Candidate &p)
Definition: HICaloUtil.cc:24
virtual float eta() const =0
momentum pseudorapidity
T sqrt(T t)
Definition: SSEVec.h:48
virtual const Point & vertex() const =0
vertex position
double HICaloUtil::EcalEta ( double  EtaParticle,
double  Zvertex,
double  plane_Radius 
)
static

Definition at line 35 of file HICaloUtil.cc.

References ETA, fff_deleter::log, PI, and funct::tan().

36 {
37  // Calculate the eta of the particle wrt center of detector.
38 
39  if(EtaParticle != 0.) {
40  double Theta = 0.0 ;
41  double ZEcal = (kER_ECAL-plane_Radius)*sinh(EtaParticle)+Zvertex;
42 
43  if(ZEcal != 0.0) Theta = atan(kER_ECAL/ZEcal);
44  if(Theta<0.0) Theta = Theta+PI ;
45 
46  double ETA = - log(tan(0.5*Theta));
47 
48  if( fabs(ETA) > kEEtaBarrelEndcap )
49  {
50  double Zend = kEZ_Endcap ;
51  if(EtaParticle<0.0 ) Zend = -Zend ;
52  double Zlen = Zend - Zvertex ;
53  double RR = Zlen/sinh(EtaParticle);
54  Theta = atan((RR+plane_Radius)/Zend);
55  if(Theta<0.0) Theta = Theta+PI ;
56  ETA = - log(tan(0.5*Theta));
57  }
58 
59  return ETA;
60  }
61 
62  return EtaParticle;
63 }
static const double kER_ECAL
Definition: HICaloUtil.h:26
#define PI
#define ETA
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
static const double kEEtaBarrelEndcap
Definition: HICaloUtil.h:25
static const double kEZ_Endcap
Definition: HICaloUtil.h:27
double HICaloUtil::EcalPhi ( const reco::Candidate p)
static

Definition at line 67 of file HICaloUtil.cc.

References reco::Candidate::charge(), reco::Candidate::eta(), phi, reco::Candidate::phi(), reco::Candidate::pt(), mathSSE::sqrt(), and reco::Candidate::vertex().

68 {
69  // Calculate the phi of the particle transported to the Ecal.
70 
71  int charge = p.charge();
72  double phi = p.phi();
73  double rhovert = sqrt(p.vertex().X()*p.vertex().X()+p.vertex().Y()*p.vertex().Y());
74  if(charge) phi = EcalPhi(p.pt(), p.eta(), phi, charge, rhovert);
75  return phi;
76 }
virtual float eta() const =0
momentum pseudorapidity
virtual float phi() const =0
momentum azimuthal angle
virtual float pt() const =0
transverse momentum
T sqrt(T t)
Definition: SSEVec.h:48
virtual const Point & vertex() const =0
vertex position
virtual int charge() const =0
electric charge
static double EcalPhi(const reco::Candidate &p)
Definition: HICaloUtil.cc:67
Definition: DDAxes.h:10
double HICaloUtil::EcalPhi ( double  PtParticle,
double  EtaParticle,
double  PhiParticle,
int  ChargeParticle,
double  Rstart 
)
static

Definition at line 80 of file HICaloUtil.cc.

References PHI, and PI.

82 {
83  // Calculate the phi of the particle transported to the Ecal.
84 
85  double Rbend = kERBARM-(Rstart/100.); //assumed rstart in cm
86  double Bend = 0.3 * 4. * Rbend/ 2.0 ;
87 
88  //PHI correction
89  double PHI = 0.0 ;
90  if( fabs(EtaParticle) <= kEEtaBarrelEndcap) {
91  if (fabs(Bend/PtParticle)<=1.) {
92  PHI = PhiParticle - asin(Bend/PtParticle)*ChargeParticle;
93  if(PHI > PI) {PHI = PHI - 2*PI;}
94  if(PHI < -PI) {PHI = PHI + 2*PI;}
95  return PHI;
96  }
97  } else {
98  double Rhit = 0.0 ;
99  Rhit = kEZENDM / sinh(fabs(EtaParticle));
100  if (fabs(((Rhit-(Rstart/100.))/Rbend)*Bend/PtParticle)<=1.) {
101  PHI = PhiParticle - asin(((Rhit-(Rstart/100.))/Rbend)*Bend/PtParticle)*ChargeParticle;
102  if(PHI > PI) {PHI = PHI - 2*PI;}
103  if(PHI < -PI) {PHI = PHI + 2*PI;}
104  return PHI;
105  } else {
106  return PhiParticle;
107  }
108  }
109 
110  return PhiParticle;
111 }
#define PI
#define PHI
static const double kEEtaBarrelEndcap
Definition: HICaloUtil.h:25
static const double kERBARM
Definition: HICaloUtil.h:28
static const double kEZENDM
Definition: HICaloUtil.h:29

Member Data Documentation

const double HICaloUtil::kEEtaBarrelEndcap = 1.479
static

Definition at line 25 of file HICaloUtil.h.

const double HICaloUtil::kER_ECAL = 136.5
static

Definition at line 26 of file HICaloUtil.h.

const double HICaloUtil::kERBARM = 1.357
static

Definition at line 28 of file HICaloUtil.h.

const double HICaloUtil::kEZ_Endcap = 328.0
static

Definition at line 27 of file HICaloUtil.h.

const double HICaloUtil::kEZENDM = 3.186
static

Definition at line 29 of file HICaloUtil.h.