File indexing completed on 2025-11-08 09:20:37
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include <boost/test/unit_test.hpp>
0010
0011 #include "Acts/Definitions/Algebra.hpp"
0012 #include "Acts/Definitions/TrackParametrization.hpp"
0013 #include "Acts/Definitions/Units.hpp"
0014 #include "Acts/EventData/GenericBoundTrackParameters.hpp"
0015 #include "Acts/EventData/TrackParameters.hpp"
0016 #include "Acts/Geometry/GeometryContext.hpp"
0017 #include "Acts/Geometry/GeometryIdentifier.hpp"
0018 #include "Acts/MagneticField/ConstantBField.hpp"
0019 #include "Acts/MagneticField/MagneticFieldContext.hpp"
0020 #include "Acts/Propagator/EigenStepper.hpp"
0021 #include "Acts/Propagator/Propagator.hpp"
0022 #include "Acts/Surfaces/PerigeeSurface.hpp"
0023 #include "Acts/Surfaces/Surface.hpp"
0024 #include "Acts/Utilities/AnnealingUtility.hpp"
0025 #include "Acts/Utilities/Logger.hpp"
0026 #include "Acts/Utilities/Result.hpp"
0027 #include "Acts/Vertexing/AMVFInfo.hpp"
0028 #include "Acts/Vertexing/AdaptiveMultiVertexFitter.hpp"
0029 #include "Acts/Vertexing/HelicalTrackLinearizer.hpp"
0030 #include "Acts/Vertexing/ImpactPointEstimator.hpp"
0031 #include "Acts/Vertexing/TrackAtVertex.hpp"
0032 #include "Acts/Vertexing/Vertex.hpp"
0033 #include "Acts/Vertexing/VertexingOptions.hpp"
0034 #include "ActsTests/CommonHelpers/FloatComparisons.hpp"
0035
0036 #include <iostream>
0037 #include <map>
0038 #include <memory>
0039 #include <numbers>
0040 #include <random>
0041 #include <utility>
0042 #include <vector>
0043
0044 using namespace Acts;
0045 using namespace Acts::UnitLiterals;
0046
0047 namespace ActsTests {
0048
0049 using Acts::VectorHelpers::makeVector4;
0050
0051
0052 ACTS_LOCAL_LOGGER(getDefaultLogger("AMVFitterTests", Logging::INFO))
0053
0054 using Covariance = BoundSquareMatrix;
0055 using Propagator = Acts::Propagator<EigenStepper<>>;
0056 using Linearizer = HelicalTrackLinearizer;
0057
0058
0059 GeometryContext geoContext = GeometryContext();
0060 MagneticFieldContext magFieldContext = MagneticFieldContext();
0061
0062
0063 std::uniform_real_distribution<double> vXYDist(-0.1_mm, 0.1_mm);
0064
0065 std::uniform_real_distribution<double> vZDist(-20_mm, 20_mm);
0066
0067 std::uniform_real_distribution<double> d0Dist(-0.01_mm, 0.01_mm);
0068
0069 std::uniform_real_distribution<double> z0Dist(-0.2_mm, 0.2_mm);
0070
0071 std::uniform_real_distribution<double> pTDist(1._GeV, 30._GeV);
0072
0073 std::uniform_real_distribution<double> phiDist(-std::numbers::pi,
0074 std::numbers::pi);
0075
0076 std::uniform_real_distribution<double> thetaDist(1., std::numbers::pi - 1.);
0077
0078 std::uniform_real_distribution<double> qDist(-1, 1);
0079
0080
0081 std::uniform_real_distribution<double> relTDist(-4_ps, 4_ps);
0082
0083 std::uniform_real_distribution<double> resIPDist(0., 100._um);
0084
0085 std::uniform_real_distribution<double> resAngDist(0., 0.1);
0086
0087 std::uniform_real_distribution<double> resQoPDist(-0.1, 0.1);
0088
0089
0090 std::uniform_real_distribution<double> resTDist(0_ps, 8_ps);
0091
0092 BOOST_AUTO_TEST_SUITE(VertexingSuite)
0093
0094
0095
0096 BOOST_AUTO_TEST_CASE(adaptive_multi_vertex_fitter_test) {
0097
0098 int mySeed = 31415;
0099 std::mt19937 gen(mySeed);
0100
0101
0102 auto bField = std::make_shared<ConstantBField>(Vector3{0.0, 0.0, 1_T});
0103
0104
0105 EigenStepper<> stepper(bField);
0106
0107
0108 auto propagator = std::make_shared<Propagator>(stepper);
0109
0110 VertexingOptions vertexingOptions(geoContext, magFieldContext);
0111
0112
0113 ImpactPointEstimator::Config ip3dEstCfg(bField, propagator);
0114 ImpactPointEstimator ip3dEst(ip3dEstCfg);
0115
0116
0117 Linearizer::Config ltConfig;
0118 ltConfig.bField = bField;
0119 ltConfig.propagator = propagator;
0120 Linearizer linearizer(ltConfig);
0121
0122 AdaptiveMultiVertexFitter::Config fitterCfg(ip3dEst);
0123 fitterCfg.trackLinearizer.connect<&Linearizer::linearizeTrack>(&linearizer);
0124
0125
0126 fitterCfg.doSmoothing = true;
0127 fitterCfg.extractParameters.connect<&InputTrack::extractParameters>();
0128
0129 AdaptiveMultiVertexFitter fitter(std::move(fitterCfg));
0130
0131
0132
0133 Vector3 vtxPos1(-0.15_mm, -0.1_mm, -1.5_mm);
0134 Vector3 vtxPos2(-0.1_mm, -0.15_mm, -3._mm);
0135 Vector3 vtxPos3(0.2_mm, 0.2_mm, 10._mm);
0136
0137 std::vector<Vector3> vtxPosVec{vtxPos1, vtxPos2, vtxPos3};
0138
0139
0140 double resD0 = resIPDist(gen);
0141 double resZ0 = resIPDist(gen);
0142 double resPh = resAngDist(gen);
0143 double resTh = resAngDist(gen);
0144 double resQp = resQoPDist(gen);
0145
0146 std::vector<Vertex> vtxList;
0147 for (auto& vtxPos : vtxPosVec) {
0148 Vertex vtx(vtxPos);
0149
0150 SquareMatrix4 posCovariance(SquareMatrix4::Identity());
0151 vtx.setFullCovariance(posCovariance);
0152
0153 vtxList.push_back(vtx);
0154 }
0155
0156 std::vector<Vertex*> vtxPtrList;
0157 ACTS_DEBUG("All vertices in test case:");
0158 int cv = 0;
0159 for (auto& vtx : vtxList) {
0160 cv++;
0161 ACTS_DEBUG("\t" << cv << ". vertex ptr: " << &vtx);
0162 vtxPtrList.push_back(&vtx);
0163 }
0164
0165 std::vector<BoundTrackParameters> allTracks;
0166
0167 unsigned int nTracksPerVtx = 4;
0168
0169
0170 for (unsigned int iTrack = 0; iTrack < nTracksPerVtx * vtxPosVec.size();
0171 iTrack++) {
0172
0173 double q = qDist(gen) < 0 ? -1. : 1.;
0174
0175
0176 Covariance covMat;
0177
0178 covMat << resD0 * resD0, 0., 0., 0., 0., 0., 0., resZ0 * resZ0, 0., 0., 0.,
0179 0., 0., 0., resPh * resPh, 0., 0., 0., 0., 0., 0., resTh * resTh, 0.,
0180 0., 0., 0., 0., 0., resQp * resQp, 0., 0., 0., 0., 0., 0., 1.;
0181
0182
0183 int vtxIdx = static_cast<int>(iTrack / nTracksPerVtx);
0184
0185
0186 BoundTrackParameters::ParametersVector paramVec;
0187 paramVec << d0Dist(gen), z0Dist(gen), phiDist(gen), thetaDist(gen),
0188 q / pTDist(gen), 0.;
0189
0190 std::shared_ptr<PerigeeSurface> perigeeSurface =
0191 Surface::makeShared<PerigeeSurface>(vtxPosVec[vtxIdx]);
0192
0193 allTracks.emplace_back(perigeeSurface, paramVec, std::move(covMat),
0194 ParticleHypothesis::pion());
0195 }
0196
0197 int ct = 0;
0198 ACTS_DEBUG("All tracks in test case:");
0199 for (auto& trk : allTracks) {
0200 ct++;
0201 ACTS_DEBUG("\t" << ct << ". track ptr: " << &trk);
0202 }
0203
0204 AdaptiveMultiVertexFitter::State state(*bField, magFieldContext);
0205
0206 for (unsigned int iTrack = 0; iTrack < nTracksPerVtx * vtxPosVec.size();
0207 iTrack++) {
0208
0209 int vtxIdx = static_cast<int>(iTrack / nTracksPerVtx);
0210
0211 InputTrack inputTrack{&allTracks[iTrack]};
0212
0213 state.vtxInfoMap[&(vtxList[vtxIdx])].trackLinks.push_back(inputTrack);
0214 state.tracksAtVerticesMap.insert(
0215 std::make_pair(std::make_pair(inputTrack, &(vtxList[vtxIdx])),
0216 TrackAtVertex(1., allTracks[iTrack], inputTrack)));
0217
0218
0219
0220 if (iTrack == 0) {
0221 state.vtxInfoMap[&(vtxList.at(1))].trackLinks.push_back(inputTrack);
0222 state.tracksAtVerticesMap.insert(
0223 std::make_pair(std::make_pair(inputTrack, &(vtxList.at(1))),
0224 TrackAtVertex(1., allTracks[iTrack], inputTrack)));
0225 }
0226 }
0227
0228 for (auto& vtx : vtxPtrList) {
0229 state.addVertexToMultiMap(*vtx);
0230 ACTS_DEBUG("Vertex, with ptr: " << vtx);
0231 for (auto& trk : state.vtxInfoMap[vtx].trackLinks) {
0232 ACTS_DEBUG("\t track ptr: " << trk);
0233 }
0234 }
0235
0236 ACTS_DEBUG("Checking all vertices linked to a single track:");
0237 for (auto& trk : allTracks) {
0238 ACTS_DEBUG("Track with ptr: " << &trk);
0239 auto range = state.trackToVerticesMultiMap.equal_range(InputTrack{&trk});
0240 for (auto vtxIter = range.first; vtxIter != range.second; ++vtxIter) {
0241 ACTS_DEBUG("\t used by vertex: " << vtxIter->second);
0242 }
0243 }
0244
0245
0246
0247 std::vector<Vertex> seedListCopy = vtxList;
0248
0249 std::vector<Vertex*> vtxFitPtr = {&vtxList.at(0)};
0250 auto res1 = fitter.addVtxToFit(state, vtxFitPtr, vertexingOptions);
0251 ACTS_DEBUG("Tracks linked to each vertex AFTER fit:");
0252 int c = 0;
0253 for (auto& vtx : vtxPtrList) {
0254 c++;
0255 ACTS_DEBUG(c << ". vertex, with ptr: " << vtx);
0256 for (const auto& trk : state.vtxInfoMap[vtx].trackLinks) {
0257 ACTS_DEBUG("\t track ptr: " << trk);
0258 }
0259 }
0260
0261 ACTS_DEBUG("Checking all vertices linked to a single track AFTER fit:");
0262 for (auto& trk : allTracks) {
0263 ACTS_DEBUG("Track with ptr: " << &trk);
0264 auto range = state.trackToVerticesMultiMap.equal_range(InputTrack{&trk});
0265 for (auto vtxIter = range.first; vtxIter != range.second; ++vtxIter) {
0266 ACTS_DEBUG("\t used by vertex: " << vtxIter->second);
0267 }
0268 }
0269
0270 BOOST_CHECK(res1.ok());
0271
0272 ACTS_DEBUG("Vertex positions after fit of vertex 1 and 2:");
0273 for (std::size_t vtxIter = 0; vtxIter < 3; vtxIter++) {
0274 ACTS_DEBUG("Vtx " << vtxIter + 1 << ", seed position:\n "
0275 << seedListCopy.at(vtxIter).fullPosition()
0276 << "\nFitted position:\n "
0277 << vtxList.at(vtxIter).fullPosition());
0278 }
0279
0280
0281
0282 BOOST_CHECK_NE(vtxList.at(0).fullPosition(),
0283 seedListCopy.at(0).fullPosition());
0284 BOOST_CHECK_NE(vtxList.at(1).fullPosition(),
0285 seedListCopy.at(1).fullPosition());
0286 BOOST_CHECK_EQUAL(vtxList.at(2).fullPosition(),
0287 seedListCopy.at(2).fullPosition());
0288
0289 CHECK_CLOSE_ABS(vtxList.at(0).fullPosition(),
0290 seedListCopy.at(0).fullPosition(), 1_mm);
0291 CHECK_CLOSE_ABS(vtxList.at(1).fullPosition(),
0292 seedListCopy.at(1).fullPosition(), 1_mm);
0293
0294 vtxFitPtr = {&vtxList.at(2)};
0295 auto res2 = fitter.addVtxToFit(state, vtxFitPtr, vertexingOptions);
0296 BOOST_CHECK(res2.ok());
0297
0298
0299 BOOST_CHECK_NE(vtxList.at(2).fullPosition(),
0300 seedListCopy.at(2).fullPosition());
0301 CHECK_CLOSE_ABS(vtxList.at(2).fullPosition(),
0302 seedListCopy.at(2).fullPosition(), 1_mm);
0303
0304 ACTS_DEBUG("Vertex positions after fit of vertex 3:");
0305 ACTS_DEBUG("Vtx 1, seed position:\n " << seedListCopy.at(0).fullPosition()
0306 << "\nFitted position:\n "
0307 << vtxList.at(0).fullPosition());
0308 ACTS_DEBUG("Vtx 2, seed position:\n " << seedListCopy.at(1).fullPosition()
0309 << "\nFitted position:\n "
0310 << vtxList.at(1).fullPosition());
0311 ACTS_DEBUG("Vtx 3, seed position:\n " << seedListCopy.at(2).fullPosition()
0312 << "\nFitted position:\n "
0313 << vtxList.at(2).fullPosition());
0314 }
0315
0316
0317
0318 BOOST_AUTO_TEST_CASE(time_fitting) {
0319
0320 int mySeed = 31415;
0321 std::mt19937 gen(mySeed);
0322
0323
0324 auto bField = std::make_shared<ConstantBField>(Vector3{0.0, 0.0, 1_T});
0325
0326
0327 EigenStepper<> stepper(bField);
0328
0329
0330 auto propagator = std::make_shared<Propagator>(stepper);
0331
0332 VertexingOptions vertexingOptions(geoContext, magFieldContext);
0333
0334 ImpactPointEstimator::Config ip3dEstCfg(bField, propagator);
0335 ImpactPointEstimator ip3dEst(ip3dEstCfg);
0336
0337 AdaptiveMultiVertexFitter::Config fitterCfg(ip3dEst);
0338
0339
0340 Linearizer::Config ltConfig;
0341 ltConfig.bField = bField;
0342 ltConfig.propagator = propagator;
0343 Linearizer linearizer(ltConfig);
0344
0345
0346 fitterCfg.doSmoothing = true;
0347
0348 fitterCfg.useTime = true;
0349 fitterCfg.extractParameters.connect<&InputTrack::extractParameters>();
0350 fitterCfg.trackLinearizer.connect<&Linearizer::linearizeTrack>(&linearizer);
0351
0352 AdaptiveMultiVertexFitter fitter(std::move(fitterCfg));
0353
0354
0355 double trueVtxTime = 40.0_ps;
0356 Vector3 trueVtxPos(-0.15_mm, -0.1_mm, -1.5_mm);
0357
0358
0359 Vector4 vtxSeedPos(0.0_mm, 0.0_mm, -1.4_mm, 0.0_ps);
0360
0361 Vertex vtx(vtxSeedPos);
0362
0363 SquareMatrix4 initialCovariance(SquareMatrix4::Identity() * 1e+8);
0364 vtx.setFullCovariance(initialCovariance);
0365
0366
0367 std::vector<BoundTrackParameters> trks;
0368
0369 unsigned int nTracks = 4;
0370 for (unsigned int _ = 0; _ < nTracks; _++) {
0371
0372 double q = qDist(gen) < 0 ? -1. : 1.;
0373
0374
0375 double resD0 = resIPDist(gen);
0376 double resZ0 = resIPDist(gen);
0377 double resPh = resAngDist(gen);
0378 double resTh = resAngDist(gen);
0379 double resQp = resQoPDist(gen);
0380 double resT = resTDist(gen);
0381
0382
0383 Covariance covMat;
0384
0385
0386 covMat <<
0387 resD0 * resD0, 0., 0., 0., 0., 0.,
0388 0., resZ0 * resZ0, 0., 0., 0., 0.,
0389 0., 0., resPh * resPh, 0., 0., 0.,
0390 0., 0., 0., resTh * resTh, 0., 0.,
0391 0., 0., 0., 0., resQp * resQp, 0.,
0392 0., 0., 0., 0., 0., resT * resT;
0393
0394
0395
0396 BoundTrackParameters::ParametersVector paramVec;
0397 paramVec << d0Dist(gen), z0Dist(gen), phiDist(gen), thetaDist(gen),
0398 q / pTDist(gen), trueVtxTime + relTDist(gen);
0399
0400 std::shared_ptr<PerigeeSurface> perigeeSurface =
0401 Surface::makeShared<PerigeeSurface>(trueVtxPos);
0402
0403 trks.emplace_back(perigeeSurface, paramVec, std::move(covMat),
0404 ParticleHypothesis::pion());
0405 }
0406
0407 std::vector<const BoundTrackParameters*> trksPtr;
0408 for (const auto& trk : trks) {
0409 trksPtr.push_back(&trk);
0410 }
0411
0412
0413 AdaptiveMultiVertexFitter::State state(*bField, magFieldContext);
0414
0415 for (const auto& trk : trks) {
0416 ACTS_DEBUG("Track parameters:\n" << trk);
0417
0418 state.vtxInfoMap[&vtx].trackLinks.push_back(InputTrack{&trk});
0419 state.tracksAtVerticesMap.insert(
0420 std::make_pair(std::make_pair(InputTrack{&trk}, &vtx),
0421 TrackAtVertex(1., trk, InputTrack{&trk})));
0422 }
0423
0424 state.addVertexToMultiMap(vtx);
0425
0426 std::vector<Vertex*> vtxFitPtr = {&vtx};
0427 auto res = fitter.addVtxToFit(state, vtxFitPtr, vertexingOptions);
0428
0429 BOOST_CHECK(res.ok());
0430
0431 ACTS_DEBUG("Truth vertex position: " << trueVtxPos.transpose());
0432 ACTS_DEBUG("Fitted vertex position: " << vtx.position().transpose());
0433
0434 ACTS_DEBUG("Truth vertex time: " << trueVtxTime);
0435 ACTS_DEBUG("Fitted vertex time: " << vtx.time());
0436
0437
0438 CHECK_CLOSE_ABS(trueVtxPos, vtx.position(), 60_um);
0439 CHECK_CLOSE_ABS(trueVtxTime, vtx.time(), 2_ps);
0440
0441 const SquareMatrix4& vtxCov = vtx.fullCovariance();
0442
0443 ACTS_DEBUG("Vertex 4D covariance after the fit:\n" << vtxCov);
0444
0445
0446 for (std::size_t i = 0; i <= 3; i++) {
0447 BOOST_CHECK_GT(vtxCov(i, i), 0.);
0448 }
0449
0450
0451 CHECK_CLOSE_ABS(vtxCov - vtxCov.transpose(), SquareMatrix4::Zero(), 1e-5);
0452 }
0453
0454
0455
0456
0457 BOOST_AUTO_TEST_CASE(adaptive_multi_vertex_fitter_test_athena) {
0458
0459 auto bField = std::make_shared<ConstantBField>(Vector3{0.0, 0.0, 2_T});
0460
0461
0462
0463 EigenStepper<> stepper(bField);
0464
0465
0466 auto propagator = std::make_shared<Propagator>(stepper);
0467
0468 VertexingOptions vertexingOptions(geoContext, magFieldContext);
0469
0470 ImpactPointEstimator::Config ip3dEstCfg(bField, propagator);
0471 ImpactPointEstimator ip3dEst(ip3dEstCfg);
0472
0473 std::vector<double> temperatures(1, 3.);
0474 AnnealingUtility::Config annealingConfig;
0475 annealingConfig.setOfTemperatures = temperatures;
0476 AnnealingUtility annealingUtility(annealingConfig);
0477
0478 AdaptiveMultiVertexFitter::Config fitterCfg(ip3dEst);
0479
0480 fitterCfg.annealingTool = annealingUtility;
0481 fitterCfg.extractParameters.connect<&InputTrack::extractParameters>();
0482
0483
0484 Linearizer::Config ltConfig;
0485 ltConfig.bField = bField;
0486 ltConfig.propagator = propagator;
0487 Linearizer linearizer(ltConfig);
0488
0489 fitterCfg.trackLinearizer.connect<&Linearizer::linearizeTrack>(&linearizer);
0490
0491
0492
0493
0494 AdaptiveMultiVertexFitter fitter(std::move(fitterCfg));
0495
0496
0497 Vector3 pos1a(0.5_mm, -0.5_mm, 2.4_mm);
0498 Vector3 mom1a(1000_MeV, 0_MeV, -500_MeV);
0499 Vector3 pos1b(0.5_mm, -0.5_mm, 3.5_mm);
0500 Vector3 mom1b(0_MeV, 1000_MeV, 500_MeV);
0501 Vector3 pos1c(-0.2_mm, 0.1_mm, 3.4_mm);
0502 Vector3 mom1c(-50_MeV, 180_MeV, 300_MeV);
0503
0504 Vector3 pos1d(-0.1_mm, 0.3_mm, 3.0_mm);
0505 Vector3 mom1d(-80_MeV, 480_MeV, -100_MeV);
0506 Vector3 pos1e(-0.01_mm, 0.01_mm, 2.9_mm);
0507 Vector3 mom1e(-600_MeV, 10_MeV, 210_MeV);
0508
0509 Vector3 pos1f(-0.07_mm, 0.03_mm, 2.5_mm);
0510 Vector3 mom1f(240_MeV, 110_MeV, 150_MeV);
0511
0512
0513 Covariance covMat1;
0514 covMat1 << 1_mm * 1_mm, 0, 0., 0, 0., 0, 0, 1_mm * 1_mm, 0, 0., 0, 0, 0., 0,
0515 0.1, 0, 0, 0, 0, 0., 0, 0.1, 0, 0, 0., 0, 0, 0, 1. / (10_GeV * 10_GeV), 0,
0516 0, 0, 0, 0, 0, 1_ns;
0517
0518 std::vector<BoundTrackParameters> params1 = {
0519 BoundTrackParameters::create(
0520 geoContext, Surface::makeShared<PerigeeSurface>(pos1a),
0521 makeVector4(pos1a, 0), mom1a.normalized(), 1_e / mom1a.norm(),
0522 covMat1, ParticleHypothesis::pion())
0523 .value(),
0524 BoundTrackParameters::create(
0525 geoContext, Surface::makeShared<PerigeeSurface>(pos1b),
0526 makeVector4(pos1b, 0), mom1b.normalized(), -1_e / mom1b.norm(),
0527 covMat1, ParticleHypothesis::pion())
0528 .value(),
0529 BoundTrackParameters::create(
0530 geoContext, Surface::makeShared<PerigeeSurface>(pos1c),
0531 makeVector4(pos1c, 0), mom1c.normalized(), 1_e / mom1c.norm(),
0532 covMat1, ParticleHypothesis::pion())
0533 .value(),
0534 BoundTrackParameters::create(
0535 geoContext, Surface::makeShared<PerigeeSurface>(pos1d),
0536 makeVector4(pos1d, 0), mom1d.normalized(), -1_e / mom1d.norm(),
0537 covMat1, ParticleHypothesis::pion())
0538 .value(),
0539 BoundTrackParameters::create(
0540 geoContext, Surface::makeShared<PerigeeSurface>(pos1e),
0541 makeVector4(pos1e, 0), mom1e.normalized(), 1_e / mom1e.norm(),
0542 covMat1, ParticleHypothesis::pion())
0543 .value(),
0544 BoundTrackParameters::create(
0545 geoContext, Surface::makeShared<PerigeeSurface>(pos1f),
0546 makeVector4(pos1f, 0), mom1f.normalized(), -1_e / mom1f.norm(),
0547 covMat1, ParticleHypothesis::pion())
0548 .value(),
0549 };
0550
0551
0552 Vector3 pos2a(0.2_mm, 0_mm, -4.9_mm);
0553 Vector3 mom2a(5000_MeV, 30_MeV, 200_MeV);
0554 Vector3 pos2b(-0.5_mm, 0.1_mm, -5.1_mm);
0555 Vector3 mom2b(800_MeV, 1200_MeV, 200_MeV);
0556 Vector3 pos2c(0.05_mm, -0.5_mm, -4.7_mm);
0557 Vector3 mom2c(400_MeV, -300_MeV, -200_MeV);
0558
0559
0560 Covariance covMat2 = covMat1;
0561
0562 std::vector<BoundTrackParameters> params2 = {
0563 BoundTrackParameters::create(
0564 geoContext, Surface::makeShared<PerigeeSurface>(pos2a),
0565 makeVector4(pos2a, 0), mom2a.normalized(), 1_e / mom2a.norm(),
0566 covMat2, ParticleHypothesis::pion())
0567 .value(),
0568 BoundTrackParameters::create(
0569 geoContext, Surface::makeShared<PerigeeSurface>(pos2b),
0570 makeVector4(pos2b, 0), mom2b.normalized(), -1_e / mom2b.norm(),
0571 covMat2, ParticleHypothesis::pion())
0572 .value(),
0573 BoundTrackParameters::create(
0574 geoContext, Surface::makeShared<PerigeeSurface>(pos2c),
0575 makeVector4(pos2c, 0), mom2c.normalized(), -1_e / mom2c.norm(),
0576 covMat2, ParticleHypothesis::pion())
0577 .value(),
0578 };
0579
0580 AdaptiveMultiVertexFitter::State state(*bField, magFieldContext);
0581
0582
0583 SquareMatrix4 covConstr(SquareMatrix4::Identity());
0584 covConstr = covConstr * 1e+8;
0585 covConstr(3, 3) = 0.;
0586
0587
0588 Vector3 vtxPos1(0.15_mm, 0.15_mm, 2.9_mm);
0589 Vertex vtx1(vtxPos1);
0590
0591
0592 state.vertexCollection.push_back(&vtx1);
0593
0594
0595 Vertex vtx1Constr(vtxPos1);
0596 vtx1Constr.setFullCovariance(covConstr);
0597 vtx1Constr.setFitQuality(0, -3);
0598
0599
0600 VertexInfo vtxInfo1;
0601 vtxInfo1.seedPosition = vtxInfo1.linPoint;
0602 vtxInfo1.linPoint.setZero();
0603 vtxInfo1.linPoint.head<3>() = vtxPos1;
0604 vtxInfo1.constraint = std::move(vtx1Constr);
0605 vtxInfo1.oldPosition = vtxInfo1.linPoint;
0606
0607 for (const auto& trk : params1) {
0608 vtxInfo1.trackLinks.push_back(InputTrack{&trk});
0609 state.tracksAtVerticesMap.insert(
0610 std::make_pair(std::make_pair(InputTrack{&trk}, &vtx1),
0611 TrackAtVertex(1.5, trk, InputTrack{&trk})));
0612 }
0613
0614
0615 Vector3 vtxPos2(0.3_mm, -0.2_mm, -4.8_mm);
0616 Vertex vtx2(vtxPos2);
0617
0618
0619 state.vertexCollection.push_back(&vtx2);
0620
0621
0622 Vertex vtx2Constr(vtxPos2);
0623 vtx2Constr.setFullCovariance(covConstr);
0624 vtx2Constr.setFitQuality(0, -3);
0625
0626
0627 VertexInfo vtxInfo2;
0628 vtxInfo2.linPoint.setZero();
0629 vtxInfo2.linPoint.head<3>() = vtxPos2;
0630 vtxInfo2.constraint = std::move(vtx2Constr);
0631 vtxInfo2.oldPosition = vtxInfo2.linPoint;
0632 vtxInfo2.seedPosition = vtxInfo2.linPoint;
0633
0634 for (const auto& trk : params2) {
0635 vtxInfo2.trackLinks.push_back(InputTrack{&trk});
0636 state.tracksAtVerticesMap.insert(
0637 std::make_pair(std::make_pair(InputTrack{&trk}, &vtx2),
0638 TrackAtVertex(1.5, trk, InputTrack{&trk})));
0639 }
0640
0641 state.vtxInfoMap[&vtx1] = std::move(vtxInfo1);
0642 state.vtxInfoMap[&vtx2] = std::move(vtxInfo2);
0643
0644 state.addVertexToMultiMap(vtx1);
0645 state.addVertexToMultiMap(vtx2);
0646
0647
0648 fitter.fit(state, vertexingOptions);
0649
0650 auto vtx1Fitted = state.vertexCollection.at(0);
0651 auto vtx1PosFitted = vtx1Fitted->position();
0652 auto vtx1CovFitted = vtx1Fitted->covariance();
0653 auto trks1 = state.vtxInfoMap.at(vtx1Fitted).trackLinks;
0654 auto vtx1FQ = vtx1Fitted->fitQuality();
0655
0656 auto vtx2Fitted = state.vertexCollection.at(1);
0657 auto vtx2PosFitted = vtx2Fitted->position();
0658 auto vtx2CovFitted = vtx2Fitted->covariance();
0659 auto trks2 = state.vtxInfoMap.at(vtx2Fitted).trackLinks;
0660 auto vtx2FQ = vtx2Fitted->fitQuality();
0661
0662
0663 ACTS_DEBUG("Vertex 1, position: " << vtx1PosFitted);
0664 ACTS_DEBUG("Vertex 1, covariance: " << vtx1CovFitted);
0665 for (const auto& trk : trks1) {
0666 auto& trkAtVtx =
0667 state.tracksAtVerticesMap.at(std::make_pair(trk, vtx1Fitted));
0668 ACTS_DEBUG("\tTrack weight:" << trkAtVtx.trackWeight);
0669 }
0670 ACTS_DEBUG("Vertex 1, chi2: " << vtx1FQ.first);
0671 ACTS_DEBUG("Vertex 1, ndf: " << vtx1FQ.second);
0672
0673
0674 ACTS_DEBUG("Vertex 2, position: " << vtx2PosFitted);
0675 ACTS_DEBUG("Vertex 2, covariance: " << vtx2CovFitted);
0676 for (const auto& trk : trks2) {
0677 auto& trkAtVtx =
0678 state.tracksAtVerticesMap.at(std::make_pair(trk, vtx2Fitted));
0679 ACTS_DEBUG("\tTrack weight:" << trkAtVtx.trackWeight);
0680 }
0681 ACTS_DEBUG("Vertex 2, chi2: " << vtx2FQ.first);
0682 ACTS_DEBUG("Vertex 2, ndf: " << vtx2FQ.second);
0683
0684
0685
0686 const Vector3 expVtx1Pos(0.077_mm, -0.189_mm, 2.924_mm);
0687
0688
0689 SquareMatrix3 expVtx1Cov;
0690 expVtx1Cov << 0.329, 0.016, -0.035, 0.016, 0.250, 0.085, -0.035, 0.085, 0.242;
0691
0692 ActsVector<6> expVtx1TrkWeights;
0693 expVtx1TrkWeights << 0.8128, 0.7994, 0.8164, 0.8165, 0.8165, 0.8119;
0694 const double expVtx1chi2 = 0.9812;
0695 const double expVtx1ndf = 6.7474;
0696
0697
0698 const Vector3 expVtx2Pos(-0.443_mm, -0.044_mm, -4.829_mm);
0699
0700 SquareMatrix3 expVtx2Cov;
0701 expVtx2Cov << 1.088, 0.028, -0.066, 0.028, 0.643, 0.073, -0.066, 0.073, 0.435;
0702
0703 const Vector3 expVtx2TrkWeights(0.8172, 0.8150, 0.8137);
0704 const double expVtx2chi2 = 0.2114;
0705 const double expVtx2ndf = 1.8920;
0706
0707
0708
0709 CHECK_CLOSE_ABS(vtx1PosFitted, expVtx1Pos, 0.001_mm);
0710 CHECK_CLOSE_ABS(vtx1CovFitted, expVtx1Cov, 0.001_mm);
0711 int trkCount = 0;
0712 for (const auto& trk : trks1) {
0713 auto& trkAtVtx =
0714 state.tracksAtVerticesMap.at(std::make_pair(trk, vtx1Fitted));
0715 CHECK_CLOSE_ABS(trkAtVtx.trackWeight, expVtx1TrkWeights[trkCount], 0.001);
0716 trkCount++;
0717 }
0718 CHECK_CLOSE_ABS(vtx1FQ.first, expVtx1chi2, 0.001);
0719 CHECK_CLOSE_ABS(vtx1FQ.second, expVtx1ndf, 0.001);
0720
0721
0722 CHECK_CLOSE_ABS(vtx2PosFitted, expVtx2Pos, 0.001_mm);
0723 CHECK_CLOSE_ABS(vtx2CovFitted, expVtx2Cov, 0.001_mm);
0724 trkCount = 0;
0725 for (const auto& trk : trks2) {
0726 auto& trkAtVtx =
0727 state.tracksAtVerticesMap.at(std::make_pair(trk, vtx2Fitted));
0728 CHECK_CLOSE_ABS(trkAtVtx.trackWeight, expVtx2TrkWeights[trkCount], 0.001);
0729 trkCount++;
0730 }
0731 CHECK_CLOSE_ABS(vtx2FQ.first, expVtx2chi2, 0.001);
0732 CHECK_CLOSE_ABS(vtx2FQ.second, expVtx2ndf, 0.001);
0733 }
0734
0735 BOOST_AUTO_TEST_SUITE_END()
0736
0737 }