CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions
TrackPropagation Namespace Reference

Functions

G4ErrorTrajErr algebraicSymMatrix55ToG4ErrorTrajErr (const AlgebraicSymMatrix55 &e, const int q)
 
AlgebraicSymMatrix55 g4ErrorTrajErrToAlgebraicSymMatrix55 (const G4ErrorTrajErr &e, const int q)
 
CLHEP::Hep3Vector globalPointToHep3Vector (const GlobalPoint &r)
 
HepGeom::Point3D< double > globalPointToHepPoint3D (const GlobalPoint &r)
 
CLHEP::Hep3Vector globalVectorToHep3Vector (const GlobalVector &p)
 
HepGeom::Normal3D< double > globalVectorToHepNormal3D (const GlobalVector &p)
 
GlobalPoint hep3VectorToGlobalPoint (const CLHEP::Hep3Vector &v)
 
GlobalVector hep3VectorToGlobalVector (const CLHEP::Hep3Vector &p)
 
GlobalVector hepNormal3DToGlobalVector (const HepGeom::Normal3D< double > &p)
 
GlobalPoint hepPoint3DToGlobalPoint (const HepGeom::Point3D< double > &r)
 
TkRotation< float > hepRotationToTkRotationF (const CLHEP::HepRotation &r)
 
CLHEP::HepRotation tkRotationFToHepRotation (const TkRotation< float > &tkr)
 

Detailed Description

Utilities to convert among CLHEP and CMS points and vectors

Function Documentation

G4ErrorTrajErr TrackPropagation::algebraicSymMatrix55ToG4ErrorTrajErr ( const AlgebraicSymMatrix55 e,
const int  q 
)

Convert a CMS Algebraic Sym Matrix (for curv error) to a G4 Trajectory Error Matrix

Definition at line 141 of file ConvertFromToCLHEP.h.

References alignCSCRings::e, i, and j.

141  {
142  G4ErrorTrajErr g4err(5,1);
143  for (unsigned int i = 0; i < 5; i++)
144  for (unsigned int j = 0; j < 5; j++) {
145  g4err(i+1, j+1) = e(i,j);
146  if(i==0) g4err(i+1,j+1) = q*g4err(i+1,j+1);
147  if(j==0) g4err(i+1,j+1) = q*g4err(i+1,j+1);
148  }
149  return g4err;
150  }
int i
Definition: DBlmapReader.cc:9
int j
Definition: DBlmapReader.cc:9
AlgebraicSymMatrix55 TrackPropagation::g4ErrorTrajErrToAlgebraicSymMatrix55 ( const G4ErrorTrajErr &  e,
const int  q 
)

Convert a G4 Trajectory Error Matrix to the CMS Algebraic Sym Matrix CMS uses q/p as first parameter, G4 uses 1/p

Definition at line 125 of file ConvertFromToCLHEP.h.

References alignCSCRings::e, i, and j.

125  {
126  //From DataFormats/CLHEP/interface/Migration.h
127  //typedef ROOT::Math::SMatrix<double,5,5,ROOT::Math::MatRepSym<double,5> > AlgebraicSymMatrix55;
129  for (unsigned int i = 0; i < 5; i++)
130  for (unsigned int j = 0; j < 5; j++) {
131  m55(i, j) = e(i+1,j+1);
132  if(i==0) m55(i,j) = q*m55(i,j);
133  if(j==0) m55(i,j) = q*m55(i,j);
134  }
135  return m55;
136  }
int i
Definition: DBlmapReader.cc:9
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
int j
Definition: DBlmapReader.cc:9
CLHEP::Hep3Vector TrackPropagation::globalPointToHep3Vector ( const GlobalPoint r)

Convert a CMS GlobalPoint to a CLHEP CLHEP::Hep3Vector CMS uses cm while Geant4 uses mm. This is taken into account in the conversion.

Definition at line 82 of file ConvertFromToCLHEP.h.

References PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

82  {
83  return CLHEP::Hep3Vector(r.x()*cm, r.y()*cm, r.z()*cm);
84  }
T y() const
Definition: PV3DBase.h:63
T z() const
Definition: PV3DBase.h:64
T x() const
Definition: PV3DBase.h:62
HepGeom::Point3D<double> TrackPropagation::globalPointToHepPoint3D ( const GlobalPoint r)

Convert a CMS GlobalPoint to a CLHEP HepGeom::Point3D<double> CMS uses cm while Geant4 uses mm. This is taken into account in the conversion.

Definition at line 31 of file ConvertFromToCLHEP.h.

References PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

