1 |
#include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistanceLineLine.h"
|
2 |
#include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
|
3 |
#include "DataFormats/GeometryVector/interface/GlobalVector.h"
|
4 |
#include "FWCore/MessageLogger/interface/MessageLogger.h"
|
5 |
#include "MagneticField/Engine/interface/MagneticField.h"
|
6 |
|
7 |
|
8 |
bool TwoTrackMinimumDistanceLineLine::calculate(
|
9 |
const GlobalTrajectoryParameters & theG,
|
10 |
const GlobalTrajectoryParameters & theH)
|
11 |
{
|
12 |
GlobalPoint gOrig = theG.position();
|
13 |
GlobalPoint hOrig = theH.position();
|
14 |
if ( ( ( theH.charge() != 0. || theG.charge() != 0. ) ) &&
|
15 |
((theG.magneticField().inTesla(gOrig).z() != 0.)||
|
16 |
(theH.magneticField().inTesla(hOrig).z() != 0.)) )
|
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 |
}
|
64 |
|
65 |
std::pair<GlobalPoint, GlobalPoint> TwoTrackMinimumDistanceLineLine::points() const
|
66 |
{
|
67 |
return std::pair<GlobalPoint, GlobalPoint > ( gPos, hPos );
|
68 |
}
|
69 |
|
70 |
std::pair<double, double> TwoTrackMinimumDistanceLineLine::pathLength() const
|
71 |
{
|
72 |
return std::pair<double, double> ( pathG, pathH);
|
73 |
}
|