CMS 3D CMS Logo

Functions
TrackPropagation Namespace Reference

Functions

G4ErrorTrajErr algebraicSymMatrix55ToG4ErrorTrajErr (const AlgebraicSymMatrix55 &e, const int q)
 
double g4doubleToCmsDouble (const G4double &d)
 
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 
)
inline

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

Definition at line 146 of file ConvertFromToCLHEP.h.

References MillePedeFileConverter_cfg::e, i, and j.

Referenced by Geant4ePropagator::propagateGeneric().

148 {
149  assert(q != 0);
150  G4ErrorTrajErr g4err(5, 1);
151  for (unsigned int i = 0; i < 5; i++)
152  for (unsigned int j = 0; j < 5; j++)
153  {
154  g4err(i + 1, j + 1) = e(i, j);
155 
156  if (i == 0)
157  g4err(i + 1, j + 1) = g4err(i + 1, j + 1) * double(q);
158  if (j == 0)
159  g4err(i + 1, j + 1) = g4err(i + 1, j + 1) * double(q);
160 
161  }
162  return g4err;
163 }
int i
Definition: DBlmapReader.cc:9
int j
Definition: DBlmapReader.cc:9
double TrackPropagation::g4doubleToCmsDouble ( const G4double &  d)
inline

Convert a G4double representing a scaler measure ( e.g. Track length ) to the CMS convention, which is using mm

Definition at line 48 of file ConvertFromToCLHEP.h.

Referenced by Geant4ePropagator::propagateGeneric().

49 {
50  return d / cm;
51 }
AlgebraicSymMatrix55 TrackPropagation::g4ErrorTrajErrToAlgebraicSymMatrix55 ( const G4ErrorTrajErr &  e,
const int  q 
)
inline

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 124 of file ConvertFromToCLHEP.h.

References MillePedeFileConverter_cfg::e, i, and j.

Referenced by Geant4ePropagator::propagateGeneric().

126 {
127  assert(q != 0);
128  //From DataFormats/CLHEP/interface/Migration.h
129  //typedef ROOT::Math::SMatrix<double,5,5,ROOT::Math::MatRepSym<double,5> > AlgebraicSymMatrix55;
131  for (unsigned int i = 0; i < 5; i++)
132  for (unsigned int j = 0; j < 5; j++)
133  {
134  m55(i, j) = e(i + 1, j + 1);
135  if (i == 0)
136  m55(i, j) = double(q) * m55(i, j);
137  if (j == 0)
138  m55(i, j) = double(q) * m55(i, j);
139 
140  }
141  return m55;
142 }
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)
inline

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 89 of file ConvertFromToCLHEP.h.

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

Referenced by Geant4ePropagator::propagateGeneric(), and Geant4ePropagator::transformToG4SurfaceTarget().

90 {
91  return CLHEP::Hep3Vector(r.x() * cm, r.y() * cm, r.z() * cm);
92 }
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)
inline

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().

Referenced by Geant4ePropagator::transformToG4SurfaceTarget().

32 {
33  return HepGeom::Point3D<double>(r.x() * cm, r.y() * cm, r.z() * cm);
34 }
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)
inline

Convert a CMS GlobalVector to a CLHEP CLHEP::Hep3Vector

Definition at line 73 of file ConvertFromToCLHEP.h.

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

Referenced by Geant4ePropagator::propagateGeneric().

74 {
75  return CLHEP::Hep3Vector(p.x(), p.y(), p.z());
76 }
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)
inline

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

Definition at line 56 of file ConvertFromToCLHEP.h.

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

Referenced by Geant4ePropagator::transformToG4SurfaceTarget().

58 {
59  return HepGeom::Normal3D<double>(p.x(), p.y(), p.z());
60 }
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)
inline

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 98 of file ConvertFromToCLHEP.h.

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

Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalVector

Definition at line 80 of file ConvertFromToCLHEP.h.

Referenced by Geant4ePropagator::propagateGeneric().

81 {
82  return GlobalVector(p.x(), p.y(), p.z());
83 }
Global3DVector GlobalVector
Definition: GlobalVector.h:10
GlobalVector TrackPropagation::hepNormal3DToGlobalVector ( const HepGeom::Normal3D< double > &  p)
inline

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

Definition at line 65 of file ConvertFromToCLHEP.h.

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

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 40 of file ConvertFromToCLHEP.h.

Referenced by Geant4ePropagator::propagateGeneric().

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)
inline

Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalPoint

Definition at line 114 of file ConvertFromToCLHEP.h.

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

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

Definition at line 105 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().

Referenced by Geant4ePropagator::transformToG4SurfaceTarget().

106 {
107  return CLHEP::HepRotation(CLHEP::Hep3Vector(tkr.xx(), tkr.yx(), tkr.zx()),
108  CLHEP::Hep3Vector(tkr.xy(), tkr.yy(), tkr.zy()),
109  CLHEP::Hep3Vector(tkr.xz(), tkr.yz(), tkr.zz()));
110 }
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