31  {
32  return HepGeom::Point3D<double> (r.x()*cm, r.y()*cm, r.z()*cm);
33  }
T y() const
Definition: PV3DBase.h:63
T z() const
Definition: PV3DBase.h:64
T x() const
Definition: PV3DBase.h:62
CLHEP::Hep3Vector TrackPropagation::globalVectorToHep3Vector ( const GlobalVector p)

Convert a CMS GlobalVector to a CLHEP CLHEP::Hep3Vector

Definition at line 65 of file ConvertFromToCLHEP.h.

References PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

65  {
66  return CLHEP::Hep3Vector(p.x(), p.y(), p.z());
67  }
T y() const
Definition: PV3DBase.h:63
T z() const
Definition: PV3DBase.h:64
T x() const
Definition: PV3DBase.h:62
HepGeom::Normal3D<double> TrackPropagation::globalVectorToHepNormal3D ( const GlobalVector p)

Convert a CMS GlobalVector to a CLHEP HepGeom::Normal3D<double> CMS uses GeV while G4 uses MeV

Definition at line 49 of file ConvertFromToCLHEP.h.

References PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

49  {
50  return HepGeom::Normal3D<double> (p.x(), p.y(), p.z());
51  }
T y() const
Definition: PV3DBase.h:63
T z() const
Definition: PV3DBase.h:64
T x() const
Definition: PV3DBase.h:62
GlobalPoint TrackPropagation::hep3VectorToGlobalPoint ( const CLHEP::Hep3Vector &  v)

Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalPoint CMS uses cm while Geant4 uses mm. This is taken into account in the conversion.

Definition at line 90 of file ConvertFromToCLHEP.h.

90  {
91  return GlobalPoint(v.x()/cm, v.y()/cm, v.z()/cm);
92  }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
GlobalVector TrackPropagation::hep3VectorToGlobalVector ( const CLHEP::Hep3Vector &  p)

Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalVector

Definition at line 71 of file ConvertFromToCLHEP.h.

71  {
72  return GlobalVector(p.x(), p.y(), p.z());
73  }
Global3DVector GlobalVector
Definition: GlobalVector.h:10
GlobalVector TrackPropagation::hepNormal3DToGlobalVector ( const HepGeom::Normal3D< double > &  p)

Convert a CLHEP HepGeom::Normal3D<double> to a CMS GlobalVector CMS uses GeV while G4 uses MeV

Definition at line 56 of file ConvertFromToCLHEP.h.

56  {
57  return GlobalVector(p.x(), p.y(), p.z());
58  }
Global3DVector GlobalVector
Definition: GlobalVector.h:10
GlobalPoint TrackPropagation::hepPoint3DToGlobalPoint ( const HepGeom::Point3D< double > &  r)

Convert a CLHEP HepGeom::Point3D<double> to a CMS GlobalPoint CMS uses cms while Geant4 uses mm. This is taken into account in the conversion.

Definition at line 41 of file ConvertFromToCLHEP.h.

41  {
42  return GlobalPoint(r.x()/cm, r.y()/cm, r.z()/cm);
43  }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
TkRotation<float> TrackPropagation::hepRotationToTkRotationF ( const CLHEP::HepRotation &  r)

Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalPoint

Definition at line 112 of file ConvertFromToCLHEP.h.

112  {
113  return TkRotation<float>(r.xx(), r.xy(), r.xz(),
114  r.yx(), r.yy(), r.yz(),
115  r.zx(), r.zy(), r.zz());
116  }
CLHEP::HepRotation TrackPropagation::tkRotationFToHepRotation ( const TkRotation< float > &  tkr)

Convert a CMS TkRotation<float> to a CLHEP CLHEP::HepRotation=G4RotationMatrix

Definition at line 104 of file ConvertFromToCLHEP.h.

References TkRotation< T >::xx(), TkRotation< T >::xy(), TkRotation< T >::xz(), TkRotation< T >::yx(), TkRotation< T >::yy(), TkRotation< T >::yz(), TkRotation< T >::zx(), TkRotation< T >::zy(), and TkRotation< T >::zz().

104  {
105  return CLHEP::HepRotation(CLHEP::Hep3Vector(tkr.xx(),tkr.yx(), tkr.zx()),
106  CLHEP::Hep3Vector(tkr.xy(),tkr.yy(), tkr.zy()),
107  CLHEP::Hep3Vector(tkr.xz(),tkr.yz(), tkr.zz()));
108  }
T xx() const
T yx() const
T zx() const
T xy() const
T zz() const
T zy() const
T yy() const
T xz() const
T yz() const