CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CaloVectors.h
Go to the documentation of this file.
2 
4 {
5  double pt = sc->energy()/cosh(sc->eta());
6  math::XYZTLorentzVector detVec(pt*cos(sc->phi()), pt*sin(sc->phi()), pt*sinh(sc->eta()), sc->energy());
7  return detVec;
8 }
10 {
11  double pt = sc.energy()/cosh(sc.eta());
12  math::XYZTLorentzVector detVec(pt*cos(sc.phi()), pt*sin(sc.phi()), pt*sinh(sc.eta()), sc.energy());
13  return detVec;
14 }
16 {
17  math::XYZPoint hitPos(pos.x(), pos.y(), pos.z());
18  math::XYZVector Vec = hitPos - vertex;
19  double eta = Vec.Eta();
20  double phi = Vec.Phi();
21  double pt = energy/cosh(eta);
22  math::XYZTLorentzVector detVec(pt*cos(phi), pt*sin(phi), pt*sinh(eta), energy);
23  return detVec;
24 }
26 {
27  math::XYZVector Vec = sc.position() - vertex;
28  double eta = Vec.Eta();
29  double phi = Vec.Phi();
30  double pt = sc.energy()/cosh(eta);
31  math::XYZTLorentzVector probe(pt*cos(phi), pt*sin(phi), pt*sinh(eta), sc.energy());
32  return probe;
33 }
35 {
36  math::XYZVector Vec = sc.position() - vertex;
37  double eta = Vec.Eta();
38  double phi = Vec.Phi();
39  double pt = sc.rawEnergy()/cosh(eta);
40  math::XYZTLorentzVector probe(pt*cos(phi), pt*sin(phi), pt*sinh(eta), sc.rawEnergy());
41  return probe;
42 }
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:123
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T y() const
Definition: PV3DBase.h:63
math::XYZTLorentzVector DetectorVector(const reco::SuperClusterRef &sc)
Definition: CaloVectors.h:3
T eta() const
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:160
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double energy() const
cluster energy
Definition: CaloCluster.h:120
math::XYZTLorentzVector PhysicsVectorRaw(const math::XYZPoint &vertex, const reco::SuperCluster &sc)
Definition: CaloVectors.h:34
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:49
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:163
T x() const
Definition: PV3DBase.h:62
math::XYZTLorentzVector PhysicsVector(const math::XYZPoint &vertex, const reco::SuperCluster &sc)
Definition: CaloVectors.h:25
Definition: DDAxes.h:10