Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-10 02:21:30

0001 #ifndef TrackPropagation_ConvertFromToCLHEP_h
0002 #define TrackPropagation_ConvertFromToCLHEP_h
0003 
0004 // CLHEP
0005 #include <CLHEP/Geometry/Normal3D.h>
0006 #include <CLHEP/Geometry/Point3D.h>
0007 #include <CLHEP/Units/SystemOfUnits.h>
0008 #include <CLHEP/Vector/Rotation.h>
0009 #include <CLHEP/Vector/ThreeVector.h>
0010 
0011 // CMS
0012 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
0013 #include "DataFormats/GeometrySurface/interface/TkRotation.h"
0014 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0015 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0016 
0017 // Geant4
0018 #include "G4ErrorFreeTrajState.hh"
0019 
0020 /** Utilities to convert among CLHEP and CMS points and vectors
0021  */
0022 
0023 namespace TrackPropagation {
0024   using CLHEP::cm;
0025   /**
0026  Convert a CMS GlobalPoint to a CLHEP HepGeom::Point3D<double>
0027  CMS uses cm while Geant4 uses mm. This is taken into account in the
0028  conversion.
0029  */
0030 
0031   inline HepGeom::Point3D<double> globalPointToHepPoint3D(const GlobalPoint &r) {
0032     return HepGeom::Point3D<double>(r.x() * cm, r.y() * cm, r.z() * cm);
0033   }
0034 
0035   /** Convert a CLHEP HepGeom::Point3D<double>  to a CMS GlobalPoint
0036  CMS uses cms while Geant4 uses mm. This is taken into account in the
0037  conversion.
0038  */
0039   inline GlobalPoint hepPoint3DToGlobalPoint(const HepGeom::Point3D<double> &r) {
0040     return GlobalPoint(r.x() / cm, r.y() / cm, r.z() / cm);
0041   }
0042 
0043   /** Convert a G4double representing a scaler measure ( e.g. Track length ) to
0044  the CMS convention, which is using mm
0045  */
0046   inline double g4doubleToCmsDouble(const G4double &d) { return d / cm; }
0047 
0048   /** Convert a CMS GlobalVector to a CLHEP HepGeom::Normal3D<double>
0049  CMS uses GeV while G4 uses MeV
0050  */
0051   inline HepGeom::Normal3D<double> globalVectorToHepNormal3D(const GlobalVector &p) {
0052     return HepGeom::Normal3D<double>(p.x(), p.y(), p.z());
0053   }
0054 
0055   /** Convert a CLHEP HepGeom::Normal3D<double>  to a CMS GlobalVector
0056  CMS uses GeV while G4 uses MeV
0057  */
0058   inline GlobalVector hepNormal3DToGlobalVector(const HepGeom::Normal3D<double> &p) {
0059     return GlobalVector(p.x(), p.y(), p.z());
0060   }
0061 
0062   /** Convert a CMS GlobalVector to a CLHEP CLHEP::Hep3Vector
0063  */
0064   inline CLHEP::Hep3Vector globalVectorToHep3Vector(const GlobalVector &p) {
0065     return CLHEP::Hep3Vector(p.x(), p.y(), p.z());
0066   }
0067 
0068   /** Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalVector
0069  */
0070   inline GlobalVector hep3VectorToGlobalVector(const CLHEP::Hep3Vector &p) { return GlobalVector(p.x(), p.y(), p.z()); }
0071 
0072   /** Convert a CMS GlobalPoint to a CLHEP CLHEP::Hep3Vector
0073  CMS uses cm while Geant4 uses mm. This is taken into account in the
0074  conversion.
0075  */
0076   inline CLHEP::Hep3Vector globalPointToHep3Vector(const GlobalPoint &r) {
0077     return CLHEP::Hep3Vector(r.x() * cm, r.y() * cm, r.z() * cm);
0078   }
0079 
0080   /** Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalPoint
0081  CMS uses cm while Geant4 uses mm. This is taken into account in the
0082  conversion.
0083  */
0084   inline GlobalPoint hep3VectorToGlobalPoint(const CLHEP::Hep3Vector &v) {
0085     return GlobalPoint(v.x() / cm, v.y() / cm, v.z() / cm);
0086   }
0087 
0088   /** Convert a CMS TkRotation<float> to a CLHEP
0089  * CLHEP::HepRotation=G4RotationMatrix
0090  */
0091   inline CLHEP::HepRotation tkRotationFToHepRotation(const TkRotation<float> &tkr) {
0092     return CLHEP::HepRotation(CLHEP::Hep3Vector(tkr.xx(), tkr.yx(), tkr.zx()),
0093                               CLHEP::Hep3Vector(tkr.xy(), tkr.yy(), tkr.zy()),
0094                               CLHEP::Hep3Vector(tkr.xz(), tkr.yz(), tkr.zz()));
0095   }
0096 
0097   /** Convert a CLHEP CLHEP::Hep3Vector to a CMS GlobalPoint
0098  */
0099   inline TkRotation<float> hepRotationToTkRotationF(const CLHEP::HepRotation &r) {
0100     return TkRotation<float>(r.xx(), r.xy(), r.xz(), r.yx(), r.yy(), r.yz(), r.zx(), r.zy(), r.zz());
0101   }
0102 
0103   /** Convert a G4 Trajectory Error Matrix to the CMS Algebraic Sym Matrix
0104  CMS uses q/p as first parameter, G4 uses 1/p
0105  */
0106 
0107   inline AlgebraicSymMatrix55 g4ErrorTrajErrToAlgebraicSymMatrix55(const G4ErrorTrajErr &e, const int q) {
0108     assert(q != 0);
0109     // From DataFormats/CLHEP/interface/Migration.h
0110     // typedef ROOT::Math::SMatrix<double,5,5,ROOT::Math::MatRepSym<double,5> >
0111     // AlgebraicSymMatrix55;
0112     AlgebraicSymMatrix55 m55;
0113     for (unsigned int i = 0; i < 5; i++)
0114       for (unsigned int j = 0; j < 5; j++) {
0115         m55(i, j) = e(i + 1, j + 1);
0116         if (i == 0)
0117           m55(i, j) = double(q) * m55(i, j);
0118         if (j == 0)
0119           m55(i, j) = double(q) * m55(i, j);
0120       }
0121     return m55;
0122   }
0123 
0124   /** Convert a CMS Algebraic Sym Matrix (for curv error) to a G4 Trajectory Error
0125  * Matrix
0126  */
0127   inline G4ErrorTrajErr algebraicSymMatrix55ToG4ErrorTrajErr(const AlgebraicSymMatrix55 &e, const int q) {
0128     assert(q != 0);
0129     G4ErrorTrajErr g4err(5, 1);
0130     for (unsigned int i = 0; i < 5; i++)
0131       for (unsigned int j = 0; j < 5; j++) {
0132         g4err(i + 1, j + 1) = e(i, j);
0133 
0134         if (i == 0)
0135           g4err(i + 1, j + 1) = g4err(i + 1, j + 1) * double(q);
0136         if (j == 0)
0137           g4err(i + 1, j + 1) = g4err(i + 1, j + 1) * double(q);
0138       }
0139     return g4err;
0140   }
0141 
0142 }  // namespace TrackPropagation
0143 
0144 #endif