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