ATLAS Offline Software
Loading...
Searching...
No Matches
AdaptiveMultiVertexFitterTestAlg.cxx
Go to the documentation of this file.
1/*
2 * Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration.
3 */
10
11
13#undef NDEBUG
15#include "xAODTracking/Vertex.h"
20#include "TrkTrack/Track.h"
24#include "GaudiKernel/SystemOfUnits.h"
25#include <cassert>
26#include <cmath>
27
28#include "CLHEP/Vector/LorentzVector.h"
29
30using Gaudi::Units::mm;
31using Gaudi::Units::MeV;
32using Gaudi::Units::GeV;
33
34
35namespace {
36
37
38AmgSymMatrix(5) cov5()
39{
40 AmgSymMatrix(5) m;
41 m.setIdentity();
42 return m;
43}
44
45
46AmgSymMatrix(5) cov5a()
47{
48 AmgSymMatrix(5) m;
49 m.setZero();
50 for (int i=0; i < 5; i++) {
51 (m)(i,i) = 1e-2;
52 }
53 (m)(1,1)=1;
54 return m;
55}
56
57
58using PerigeeUVec_t = std::vector<std::unique_ptr<Trk::Perigee> >;
59PerigeeUVec_t makePerigees1()
60{
61 Amg::Vector3D pos0 { 0, 0, 0 };
62
63 Amg::Vector3D pos1a { 2*mm, 1*mm, -10*mm };
64 Amg::Vector3D mom1a { 400*MeV, 600*MeV, 200*MeV };
65 Amg::Vector3D pos1b { 1*mm, 2*mm, -3*mm };
66 Amg::Vector3D mom1b { 600*MeV, 400*MeV, -200*MeV };
67 Amg::Vector3D pos1c { 1.2*mm, 1.3*mm, -7*mm };
68 Amg::Vector3D mom1c { 300*MeV, 1000*MeV, 100*MeV };
69
70 PerigeeUVec_t ret;
71
72 ret.emplace_back (std::make_unique<Trk::Perigee>(pos1a, mom1a, 1, pos1a, cov5()).release());
73 ret.emplace_back (std::make_unique<Trk::Perigee>(pos1b, mom1b, -1, pos1b, cov5()).release());
74 ret.emplace_back (std::make_unique<Trk::Perigee>(pos1c, mom1c, -1, pos1c, cov5()).release());
75
76 return ret;
77}
78
79
80using PerigeeUVec_t = std::vector<std::unique_ptr<Trk::Perigee> >;
81PerigeeUVec_t makePerigees2()
82{
83 Amg::Vector3D pos1a { 10*mm, 0*mm, -5*mm };
84 Amg::Vector3D mom1a { 1000*MeV, 0*MeV, 0*MeV };
85 Amg::Vector3D pos1b { 10.5*mm, 0.5*mm, -5.5*mm };
86 Amg::Vector3D mom1b { 800*MeV, 200*MeV, 200*MeV };
87 Amg::Vector3D pos1c { 9.5*mm, -0.5*mm, -4.5*mm };
88 Amg::Vector3D mom1c { 700*MeV, -300*MeV, -200*MeV };
89
90 PerigeeUVec_t ret;
91
92 ret.emplace_back (std::make_unique<Trk::Perigee>(pos1a, mom1a, 1, pos1a, cov5a()).release());
93 ret.emplace_back (std::make_unique<Trk::Perigee>(pos1b, mom1b, -1, pos1a, cov5a()).release());
94 ret.emplace_back (std::make_unique<Trk::Perigee>(pos1c, mom1c, -1, pos1a, cov5a()).release());
95
96 return ret;
97}
98
99
100using TrackUVec_t = std::vector<std::unique_ptr<Trk::Track> >;
101TrackUVec_t makeTracks (const PerigeeUVec_t& perigees)
102{
103 TrackUVec_t tracks;
104
105 for (const std::unique_ptr<Trk::Perigee>& p : perigees) {
107 auto fqual = std::make_unique<Trk::FitQuality> (0, 0);
108 auto tsos = std::make_unique<Trk::TrackStates>();
109 std::bitset<Trk::TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes> typePattern(0);
110 typePattern.set(Trk::TrackStateOnSurface::Perigee);
111 tsos->push_back(std::make_unique<Trk::TrackStateOnSurface>(
112 nullptr, p->uniqueClone(), nullptr, typePattern));
113 tracks.push_back(
114 std::make_unique<Trk::Track>(info, std::move(tsos), std::move(fqual)));
115 }
116
117 return tracks;
118}
119
120
121void dumpCovariance (const AmgSymMatrix(5)& m)
122{
123 for (int i=0; i < 5; i++) {
124 for (int j=0; j < 5; j++) {
125 std::cout << m(i,j) << ", ";
126 }
127 }
128}
129
130
131[[maybe_unused]]
132void dumpVertex (const xAOD::Vertex& v)
133{
134 std::cout << "vvv\n";
135 std::cout << v.x() << ", " << v.y() << ", " << v.z() << "\n";
136 static const SG::ConstAccessor<short> vertexTypeAcc ("vertexType");
137 if (vertexTypeAcc.isAvailable(v)) {
138 std::cout << "vertexType " << v.vertexType() << "\n";
139 }
140 std::cout << "chi2/ndof " << v.chiSquared() << ", " << v.numberDoF() << "\n";
141
142 std::cout << "cov ";
143 for (float f : v.covariance()) {
144 std::cout << f << ", ";
145 }
146 std::cout << "\n";
147
148 static const SG::ConstAccessor<std::vector<ElementLink<xAOD::TrackParticleContainer> > > trackParticleLinksAcc ("trackParticleLinks");
149 if (trackParticleLinksAcc.isAvailable(v)) {
150 std::cout << "tplinks ";
151 for (const ElementLink< xAOD::TrackParticleContainer >& l : v.trackParticleLinks()) {
152 std::cout << l.dataID() << "/" << l.index() << " ";
153 }
154 std::cout << "\n";
155 }
156
157 static const SG::ConstAccessor<float> trackWeightsAcc ("trackWeights");
158 if (trackWeightsAcc.isAvailable(v)) {
159 std::cout << "wt ";
160 for (float f : v.trackWeights()) {
161 std::cout << f << " ";
162 }
163 std::cout << "\n";
164 }
165
166 static const SG::ConstAccessor<std::vector<ElementLink<xAOD::NeutralParticleContainer> > > neutralParticleLinksAcc ("neutralParticleLinks");
167 if (neutralParticleLinksAcc.isAvailable(v)) {
168 std::cout << "nplinks ";
169 for (const ElementLink< xAOD::NeutralParticleContainer >& l : v.neutralParticleLinks()) {
170 std::cout << l.dataID() << "/" << l.index() << " ";
171 }
172 std::cout << "\n";
173 }
174
175 static const SG::ConstAccessor<float> neutralWeightsAcc ("neutralWeights");
176 if (neutralWeightsAcc.isAvailable(v)) {
177 std::cout << "wt ";
178 for (float f : v.neutralWeights()) {
179 std::cout << f << " ";
180 }
181 std::cout << "\n";
182 }
183
184 std::cout << v.vxTrackAtVertexAvailable() << "\n";
185 for (const Trk::VxTrackAtVertex& vv : v.vxTrackAtVertex()) {
186 vv.dump (std::cout);
187 std::cout << "cov ";
188 if (vv.perigeeAtVertex()) {
189 dumpCovariance (*vv.perigeeAtVertex()->covariance());
190 }
191 else {
192 std::cout << "(null)";
193 }
194 std::cout << "\n";
195 }
196}
197
198
199void assertVec3D (const char* which,
200 const Amg::Vector3D& a,
201 const Amg::Vector3D& b,
202 double thresh = 2e-5)
203{
204 if ( ! Athena_test::isEqual (a.x(), b.x(), thresh) ||
205 ! Athena_test::isEqual (a.y(), b.y(), thresh) ||
206 ! Athena_test::isEqual (a.z(), b.z(), thresh) )
207 {
208 std::cerr << "TrkVKalVrtFitterTestAlg::assertVec3D mismatch " << which
209 << ": ["
210 << a.x() << ", "
211 << a.y() << ", "
212 << a.z() << "] / ["
213 << b.x() << ", "
214 << b.y() << ", "
215 << b.z() << "]\n";
216 std::abort();
217 }
218}
219
220
221void compareVertex (const xAOD::Vertex& a, const xAOD::Vertex& b)
222{
223 assertVec3D ("vertex pos", a.position(), b.position(), 5e-4);
224 assert (Athena_test::isEqual (a.chiSquared(), b.chiSquared(), 5e-3) );
225 assert (Athena_test::isEqual (a.numberDoF(), b.numberDoF(), 5e-3) );
226
227 assert (a.covariance().size() == b.covariance().size());
228 for (unsigned int i = 0; i < a.covariance().size(); i++) {
229 if (std::isinf(a.covariance()[i]) && std::isinf(b.covariance()[i])) continue;
230 assert (Athena_test::isEqual (a.covariance()[i], b.covariance()[i], 2e-2) );
231 }
232}
233
234
235void setInitialPerigee (xAOD::Vertex& v, unsigned i, const Trk::Perigee* p)
236{
237 std::vector< Trk::VxTrackAtVertex >& vec = v.vxTrackAtVertex();
238 if (vec.size() <= i) vec.resize(i+1);
239 vec[i].setInitialPerigee (p);
240}
241
242
243void setInitialPerigees (xAOD::Vertex& v, const TrackUVec_t& tracks)
244{
245 int i = 0;
246 for (const std::unique_ptr<Trk::Track>& t : tracks) {
247 setInitialPerigee (v, i, t->perigeeParameters());
248 ++i;
249 }
250}
251
252
253struct VertexInfo
254{
255 xAOD::Vertex v;
256 TrackUVec_t tracks;
257 PerigeeUVec_t perigees;
258 std::unique_ptr<Trk::MvfFitInfo> fi;
259 std::vector<Trk::VxTrackAtVertex*> vtracks;
260 std::vector<Trk::TrackToVtxLink> links;
261};
262
263
264void initVertex (VertexInfo& vi,
265 const Amg::Vector3D& pos,
266 PerigeeUVec_t&& perigees)
267{
268 xAOD::Vertex::Decorator< Trk::MvfFitInfo* > MvfFitInfo("MvfFitInfo");
269 xAOD::Vertex::Decorator< bool > isInitialized("isInitialized");
271
272 vi.perigees = std::move (perigees);
273
274 vi.v.makePrivateStore();
275 vi.v.setPosition (pos);
276
277 vi.tracks = makeTracks (vi.perigees);
278 setInitialPerigees (vi.v, vi.tracks);
279
280 auto cv = std::make_unique<xAOD::Vertex>();
281 cv->makePrivateStore();
282 cv->setPosition (pos);
283 Amg::MatrixX looseConstraintCovariance(3, 3);
284 looseConstraintCovariance.setIdentity();
285 looseConstraintCovariance = looseConstraintCovariance * 1e+5;
286 cv->setCovariancePosition (looseConstraintCovariance);
287 cv->setFitQuality (0, -3);
288 vi.fi = std::make_unique<Trk::MvfFitInfo> (cv.release(),
289 new Amg::Vector3D(pos),
290 new Amg::Vector3D(pos));
291 MvfFitInfo(vi.v) = vi.fi.get();
292 isInitialized(vi.v) = false;
293
294 vi.links.reserve (vi.perigees.size());
295 for (std::unique_ptr<Trk::Perigee>& p : vi.perigees) {
296 vi.links.push_back (new std::vector<xAOD::Vertex*> {&vi.v});
297 auto tav = std::make_unique<Trk::MVFVxTrackAtVertex>(1.5, p->clone(), p.get());
298 tav->setLinkToVertices (&vi.links.back());
299 vi.vtracks.push_back (tav.release());
300 }
301 VTAV(vi.v) = vi.vtracks;
302}
303
304
305} // anonymous namespace
306
307
308namespace Trk {
309
310
315{
316 ATH_CHECK( m_fitter.retrieve() );
317 return StatusCode::SUCCESS;
318}
319
320
325{
326 ATH_MSG_VERBOSE ("execute");
327
328 ATH_CHECK( test1() );
329
330 return StatusCode::SUCCESS;
331}
332
333
335{
336 VertexInfo v1;
337 initVertex (v1, {1.5*mm, 1.7*mm, -6*mm}, makePerigees1());
338
339
340 VertexInfo v2;
341 initVertex (v2, {9.8*mm, 0.2*mm, -4.8*mm}, makePerigees2());
342
343 std::vector<xAOD::Vertex*> verts {&v1.v, &v2.v};
344 m_fitter->fit (verts);
345
346 xAOD::Vertex exp_v1;
347 exp_v1.makePrivateStore();
348 exp_v1.setPosition ({8.60555, 10.276, -6.33474});
349 exp_v1.setFitQuality (0.0155243, 0.811922);
350 exp_v1.setCovariance (std::vector<float>
351 {730.916, 138.738, 926.494,
352 -7.25401, 46.1613, 66.2963});
353
354 xAOD::Vertex exp_v2;
355 exp_v2.makePrivateStore();
356 exp_v2.setPosition ({7.97405, 0.106089, -4.97332});
357 exp_v2.setFitQuality (1.74032, 1.6738);
358 exp_v2.setCovariance (std::vector<float>
359 {0.34849, -0.0194732, 0.0290884,
360 -0.00362668, 0.000706326, 0.445416});
361
362 compareVertex (v1.v, exp_v1);
363 compareVertex (v2.v, exp_v2);
364
365 return StatusCode::SUCCESS;
366}
367
368
369} // namespace Trk
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_VERBOSE(x)
std::vector< size_t > vec
#define AmgSymMatrix(dim)
static Double_t a
void makePrivateStore()
Create a new (empty) private store for this object.
Helper class to provide constant type-safe access to aux data.
SG::Decorator< T, ALLOC > Decorator
class to provide type-safe access to aux data.
Definition AuxElement.h:135
ToolHandle< Trk::AdaptiveMultiVertexFitter > m_fitter
virtual StatusCode execute() override
Execute the algorithm.
virtual StatusCode initialize() override
Standard Gaudi initialize method.
Contains information about the 'fitter' of this track.
@ Perigee
This represents a perigee, and so will contain a Perigee object only.
The VxTrackAtVertex is a common class for all present TrkVertexFitters The VxTrackAtVertex is designe...
void setCovariance(const std::vector< float > &value)
Sets the covariance matrix as a simple vector of values.
void setPosition(const Amg::Vector3D &position)
Sets the 3-position.
void setFitQuality(float chiSquared, float numberDoF)
Set the 'Fit Quality' information.
static std::string release
Definition computils.h:50
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Matrix< double, 3, 1 > Vector3D
bool isEqual(double x1, double x2, double thresh=1e-6)
Definition FLOATassert.h:26
l
Printing final latex table to .tex output file.
Ensure that the ATLAS eigen extensions are properly loaded.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
which(filename, env=os.environ)
UNIX-style which ---------------------------------------------------------—.
Definition unixtools.py:39
Vertex_v1 Vertex
Define the latest version of the vertex class.