TrackingTools
PatternTools
src
TwoTrackMinimumDistanceHelixHelix.cc
Go to the documentation of this file.
1
#include "
TrackingTools/PatternTools/interface/TwoTrackMinimumDistanceHelixHelix.h
"
2
#include "
DataFormats/GeometryVector/interface/GlobalVector.h
"
3
#include "
MagneticField/Engine/interface/MagneticField.h
"
4
#include "
TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h
"
5
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
6
#include "
FWCore/Utilities/interface/isFinite.h
"
7
8
using namespace
std
;
9
10
TwoTrackMinimumDistanceHelixHelix::TwoTrackMinimumDistanceHelixHelix
()
11
: theH(nullptr), theG(nullptr), pointsUpdated(
false
), themaxjump(20), thesingjacI(1. / 0.1), themaxiter(4) {}
12
13
TwoTrackMinimumDistanceHelixHelix::~TwoTrackMinimumDistanceHelixHelix
() {}
14
15
bool
TwoTrackMinimumDistanceHelixHelix::updateCoeffs
(
const
GlobalPoint
& gpH,
const
GlobalPoint
& gpG) {
16
const
double
Bc2kH =
theH
->
magneticField
().
inInverseGeV
(gpH).
z
();
17
const
double
Bc2kG =
theG
->
magneticField
().
inInverseGeV
(gpG).
z
();
18
const
double
Ht =
theH
->
momentum
().
perp
();
19
const
double
Gt =
theG
->
momentum
().
perp
();
20
// thelambdaH=asin ( theH->momentum().z() / Hn );
21
22
if
(Ht == 0. || Gt == 0.) {
23
edm::LogWarning
(
"TwoTrackMinimumDistanceHelixHelix"
) <<
"transverse momentum of input trajectory is zero."
;
24
return
true
;
25
};
26
27
if
(
theH
->
charge
() == 0. ||
theG
->
charge
() == 0.) {
28
edm::LogWarning
(
"TwoTrackMinimumDistanceHelixHelix"
) <<
"charge of input track is zero."
;
29
return
true
;
30
};
31
32
if
(Bc2kG == 0.) {
33
edm::LogWarning
(
"TwoTrackMinimumDistanceHelixHelix"
) <<
"magnetic field at point "
<< gpG <<
" is zero."
;
34
return
true
;
35
};
36
37
if
(Bc2kH == 0.) {
38
edm::LogWarning
(
"TwoTrackMinimumDistanceHelixHelix"
) <<
"magnetic field at point "
<< gpH <<
" is zero."
;
39
return
true
;
40
};
41
42
theh
= Ht / (
theH
->
charge
() * Bc2kH);
43
thesinpH0
= -
theH
->
momentum
().
y
() / (Ht);
44
thecospH0
= -
theH
->
momentum
().
x
() / (Ht);
45
thetanlambdaH
= -
theH
->
momentum
().
z
() / (Ht);
46
thepH0
= atan2(
thesinpH0
,
thecospH0
);
47
48
// thelambdaG=asin ( theG->momentum().z()/( Gn ) );
49
50
theg
= Gt / (
theG
->
charge
() * Bc2kG);
51
thesinpG0
= -
theG
->
momentum
().
y
() / (Gt);
52
thecospG0
= -
theG
->
momentum
().
x
() / (Gt);
53
thetanlambdaG
= -
theG
->
momentum
().
z
() / (Gt);
54
thepG0
= atan2(
thesinpG0
,
thecospG0
);
55
56
thea
=
theH
->
position
().
x
() -
theG
->
position
().
x
() +
theg
*
thesinpG0
-
theh
*
thesinpH0
;
57
theb
=
theH
->
position
().
y
() -
theG
->
position
().
y
() -
theg
*
thecospG0
+
theh
*
thecospH0
;
58
thec1
=
theh
*
thetanlambdaH
*
thetanlambdaH
;
59
thec2
= -
theg
*
thetanlambdaG
*
thetanlambdaG
;
60
thed1
= -
theg
*
thetanlambdaG
*
thetanlambdaH
;
61
thed2
=
theh
*
thetanlambdaG
*
thetanlambdaH
;
62
thee1
=
thetanlambdaH
*
63
(
theH
->
position
().
z
() -
theG
->
position
().
z
() -
theh
*
thepH0
*
thetanlambdaH
+
theg
*
thepG0
*
thetanlambdaG
);
64
thee2
=
thetanlambdaG
*
65
(
theH
->
position
().
z
() -
theG
->
position
().
z
() -
theh
*
thepH0
*
thetanlambdaH
+
theg
*
thepG0
*
thetanlambdaG
);
66
return
false
;
67
}
68
69
bool
TwoTrackMinimumDistanceHelixHelix::oneIteration
(
double
& dH,
double
& dG) {
70
thesinpH
=
sin
(
thepH
);
71
thecospH
=
cos
(
thepH
);
72
thesinpG
=
sin
(
thepG
);
73
thecospG
=
cos
(
thepG
);
74
75
const
double
A11 =
theh
* (-
thesinpH
* (
thea
-
theg
*
thesinpG
) +
thecospH
* (
theb
+
theg
*
thecospG
) +
thec1
);
76
if
(A11 < 0) {
77
return
true
;
78
};
79
const
double
A22 = -
theg
* (-
thesinpG
* (
thea
+
theh
*
thesinpH
) +
thecospG
* (
theb
-
theh
*
thecospH
) +
thec2
);
80
if
(A22 < 0) {
81
return
true
;
82
};
83
const
double
A12 =
theh
* (-
theg
*
thecospG
*
thecospH
-
theg
*
thesinpH
*
thesinpG
+
thed1
);
84
const
double
A21 = -
theg
* (
theh
*
thecospG
*
thecospH
+
theh
*
thesinpH
*
thesinpG
+
thed2
);
85
const
double
detaI = 1. / (A11 * A22 - A12 * A21);
86
const
double
z1 =
theh
* (
thecospH
* (
thea
-
theg
*
thesinpG
) +
thesinpH
* (
theb
+
theg
*
thecospG
) +
thec1
*
thepH
+
87
thed1
*
thepG
+
thee1
);
88
const
double
z2
= -
theg
* (
thecospG
* (
thea
+
theh
*
thesinpH
) +
thesinpG
* (
theb
-
theh
*
thecospH
) +
thec2
*
thepG
+
89
thed2
*
thepH
+
thee2
);
90
91
dH = (z1 * A22 -
z2
* A12) * detaI;
92
dG = (
z2
* A11 - z1 * A21) * detaI;
93
if
(fabs(detaI) >
thesingjacI
) {
94
return
true
;
95
};
96
97
thepH
-= dH;
98
thepG
-= dG;
99
return
false
;
100
}
101
102
/*
103
bool TwoTrackMinimumDistanceHelixHelix::parallelTracks() const {
104
return (fabs(theH->momentum().x() - theG->momentum().x()) < .00000001 )
105
&& (fabs(theH->momentum().y() - theG->momentum().y()) < .00000001 )
106
&& (fabs(theH->momentum().z() - theG->momentum().z()) < .00000001 )
107
&& (theH->charge()==theG->charge())
108
;
109
}
110
*/
111
112
bool
TwoTrackMinimumDistanceHelixHelix::calculate
(
const
GlobalTrajectoryParameters
&
G
,
113
const
GlobalTrajectoryParameters
&
H
,
114
const
float
qual) {
115
pointsUpdated
=
false
;
116
theG
= &
G
;
117
theH
= &
H
;
118
bool
retval =
false
;
119
120
if
(
updateCoeffs
(
theG
->
position
(),
theH
->
position
())) {
121
finalPoints
();
122
return
true
;
123
}
124
125
thepG
=
thepG0
;
126
thepH
=
thepH0
;
127
128
int
counter
= 0;
129
double
pH = 0;
130
double
pG = 0;
131
do
{
132
retval =
oneIteration
(pG, pH);
133
if
(
edm::isNotFinite
(pG) ||
edm::isNotFinite
(pH))
134
retval =
true
;
135
if
(
counter
++ >
themaxiter
)
136
retval =
true
;
137
}
while
((!retval) && (fabs(pG) > qual || fabs(pH) > qual));
138
if
(fabs(
theg
* (
thepG
-
thepG0
)) >
themaxjump
)
139
retval =
true
;
140
if
(fabs(
theh
* (
thepH
-
thepH0
)) >
themaxjump
)
141
retval =
true
;
142
143
finalPoints
();
144
145
return
retval;
146
}
147
148
void
TwoTrackMinimumDistanceHelixHelix::finalPoints
() {
149
if
(
pointsUpdated
)
150
return
;
151
GlobalVector
tmpG(
sin
(
thepG
) -
thesinpG0
, -
cos
(
thepG
) +
thecospG0
,
thetanlambdaG
* (
thepG
-
thepG0
));
152
pointG
=
theG
->
position
() +
theg
* tmpG;
153
pathG
= (
thepG
-
thepG0
) * (
theg
*
sqrt
(1 +
thetanlambdaG
*
thetanlambdaG
));
154
155
GlobalVector
tmpH(
sin
(
thepH
) -
thesinpH0
, -
cos
(
thepH
) +
thecospH0
,
thetanlambdaH
* (
thepH
-
thepH0
));
156
pointH
=
theH
->
position
() +
theh
* tmpH;
157
pathH
= (
thepH
-
thepH0
) * (
theh
*
sqrt
(1 +
thetanlambdaH
*
thetanlambdaH
));
158
159
pointsUpdated
=
true
;
160
}
Vector3DBase
Definition:
Vector3DBase.h:8
TwoTrackMinimumDistanceHelixHelix::thesinpH0
double thesinpH0
Definition:
TwoTrackMinimumDistanceHelixHelix.h:48
class-composition.H
H
Definition:
class-composition.py:31
TwoTrackMinimumDistanceHelixHelix::thed1
double thed1
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
TwoTrackMinimumDistanceHelixHelix::thetanlambdaG
double thetanlambdaG
Definition:
TwoTrackMinimumDistanceHelixHelix.h:46
counter
Definition:
counter.py:1
TwoTrackMinimumDistanceHelixHelix::thepH
double thepH
Definition:
TwoTrackMinimumDistanceHelixHelix.h:53
GlobalTrajectoryParameters.h
MessageLogger.h
funct::false
false
Definition:
Factorize.h:29
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition:
isFinite.h:9
GlobalTrajectoryParameters::position
GlobalPoint position() const
Definition:
GlobalTrajectoryParameters.h:60
PV3DBase::x
T x() const
Definition:
PV3DBase.h:59
TwoTrackMinimumDistanceHelixHelix::thesinpH
double thesinpH
Definition:
TwoTrackMinimumDistanceHelixHelix.h:54
TwoTrackMinimumDistanceHelixHelix::theG
const GlobalTrajectoryParameters * theG
Definition:
TwoTrackMinimumDistanceHelixHelix.h:42
TwoTrackMinimumDistanceHelixHelix::pathH
double pathH
Definition:
TwoTrackMinimumDistanceHelixHelix.h:57
GlobalTrajectoryParameters::charge
TrackCharge charge() const
Definition:
GlobalTrajectoryParameters.h:72
callgraph.G
G
Definition:
callgraph.py:13
TwoTrackMinimumDistanceHelixHelix::thecospG
double thecospG
Definition:
TwoTrackMinimumDistanceHelixHelix.h:55
TwoTrackMinimumDistanceHelixHelix::thec2
double thec2
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
TwoTrackMinimumDistanceHelixHelix::finalPoints
void finalPoints()
Definition:
TwoTrackMinimumDistanceHelixHelix.cc:148
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition:
MessageLogger.h:122
TwoTrackMinimumDistanceHelixHelix::~TwoTrackMinimumDistanceHelixHelix
~TwoTrackMinimumDistanceHelixHelix()
Definition:
TwoTrackMinimumDistanceHelixHelix.cc:13
funct::sin
Sin< T >::type sin(const T &t)
Definition:
Sin.h:22
testProducerWithPsetDescEmpty_cfi.z2
z2
Definition:
testProducerWithPsetDescEmpty_cfi.py:41
TwoTrackMinimumDistanceHelixHelix::TwoTrackMinimumDistanceHelixHelix
TwoTrackMinimumDistanceHelixHelix()
Definition:
TwoTrackMinimumDistanceHelixHelix.cc:10
PV3DBase::z
T z() const
Definition:
PV3DBase.h:61
TwoTrackMinimumDistanceHelixHelix::pointG
GlobalPoint pointG
Definition:
TwoTrackMinimumDistanceHelixHelix.h:56
TwoTrackMinimumDistanceHelixHelix::pointsUpdated
bool pointsUpdated
Definition:
TwoTrackMinimumDistanceHelixHelix.h:58
funct::cos
Cos< T >::type cos(const T &t)
Definition:
Cos.h:22
TwoTrackMinimumDistanceHelixHelix::thecospH
double thecospH
Definition:
TwoTrackMinimumDistanceHelixHelix.h:55
TwoTrackMinimumDistanceHelixHelix::pointH
GlobalPoint pointH
Definition:
TwoTrackMinimumDistanceHelixHelix.h:56
TwoTrackMinimumDistanceHelixHelix::thecospG0
double thecospG0
Definition:
TwoTrackMinimumDistanceHelixHelix.h:47
mathSSE::sqrt
T sqrt(T t)
Definition:
SSEVec.h:19
TwoTrackMinimumDistanceHelixHelix::thee2
double thee2
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
GlobalTrajectoryParameters
Definition:
GlobalTrajectoryParameters.h:15
Point3DBase< float, GlobalTag >
GlobalTrajectoryParameters::momentum
GlobalVector momentum() const
Definition:
GlobalTrajectoryParameters.h:65
MagneticField::inInverseGeV
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition:
MagneticField.h:36
TwoTrackMinimumDistanceHelixHelix.h
TwoTrackMinimumDistanceHelixHelix::theb
double theb
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
TwoTrackMinimumDistanceHelixHelix::thec1
double thec1
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
TwoTrackMinimumDistanceHelixHelix::thed2
double thed2
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
PV3DBase::y
T y() const
Definition:
PV3DBase.h:60
TwoTrackMinimumDistanceHelixHelix::thesingjacI
double thesingjacI
Definition:
TwoTrackMinimumDistanceHelixHelix.h:60
MagneticField.h
TwoTrackMinimumDistanceHelixHelix::calculate
bool calculate(const GlobalTrajectoryParameters &, const GlobalTrajectoryParameters &, const float qual=.001)
Definition:
TwoTrackMinimumDistanceHelixHelix.cc:112
TwoTrackMinimumDistanceHelixHelix::thea
double thea
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
TwoTrackMinimumDistanceHelixHelix::pathG
double pathG
Definition:
TwoTrackMinimumDistanceHelixHelix.h:57
TwoTrackMinimumDistanceHelixHelix::thepH0
double thepH0
Definition:
TwoTrackMinimumDistanceHelixHelix.h:49
TwoTrackMinimumDistanceHelixHelix::thepG
double thepG
Definition:
TwoTrackMinimumDistanceHelixHelix.h:53
std
Definition:
JetResolutionObject.h:76
TwoTrackMinimumDistanceHelixHelix::oneIteration
bool oneIteration(double &, double &)
Definition:
TwoTrackMinimumDistanceHelixHelix.cc:69
TwoTrackMinimumDistanceHelixHelix::theg
double theg
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
isFinite.h
TwoTrackMinimumDistanceHelixHelix::thee1
double thee1
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
TwoTrackMinimumDistanceHelixHelix::theh
double theh
Definition:
TwoTrackMinimumDistanceHelixHelix.h:44
TwoTrackMinimumDistanceHelixHelix::themaxjump
double themaxjump
Definition:
TwoTrackMinimumDistanceHelixHelix.h:60
GlobalVector.h
GlobalTrajectoryParameters::magneticField
const MagneticField & magneticField() const
Definition:
GlobalTrajectoryParameters.h:106
TwoTrackMinimumDistanceHelixHelix::thetanlambdaH
double thetanlambdaH
Definition:
TwoTrackMinimumDistanceHelixHelix.h:46
TwoTrackMinimumDistanceHelixHelix::theH
const GlobalTrajectoryParameters * theH
Definition:
TwoTrackMinimumDistanceHelixHelix.h:42
TwoTrackMinimumDistanceHelixHelix::thepG0
double thepG0
Definition:
TwoTrackMinimumDistanceHelixHelix.h:49
TwoTrackMinimumDistanceHelixHelix::updateCoeffs
bool updateCoeffs(const GlobalPoint &, const GlobalPoint &)
Definition:
TwoTrackMinimumDistanceHelixHelix.cc:15
PV3DBase::perp
T perp() const
Definition:
PV3DBase.h:69
TwoTrackMinimumDistanceHelixHelix::thecospH0
double thecospH0
Definition:
TwoTrackMinimumDistanceHelixHelix.h:48
TwoTrackMinimumDistanceHelixHelix::thesinpG
double thesinpG
Definition:
TwoTrackMinimumDistanceHelixHelix.h:54
TwoTrackMinimumDistanceHelixHelix::thesinpG0
double thesinpG0
Definition:
TwoTrackMinimumDistanceHelixHelix.h:47
TwoTrackMinimumDistanceHelixHelix::themaxiter
int themaxiter
Definition:
TwoTrackMinimumDistanceHelixHelix.h:61
Generated for CMSSW Reference Manual by
1.8.16