RecoTracker
TkSeedGenerator
interface
FastCircleFit.h
Go to the documentation of this file.
1
#ifndef RecoTracker_TkSeedGenerator_FastCircleFit_h
2
#define RecoTracker_TkSeedGenerator_FastCircleFit_h
3
4
#include "
DataFormats/GeometryVector/interface/GlobalPoint.h
"
5
#include "
DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h
"
6
#include "
CommonTools/Utils/interface/DynArray.h
"
7
8
#include <Eigen/Dense>
9
10
#ifdef MK_DEBUG
11
#include <iostream>
12
#define PRINT std::cout
13
#else
14
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
15
#define PRINT LogTrace("FastCircleFit")
16
#endif
17
18
#include <numeric>
19
#include <vector>
20
27
class
FastCircleFit
{
28
public
:
37
template
<
typename
P,
typename
E>
38
FastCircleFit
(
const
P
&
points
,
const
E&
errors
) {
39
const
auto
N
=
points
.size();
40
declareDynArray
(
float
,
N
,
x
);
41
declareDynArray
(
float
,
N
,
y
);
42
declareDynArray
(
float
,
N
,
z
);
43
declareDynArray
(
float
,
N
,
weight
);
// 1/sigma^2
44
calculate
(
points
,
errors
,
x
,
y
,
z
,
weight
);
45
}
46
50
template
<
size_t
N>
51
FastCircleFit
(
const
std::array<GlobalPoint, N>&
points
,
const
std::array<GlobalError, N>&
errors
) {
52
std::array<float, N>
x
;
53
std::array<float, N>
y
;
54
std::array<float, N>
z
;
55
std::array<float, N>
weight
;
// 1/sigma^2
56
calculate
(
points
,
errors
,
x
,
y
,
z
,
weight
);
57
}
58
59
~FastCircleFit
() =
default
;
60
61
float
x0
()
const
{
return
x0_
; }
62
float
y0
()
const
{
return
y0_
; }
63
float
rho
()
const
{
return
rho_
; }
64
65
// TODO: I'm not sure if the minimized square sum is chi2 distributed
66
float
chi2
()
const
{
return
chi2_
; }
67
68
private
:
69
template
<
typename
P,
typename
E,
typename
C>
70
void
calculate
(
const
P
&
points
,
const
E&
errors
,
C
&
x
,
C
&
y
,
C
&
z
,
C
&
weight
);
71
72
template
<
typename
T>
73
T
sqr
(
T
t
) {
74
return
t
*
t
;
75
}
76
77
float
x0_
;
78
float
y0_
;
79
float
rho_
;
80
float
chi2_
;
81
};
82
83
template
<
typename
P,
typename
E,
typename
C>
84
void
FastCircleFit::calculate
(
const
P
&
points
,
const
E&
errors
,
C
&
x
,
C
&
y
,
C
&
z
,
C
&
weight
) {
85
const
auto
N
=
points
.size();
86
87
// transform
88
for
(
size_t
i
= 0;
i
<
N
; ++
i
) {
89
const
auto
&
point
=
points
[
i
];
90
const
auto
&
p
=
point
.basicVector();
91
x
[
i
] =
p
.x();
92
y
[
i
] =
p
.y();
93
z
[
i
] =
sqr
(
p
.x()) +
sqr
(
p
.y());
94
95
// TODO: The following accounts only the hit uncertainty in phi.
96
// Albeit it "works nicely", should we account also the
97
// uncertainty in r in FPix (and the correlation between phi and r)?
98
const
float
phiErr2 =
errors
[
i
].phierr(
point
);
99
const
float
w
= 1.f / (
point
.perp2() * phiErr2);
100
weight
[
i
] =
w
;
101
102
PRINT
<<
" point "
<<
p
.x() <<
" "
<<
p
.y() <<
" transformed "
<<
x
[
i
] <<
" "
<<
y
[
i
] <<
" "
<<
z
[
i
] <<
" weight "
103
<<
w
<< std::endl;
104
}
105
const
float
invTotWeight = 1.f / std::accumulate(
weight
.begin(),
weight
.end(), 0.f);
106
PRINT
<<
" invTotWeight "
<< invTotWeight;
107
108
Eigen::Vector3f
mean
= Eigen::Vector3f::Zero();
109
for
(
size_t
i
= 0;
i
<
N
; ++
i
) {
110
const
float
w
=
weight
[
i
];
111
mean
[0] +=
w
*
x
[
i
];
112
mean
[1] +=
w
*
y
[
i
];
113
mean
[2] +=
w
*
z
[
i
];
114
}
115
mean
*= invTotWeight;
116
PRINT
<<
" mean "
<<
mean
[0] <<
" "
<<
mean
[1] <<
" "
<<
mean
[2] << std::endl;
117
118
Eigen::Matrix3f
A
= Eigen::Matrix3f::Zero();
119
for
(
size_t
i
= 0;
i
<
N
; ++
i
) {
120
const
auto
diff_x =
x
[
i
] -
mean
[0];
121
const
auto
diff_y =
y
[
i
] -
mean
[1];
122
const
auto
diff_z =
z
[
i
] -
mean
[2];
123
const
auto
w
=
weight
[
i
];
124
125
// exploit the fact that only lower triangular part of the matrix
126
// is used in the eigenvalue calculation
127
A
(0, 0) +=
w
* diff_x * diff_x;
128
A
(1, 0) +=
w
* diff_y * diff_x;
129
A
(1, 1) +=
w
* diff_y * diff_y;
130
A
(2, 0) +=
w
* diff_z * diff_x;
131
A
(2, 1) +=
w
* diff_z * diff_y;
132
A
(2, 2) +=
w
* diff_z * diff_z;
133
PRINT
<<
" i "
<<
i
<<
" A "
<<
A
(0, 0) <<
" "
<<
A
(0, 1) <<
" "
<<
A
(0, 2) << std::endl
134
<<
" "
<<
A
(1, 0) <<
" "
<<
A
(1, 1) <<
" "
<<
A
(1, 2) << std::endl
135
<<
" "
<<
A
(2, 0) <<
" "
<<
A
(2, 1) <<
" "
<<
A
(2, 2) << std::endl;
136
}
137
A
*= 1. / static_cast<float>(
N
);
138
139
PRINT
<<
" A "
<<
A
(0, 0) <<
" "
<<
A
(0, 1) <<
" "
<<
A
(0, 2) << std::endl
140
<<
" "
<<
A
(1, 0) <<
" "
<<
A
(1, 1) <<
" "
<<
A
(1, 2) << std::endl
141
<<
" "
<<
A
(2, 0) <<
" "
<<
A
(2, 1) <<
" "
<<
A
(2, 2) << std::endl;
142
143
// find eigenvalues and vectors
144
Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigen(
A
);
145
const
auto
& eigenValues = eigen.eigenvalues();
146
const
auto
& eigenVectors = eigen.eigenvectors();
147
148
// in Eigen the first eigenvalue is the smallest
149
PRINT
<<
" eigenvalues "
<< eigenValues[0] <<
" "
<< eigenValues[1] <<
" "
<< eigenValues[2] << std::endl;
150
151
PRINT
<<
" eigen "
<< eigenVectors(0, 0) <<
" "
<< eigenVectors(0, 1) <<
" "
<< eigenVectors(0, 2) << std::endl
152
<<
" "
<< eigenVectors(1, 0) <<
" "
<< eigenVectors(1, 1) <<
" "
<< eigenVectors(1, 2) << std::endl
153
<<
" "
<< eigenVectors(2, 0) <<
" "
<< eigenVectors(2, 1) <<
" "
<< eigenVectors(2, 2) << std::endl;
154
155
// eivenvector corresponding smallest eigenvalue
156
auto
n
= eigenVectors.col(0);
157
PRINT
<<
" n (unit) "
<<
n
[0] <<
" "
<<
n
[1] <<
" "
<<
n
[2] << std::endl;
158
159
const
float
c
= -1.f * (
n
[0] *
mean
[0] +
n
[1] *
mean
[1] +
n
[2] *
mean
[2]);
160
161
PRINT
<<
" c "
<<
c
<< std::endl;
162
163
// calculate fit parameters
164
const
auto
tmp
= 0.5f * 1.f /
n
[2];
165
x0_
= -
n
[0] *
tmp
;
166
y0_
= -
n
[1] *
tmp
;
167
const
float
rho2 = (1 -
sqr
(
n
[2]) - 4 *
c
*
n
[2]) *
sqr
(
tmp
);
168
rho_
=
std::sqrt
(rho2);
169
170
// calculate square sum
171
float
S
= 0;
172
for
(
size_t
i
= 0;
i
<
N
; ++
i
) {
173
const
float
incr =
sqr
(
c
+
n
[0] *
x
[
i
] +
n
[1] *
y
[
i
] +
n
[2] *
z
[
i
]) *
weight
[
i
];
174
#if defined(MK_DEBUG) || defined(EDM_ML_DEBUG)
175
const
float
distance
=
std::sqrt
(
sqr
(
x
[
i
] -
x0_
) +
sqr
(
y
[
i
] -
y0_
)) -
rho_
;
176
PRINT
<<
" i "
<<
i
<<
" chi2 incr "
<< incr <<
" d(hit, circle) "
<<
distance
<<
" sigma "
177
<< 1.f /
std::sqrt
(
weight
[
i
]) << std::endl;
178
#endif
179
S
+= incr;
180
}
181
chi2_
=
S
;
182
}
183
184
#endif
DDAxes::y
mps_fire.i
i
Definition:
mps_fire.py:428
FastCircleFit::x0
float x0() const
Definition:
FastCircleFit.h:61
MessageLogger.h
SiStripPI::mean
Definition:
SiStripPayloadInspectorHelper.h:169
dqmiodumpmetadata.n
n
Definition:
dqmiodumpmetadata.py:28
HLT_FULL_cff.points
points
Definition:
HLT_FULL_cff.py:21455
FastCircleFit::chi2_
float chi2_
Definition:
FastCircleFit.h:80
detailsBasic3DVector::z
float float float z
Definition:
extBasic3DVector.h:14
FastCircleFit::rho
float rho() const
Definition:
FastCircleFit.h:63
mps_merge.weight
weight
Definition:
mps_merge.py:88
AlCaHLTBitMon_ParallelJobs.p
p
Definition:
AlCaHLTBitMon_ParallelJobs.py:153
PRINT
#define PRINT
Definition:
FastCircleFit.h:15
DDAxes::x
createJobs.tmp
tmp
align.sh
Definition:
createJobs.py:716
FastCircleFit::sqr
T sqr(T t)
Definition:
FastCircleFit.h:73
FastCircleFit
Definition:
FastCircleFit.h:27
FastCircleFit::~FastCircleFit
~FastCircleFit()=default
errors
Definition:
errors.py:1
w
const double w
Definition:
UKUtility.cc:23
mathSSE::sqrt
T sqrt(T t)
Definition:
SSEVec.h:19
DDAxes::z
N
#define N
Definition:
blowfish.cc:9
declareDynArray
#define declareDynArray(T, n, x)
Definition:
DynArray.h:91
FastCircleFit::rho_
float rho_
Definition:
FastCircleFit.h:79
S
double S(const TLorentzVector &, const TLorentzVector &)
Definition:
Particle.cc:97
A
FastCircleFit::y0_
float y0_
Definition:
FastCircleFit.h:78
GlobalError.h
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:31
FastCircleFit::y0
float y0() const
Definition:
FastCircleFit.h:62
DynArray.h
FastCircleFit::chi2
float chi2() const
Definition:
FastCircleFit.h:66
MaterialEffects_cfi.A
A
Definition:
MaterialEffects_cfi.py:11
FastCircleFit::FastCircleFit
FastCircleFit(const P &points, const E &errors)
Definition:
FastCircleFit.h:38
gen::C
C
Definition:
PomwigHadronizer.cc:78
T
long double T
Definition:
Basic3DVectorLD.h:48
genVertex_cff.x
x
Definition:
genVertex_cff.py:12
detailsBasic3DVector::y
float float y
Definition:
extBasic3DVector.h:14
S
Definition:
CSCDBL1TPParametersExtended.h:16
FastCircleFit::x0_
float x0_
Definition:
FastCircleFit.h:77
P
std::pair< OmniClusterRef, TrackingParticleRef > P
Definition:
BDHadronTrackMonitoringAnalyzer.cc:202
point
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition:
invegas.h:5
HLT_FULL_cff.distance
distance
Definition:
HLT_FULL_cff.py:7733
submitPVValidationJobs.t
string t
Definition:
submitPVValidationJobs.py:644
GlobalPoint.h
FastCircleFit::calculate
void calculate(const P &points, const E &errors, C &x, C &y, C &z, C &weight)
Definition:
FastCircleFit.h:84
weight
Definition:
weight.py:1
FastCircleFit::FastCircleFit
FastCircleFit(const std::array< GlobalPoint, N > &points, const std::array< GlobalError, N > &errors)
Definition:
FastCircleFit.h:51
Generated for CMSSW Reference Manual by
1.8.16