CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
TwoTrackMinimumDistanceLineLine Class Reference

#include <TwoTrackMinimumDistanceLineLine.h>

Public Member Functions

bool calculate (const GlobalTrajectoryParameters &, const GlobalTrajectoryParameters &)
 
double firstAngle () const
 
std::pair< double, double > pathLength () const
 
std::pair< GlobalPoint,
GlobalPoint
points () const
 
double secondAngle () const
 

Private Attributes

GlobalPoint gPos
 
GlobalPoint hPos
 
double pathG
 
double pathH
 
double phiG
 
double phiH
 

Detailed Description

This is a helper class for TwoTrackMinimumDistance. No user should need direct access to this class. Exact solution.

Definition at line 18 of file TwoTrackMinimumDistanceLineLine.h.

Member Function Documentation

bool TwoTrackMinimumDistanceLineLine::calculate ( const GlobalTrajectoryParameters theG,
const GlobalTrajectoryParameters theH 
)

Calculates the point of closest approach on the two tracks.

Returns
false in case of success
true in case of failure. Possible failures:
  • Either of the Trajectories are charged
  • Either of the Trajectories is of zero momentum
  • Trajectories are parallel

Definition at line 8 of file TwoTrackMinimumDistanceLineLine.cc.

References GlobalTrajectoryParameters::charge(), Vector3DBase< T, FrameTag >::dot(), gPos, hPos, MagneticField::inTesla(), PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), GlobalTrajectoryParameters::momentum(), pathG, pathH, PV3DBase< T, PVType, FrameType >::phi(), phiG, phiH, GlobalTrajectoryParameters::position(), and PV3DBase< T, PVType, FrameType >::z().

11 {
12  GlobalPoint gOrig = theG.position();
13  GlobalPoint hOrig = theH.position();
14  bool isLineG = theG.magneticField().inTesla(gOrig).z() == 0 || theG.charge() == 0;
15  bool isLineH = theH.magneticField().inTesla(hOrig).z() == 0 || theH.charge() == 0;
16  if ( ! (isLineG && isLineH) )
17  {
18  edm::LogWarning ("TwoTrackMinimumDistanceLineLine")
19  << "charge of input track is not zero or field non zero"
20  << "\n positions: "<<gOrig<<" , "<<hOrig
21  << "\n Bz fields: "<<theG.magneticField().inTesla(gOrig).z()<<" , "<<theH.magneticField().inTesla(hOrig).z()
22  << "\n charges: "<<theG.charge()<<" , "<<theH.charge();
23  return true;
24  };
25 
26  GlobalVector gVec = theG.momentum();
27  const double gMag= theG.momentum().mag(); double gMag2 = gMag*gMag;
28  GlobalVector hVec = theH.momentum();
29  const double hMag= theH.momentum().mag(); double hMag2 = hMag*hMag;
30 
31  if ( gMag == 0. || hMag == 0. )
32  {
33  edm::LogWarning ("TwoTrackMinimumDistanceLineLine")
34  << "momentum of input trajectory is zero.";
35  return true;
36  };
37 
38  phiG = theG.momentum().phi();
39  phiH = theH.momentum().phi();
40 
41  double gVec_Dot_hVec = gVec.dot(hVec);
42  double norm = gVec_Dot_hVec*gVec_Dot_hVec - gMag2*hMag2;
43 
44  if ( norm == 0 )
45  {
46  edm::LogWarning ("TwoTrackMinimumDistanceLineLine")
47  << "Tracks are parallel.";
48  return true;
49  }
50 
51  GlobalVector posDiff = gOrig - hOrig;
52 
53  double tG = (posDiff.dot(gVec) * hMag2 - gVec_Dot_hVec * posDiff.dot(hVec))/ norm;
54  double tH = (gVec_Dot_hVec * posDiff.dot(gVec) - posDiff.dot(hVec) * gMag2)/ norm;
55 
56  gPos = GlobalPoint ( gOrig + tG * gVec);
57  hPos = GlobalPoint ( hOrig + tH * hVec);
58 // cout << "TwoTrackMinimumDistanceLineLine "<<gPos<<hPos<<endl;
59 
60  pathG = tG * gMag;
61  pathH = tH * hMag;
62  return false;
63 }
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:68
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
T mag() const
Definition: PV3DBase.h:66
T z() const
Definition: PV3DBase.h:63
const MagneticField & magneticField() const
double TwoTrackMinimumDistanceLineLine::firstAngle ( ) const
inline

Definition at line 36 of file TwoTrackMinimumDistanceLineLine.h.

References phiG.

std::pair< double, double > TwoTrackMinimumDistanceLineLine::pathLength ( ) const

Definition at line 70 of file TwoTrackMinimumDistanceLineLine.cc.

References pathG, and pathH.

71 {
72  return std::pair<double, double> ( pathG, pathH);
73 }
std::pair< GlobalPoint, GlobalPoint > TwoTrackMinimumDistanceLineLine::points ( ) const

Definition at line 65 of file TwoTrackMinimumDistanceLineLine.cc.

References gPos, and hPos.

66 {
67  return std::pair<GlobalPoint, GlobalPoint > ( gPos, hPos );
68 }
double TwoTrackMinimumDistanceLineLine::secondAngle ( ) const
inline

Definition at line 37 of file TwoTrackMinimumDistanceLineLine.h.

References phiH.

Member Data Documentation

GlobalPoint TwoTrackMinimumDistanceLineLine::gPos
private

Definition at line 41 of file TwoTrackMinimumDistanceLineLine.h.

Referenced by calculate(), and points().

GlobalPoint TwoTrackMinimumDistanceLineLine::hPos
private

Definition at line 41 of file TwoTrackMinimumDistanceLineLine.h.

Referenced by calculate(), and points().

double TwoTrackMinimumDistanceLineLine::pathG
private

Definition at line 40 of file TwoTrackMinimumDistanceLineLine.h.

Referenced by calculate(), and pathLength().

double TwoTrackMinimumDistanceLineLine::pathH
private

Definition at line 40 of file TwoTrackMinimumDistanceLineLine.h.

Referenced by calculate(), and pathLength().

double TwoTrackMinimumDistanceLineLine::phiG
private

Definition at line 39 of file TwoTrackMinimumDistanceLineLine.h.

Referenced by calculate(), and firstAngle().

double TwoTrackMinimumDistanceLineLine::phiH
private

Definition at line 39 of file TwoTrackMinimumDistanceLineLine.h.

Referenced by calculate(), and secondAngle().