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

◆ algebraicSymMatrix55ToG4ErrorTrajErr()

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

References cms::cuda::assert(), MillePedeFileConverter_cfg::e, mps_fire::i, dqmiolumiharvest::j, and submitPVResolutionJobs::q.

Referenced by Geant4ePropagator::propagateGeneric().

126  {
127  assert(q != 0);
128  G4ErrorTrajErr g4err(5, 1);
129  for (unsigned int i = 0; i < 5; i++)
130  for (unsigned int j = 0; j < 5; j++) {
131  g4err(i + 1, j + 1) = e(i, j);
132 
133  if (i == 0)
134  g4err(i + 1, j + 1) = g4err(i + 1, j + 1) * double(q);
135  if (j == 0)
136  g4err(i + 1, j + 1) = g4err(i + 1, j + 1) * double(q);
137  }
138  return g4err;
139  }
assert(be >=bs)

◆ g4doubleToCmsDouble()

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

References ztail::d.

Referenced by Geant4ePropagator::propagateGeneric().

45 { return d / cm; }
d
Definition: ztail.py:151

◆ g4ErrorTrajErrToAlgebraicSymMatrix55()

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

References cms::cuda::assert(), MillePedeFileConverter_cfg::e, mps_fire::i, dqmiolumiharvest::j, and submitPVResolutionJobs::q.

Referenced by Geant4ePropagator::propagateGeneric().

106  {
107  assert(q != 0);
108  // From DataFormats/CLHEP/interface/Migration.h
109  // typedef ROOT::Math::SMatrix<double,5,5,ROOT::Math::MatRepSym<double,5> >
110  // AlgebraicSymMatrix55;
112  for (unsigned int i = 0; i < 5; i++)
113  for (unsigned int j = 0; j < 5; j++) {
114  m55(i, j) = e(i + 1, j + 1);
115  if (i == 0)
116  m55(i, j) = double(q) * m55(i, j);
117  if (j == 0)
118  m55(i, j) = double(q) * m55(i, j);
119  }
120  return m55;
121  }
assert(be >=bs)
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55

◆ globalPointToHep3Vector()

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

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

75  {
76  return CLHEP::Hep3Vector(r.x() * cm, r.y() * cm, r.z() * cm);
77  }

◆ globalPointToHepPoint3D()

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

Referenced by Geant4ePropagator::transformToG4SurfaceTarget().

30  {
31  return HepGeom::Point3D<double>(r.x() * cm, r.y() * cm, r.z() * cm);
32  }

◆ globalVectorToHep3Vector()

CLHEP::Hep3Vector TrackPropagation::globalVectorToHep3Vector ( const GlobalVector p)
inline

Convert a CMS GlobalVector to a CLHEP CLHEP::Hep3Vector

Definition at line 63 of file ConvertFromToCLHEP.h.

References AlCaHLTBitMon_ParallelJobs::p.

Referenced by Geant4ePropagator::propagateGeneric().

63  {
64  return CLHEP::Hep3Vector(p.x(), p.y(), p.z());
65  }

◆ globalVectorToHepNormal3D()

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

References AlCaHLTBitMon_ParallelJobs::p.

Referenced by Geant4ePropagator::transformToG4SurfaceTarget().

50  {
51  return HepGeom::Normal3D<double>(p.x(), p.y(), p.z());
52  }

◆ hep3VectorToGlobalPoint()

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

References findQualityFiles::v.

83  {
84  return GlobalPoint(v.x() / cm, v.y() / cm, v.z() / cm);
85  }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10

◆ hep3VectorToGlobalVector()

GlobalVector TrackPropagation::hep3VectorToGlobalVector ( const CLHEP::Hep3Vector &  p)
inline

Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalVector

Definition at line 69 of file ConvertFromToCLHEP.h.

References AlCaHLTBitMon_ParallelJobs::p.

Referenced by Geant4ePropagator::propagateGeneric().

69 { return GlobalVector(p.x(), p.y(), p.z()); }
Global3DVector GlobalVector
Definition: GlobalVector.h:10

◆ hepNormal3DToGlobalVector()

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

References AlCaHLTBitMon_ParallelJobs::p.

57  {
58  return GlobalVector(p.x(), p.y(), p.z());
59  }
Global3DVector GlobalVector
Definition: GlobalVector.h:10

◆ hepPoint3DToGlobalPoint()

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

Referenced by Geant4ePropagator::propagateGeneric().

38  {
39  return GlobalPoint(r.x() / cm, r.y() / cm, r.z() / cm);
40  }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10

◆ hepRotationToTkRotationF()

TkRotation<float> TrackPropagation::hepRotationToTkRotationF ( const CLHEP::HepRotation &  r)
inline

Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalPoint

Definition at line 98 of file ConvertFromToCLHEP.h.

98  {
99  return TkRotation<float>(r.xx(), r.xy(), r.xz(), r.yx(), r.yy(), r.yz(), r.zx(), r.zy(), r.zz());
100  }

◆ tkRotationFToHepRotation()

CLHEP::HepRotation TrackPropagation::tkRotationFToHepRotation ( const TkRotation< float > &  tkr)
inline

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

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

90  {
91  return CLHEP::HepRotation(CLHEP::Hep3Vector(tkr.xx(), tkr.yx(), tkr.zx()),
92  CLHEP::Hep3Vector(tkr.xy(), tkr.yy(), tkr.zy()),
93  CLHEP::Hep3Vector(tkr.xz(), tkr.yz(), tkr.zz()));
94  }
T xx() const
T xy() const
T zz() const
T yy() const
T yz() const
T zx() const
T zy() const
T xz() const
T yx() const