File indexing completed on 2026-01-08 09:22:00
0001
0002
0003
0004
0005
0006
0007
0008
0009 #pragma once
0010
0011 #include "Acts/Definitions/Algebra.hpp"
0012 #include "Acts/EventData/MeasurementHelpers.hpp"
0013 #include "Acts/EventData/MultiTrajectory.hpp"
0014 #include "Acts/EventData/SourceLink.hpp"
0015 #include "Acts/EventData/TrackParameters.hpp"
0016 #include "Acts/EventData/TrackProxyConcept.hpp"
0017 #include "Acts/EventData/Types.hpp"
0018 #include "Acts/EventData/VectorMultiTrajectory.hpp"
0019 #include "Acts/EventData/VectorTrackContainer.hpp"
0020 #include "Acts/EventData/detail/CorrectedTransformationFreeToBound.hpp"
0021 #include "Acts/Geometry/GeometryContext.hpp"
0022 #include "Acts/Geometry/TrackingVolume.hpp"
0023 #include "Acts/MagneticField/MagneticFieldContext.hpp"
0024 #include "Acts/Material/Interactions.hpp"
0025 #include "Acts/Propagator/ActorList.hpp"
0026 #include "Acts/Propagator/DirectNavigator.hpp"
0027 #include "Acts/Propagator/PropagatorOptions.hpp"
0028 #include "Acts/Propagator/StandardAborters.hpp"
0029 #include "Acts/Propagator/detail/PointwiseMaterialInteraction.hpp"
0030 #include "Acts/TrackFitting/GlobalChiSquareFitterError.hpp"
0031 #include "Acts/TrackFitting/detail/VoidFitterComponents.hpp"
0032 #include "Acts/Utilities/CalibrationContext.hpp"
0033 #include "Acts/Utilities/Delegate.hpp"
0034 #include "Acts/Utilities/Logger.hpp"
0035 #include "Acts/Utilities/Result.hpp"
0036 #include "Acts/Utilities/TrackHelpers.hpp"
0037
0038 #include <functional>
0039 #include <limits>
0040 #include <map>
0041 #include <memory>
0042 #include <type_traits>
0043 #include <unordered_map>
0044
0045 namespace Acts::Experimental {
0046
0047 namespace Gx2fConstants {
0048 constexpr std::string_view gx2fnUpdateColumn = "Gx2fnUpdateColumn";
0049
0050
0051 constexpr TrackStatePropMask trackStateMask = TrackStatePropMask::Smoothed |
0052 TrackStatePropMask::Jacobian |
0053 TrackStatePropMask::Calibrated;
0054
0055
0056
0057 const Eigen::Matrix<double, eBoundSize, 2> phiThetaProjector = [] {
0058 Eigen::Matrix<double, eBoundSize, 2> m =
0059 Eigen::Matrix<double, eBoundSize, 2>::Zero();
0060 m(eBoundPhi, 0) = 1.0;
0061 m(eBoundTheta, 1) = 1.0;
0062 return m;
0063 }();
0064 }
0065
0066
0067 template <typename traj_t>
0068 struct Gx2FitterExtensions {
0069
0070 using TrackStateProxy = typename MultiTrajectory<traj_t>::TrackStateProxy;
0071
0072 using ConstTrackStateProxy =
0073 typename MultiTrajectory<traj_t>::ConstTrackStateProxy;
0074
0075 using Parameters = typename TrackStateProxy::Parameters;
0076
0077
0078 using Calibrator =
0079 Delegate<void(const GeometryContext&, const CalibrationContext&,
0080 const SourceLink&, TrackStateProxy)>;
0081
0082
0083
0084 using Updater = Delegate<Result<void>(const GeometryContext&, TrackStateProxy,
0085 const Logger&)>;
0086
0087
0088 using OutlierFinder = Delegate<bool(ConstTrackStateProxy)>;
0089
0090
0091
0092
0093 Calibrator calibrator;
0094
0095
0096 Updater updater;
0097
0098
0099
0100 OutlierFinder outlierFinder;
0101
0102
0103 SourceLinkSurfaceAccessor surfaceAccessor;
0104
0105
0106 Gx2FitterExtensions() {
0107 calibrator.template connect<&Acts::detail::voidFitterCalibrator<traj_t>>();
0108 updater.template connect<&Acts::detail::voidFitterUpdater<traj_t>>();
0109 outlierFinder.template connect<&Acts::detail::voidOutlierFinder<traj_t>>();
0110 surfaceAccessor.connect<&Acts::detail::voidSurfaceAccessor>();
0111 }
0112 };
0113
0114
0115
0116
0117 template <typename traj_t>
0118 struct Gx2FitterOptions {
0119
0120
0121
0122
0123
0124
0125
0126
0127
0128
0129
0130
0131
0132 Gx2FitterOptions(const GeometryContext& gctx,
0133 const MagneticFieldContext& mctx,
0134 std::reference_wrapper<const CalibrationContext> cctx,
0135 Gx2FitterExtensions<traj_t> extensions_,
0136 const PropagatorPlainOptions& pOptions,
0137 const Surface* rSurface = nullptr, bool mScattering = false,
0138 bool eLoss = false,
0139 const FreeToBoundCorrection& freeToBoundCorrection_ =
0140 FreeToBoundCorrection(false),
0141 const std::size_t nUpdateMax_ = 5,
0142 double relChi2changeCutOff_ = 1e-5)
0143 : geoContext(gctx),
0144 magFieldContext(mctx),
0145 calibrationContext(cctx),
0146 extensions(extensions_),
0147 propagatorPlainOptions(pOptions),
0148 referenceSurface(rSurface),
0149 multipleScattering(mScattering),
0150 energyLoss(eLoss),
0151 freeToBoundCorrection(freeToBoundCorrection_),
0152 nUpdateMax(nUpdateMax_),
0153 relChi2changeCutOff(relChi2changeCutOff_) {}
0154
0155
0156 Gx2FitterOptions() = delete;
0157
0158
0159 std::reference_wrapper<const GeometryContext> geoContext;
0160
0161 std::reference_wrapper<const MagneticFieldContext> magFieldContext;
0162
0163 std::reference_wrapper<const CalibrationContext> calibrationContext;
0164
0165
0166 Gx2FitterExtensions<traj_t> extensions;
0167
0168
0169 PropagatorPlainOptions propagatorPlainOptions;
0170
0171
0172 const Surface* referenceSurface = nullptr;
0173
0174
0175 bool multipleScattering = false;
0176
0177
0178 bool energyLoss = false;
0179
0180
0181
0182 FreeToBoundCorrection freeToBoundCorrection;
0183
0184
0185 std::size_t nUpdateMax = 5;
0186
0187
0188 double relChi2changeCutOff = 1e-7;
0189 };
0190
0191 template <typename traj_t>
0192 struct Gx2FitterResult {
0193
0194 traj_t* fittedStates{nullptr};
0195
0196
0197
0198
0199
0200 std::size_t lastMeasurementIndex = Acts::kTrackIndexInvalid;
0201
0202
0203
0204
0205
0206 std::size_t lastTrackIndex = Acts::kTrackIndexInvalid;
0207
0208
0209 std::optional<BoundTrackParameters> fittedParameters;
0210
0211
0212 std::size_t measurementStates = 0;
0213
0214
0215
0216
0217
0218 std::size_t measurementHoles = 0;
0219
0220
0221 std::size_t processedStates = 0;
0222
0223
0224 std::size_t processedMeasurements = 0;
0225
0226
0227 bool finished = false;
0228
0229
0230 std::vector<const Surface*> missedActiveSurfaces;
0231
0232
0233
0234 std::vector<const Surface*> passedAgainSurfaces;
0235
0236
0237 std::size_t surfaceCount = 0;
0238 };
0239
0240
0241
0242
0243
0244
0245 struct ScatteringProperties {
0246 public:
0247
0248
0249
0250
0251
0252 ScatteringProperties(const BoundVector& scatteringAngles_,
0253 const double invCovarianceMaterial_,
0254 const bool materialIsValid_)
0255 : m_scatteringAngles(scatteringAngles_),
0256 m_invCovarianceMaterial(invCovarianceMaterial_),
0257 m_materialIsValid(materialIsValid_) {}
0258
0259
0260
0261 const BoundVector& scatteringAngles() const { return m_scatteringAngles; }
0262
0263
0264
0265 BoundVector& scatteringAngles() { return m_scatteringAngles; }
0266
0267
0268
0269 double invCovarianceMaterial() const { return m_invCovarianceMaterial; }
0270
0271
0272
0273 bool materialIsValid() const { return m_materialIsValid; }
0274
0275 private:
0276
0277
0278 BoundVector m_scatteringAngles;
0279
0280
0281
0282 double m_invCovarianceMaterial;
0283
0284
0285
0286 bool m_materialIsValid;
0287 };
0288
0289
0290
0291
0292
0293 struct Gx2fSystem {
0294 public:
0295
0296
0297
0298 explicit Gx2fSystem(std::size_t nDims)
0299 : m_nDims{nDims},
0300 m_aMatrix{Eigen::MatrixXd::Zero(nDims, nDims)},
0301 m_bVector{Eigen::VectorXd::Zero(nDims)} {}
0302
0303
0304
0305 std::size_t nDims() const { return m_nDims; }
0306
0307
0308
0309 double chi2() const { return m_chi2; }
0310
0311
0312
0313 double& chi2() { return m_chi2; }
0314
0315
0316
0317 const Eigen::MatrixXd& aMatrix() const { return m_aMatrix; }
0318
0319
0320
0321 Eigen::MatrixXd& aMatrix() { return m_aMatrix; }
0322
0323
0324
0325 const Eigen::VectorXd& bVector() const { return m_bVector; }
0326
0327
0328
0329 Eigen::VectorXd& bVector() { return m_bVector; }
0330
0331
0332
0333 std::size_t ndf() const { return m_ndf; }
0334
0335
0336
0337 std::size_t& ndf() { return m_ndf; }
0338
0339
0340
0341
0342
0343
0344
0345
0346
0347
0348 std::size_t findRequiredNdf() {
0349 std::size_t ndfSystem = 0;
0350 if (m_aMatrix(4, 4) == 0) {
0351 ndfSystem = 4;
0352 } else if (m_aMatrix(5, 5) == 0) {
0353 ndfSystem = 5;
0354 } else {
0355 ndfSystem = 6;
0356 }
0357
0358 return ndfSystem;
0359 }
0360
0361
0362
0363 bool isWellDefined() { return m_ndf > findRequiredNdf(); }
0364
0365 private:
0366
0367 std::size_t m_nDims;
0368
0369
0370 double m_chi2 = 0.;
0371
0372
0373 Eigen::MatrixXd m_aMatrix;
0374
0375
0376 Eigen::VectorXd m_bVector;
0377
0378
0379 std::size_t m_ndf = 0u;
0380 };
0381
0382
0383
0384
0385
0386
0387
0388
0389
0390
0391
0392
0393
0394
0395
0396
0397
0398 void addMeasurementToGx2fSumsBackend(
0399 Gx2fSystem& extendedSystem,
0400 const std::vector<BoundMatrix>& jacobianFromStart,
0401 const Eigen::MatrixXd& covarianceMeasurement, const BoundVector& predicted,
0402 const Eigen::VectorXd& measurement, const Eigen::MatrixXd& projector,
0403 const Logger& logger);
0404
0405
0406
0407
0408
0409
0410
0411
0412
0413
0414
0415
0416
0417
0418 template <std::size_t kMeasDim, typename track_state_t>
0419 void addMeasurementToGx2fSums(Gx2fSystem& extendedSystem,
0420 const std::vector<BoundMatrix>& jacobianFromStart,
0421 const track_state_t& trackState,
0422 const Logger& logger) {
0423 const ActsSquareMatrix<kMeasDim> covarianceMeasurement =
0424 trackState.template calibratedCovariance<kMeasDim>();
0425
0426 const BoundVector predicted = trackState.smoothed();
0427
0428 const ActsVector<kMeasDim> measurement =
0429 trackState.template calibrated<kMeasDim>();
0430
0431 const ActsMatrix<kMeasDim, eBoundSize> projector =
0432 trackState.template projectorSubspaceHelper<kMeasDim>().projector();
0433
0434 addMeasurementToGx2fSumsBackend(extendedSystem, jacobianFromStart,
0435 covarianceMeasurement, predicted, measurement,
0436 projector, logger);
0437 }
0438
0439
0440
0441
0442
0443
0444
0445
0446
0447
0448
0449
0450
0451
0452 template <typename track_state_t>
0453 void addMaterialToGx2fSums(
0454 Gx2fSystem& extendedSystem, const std::size_t nMaterialsHandled,
0455 const std::unordered_map<GeometryIdentifier, ScatteringProperties>&
0456 scatteringMap,
0457 const track_state_t& trackState, const Logger& logger) {
0458
0459 const GeometryIdentifier geoId = trackState.referenceSurface().geometryId();
0460 const auto scatteringMapId = scatteringMap.find(geoId);
0461 if (scatteringMapId == scatteringMap.end()) {
0462 ACTS_ERROR("No scattering angles found for material surface " << geoId);
0463 throw std::runtime_error(
0464 "No scattering angles found for material surface.");
0465 }
0466
0467 const double sinThetaLoc = std::sin(trackState.smoothed()[eBoundTheta]);
0468
0469
0470 const std::size_t deltaPosition = eBoundSize + 2 * nMaterialsHandled;
0471
0472 const BoundVector& scatteringAngles =
0473 scatteringMapId->second.scatteringAngles();
0474
0475 const double invCov = scatteringMapId->second.invCovarianceMaterial();
0476
0477
0478 extendedSystem.aMatrix()(deltaPosition, deltaPosition) +=
0479 invCov * sinThetaLoc * sinThetaLoc;
0480 extendedSystem.bVector()(deltaPosition, 0) -=
0481 invCov * scatteringAngles[eBoundPhi] * sinThetaLoc;
0482 extendedSystem.chi2() += invCov * scatteringAngles[eBoundPhi] * sinThetaLoc *
0483 scatteringAngles[eBoundPhi] * sinThetaLoc;
0484
0485
0486 extendedSystem.aMatrix()(deltaPosition + 1, deltaPosition + 1) += invCov;
0487 extendedSystem.bVector()(deltaPosition + 1, 0) -=
0488 invCov * scatteringAngles[eBoundTheta];
0489 extendedSystem.chi2() +=
0490 invCov * scatteringAngles[eBoundTheta] * scatteringAngles[eBoundTheta];
0491
0492 ACTS_VERBOSE(
0493 "Contributions in addMaterialToGx2fSums:\n"
0494 << " invCov: " << invCov << "\n"
0495 << " sinThetaLoc: " << sinThetaLoc << "\n"
0496 << " deltaPosition: " << deltaPosition << "\n"
0497 << " Phi:\n"
0498 << " scattering angle: " << scatteringAngles[eBoundPhi] << "\n"
0499 << " aMatrix contribution: " << invCov * sinThetaLoc * sinThetaLoc
0500 << "\n"
0501 << " bVector contribution: "
0502 << invCov * scatteringAngles[eBoundPhi] * sinThetaLoc << "\n"
0503 << " chi2sum contribution: "
0504 << invCov * scatteringAngles[eBoundPhi] * sinThetaLoc *
0505 scatteringAngles[eBoundPhi] * sinThetaLoc
0506 << "\n"
0507 << " Theta:\n"
0508 << " scattering angle: " << scatteringAngles[eBoundTheta]
0509 << "\n"
0510 << " aMatrix contribution: " << invCov << "\n"
0511 << " bVector contribution: "
0512 << invCov * scatteringAngles[eBoundTheta] << "\n"
0513 << " chi2sum contribution: "
0514 << invCov * scatteringAngles[eBoundTheta] * scatteringAngles[eBoundTheta]
0515 << "\n");
0516
0517 return;
0518 }
0519
0520
0521
0522
0523
0524
0525
0526
0527
0528
0529
0530
0531
0532
0533
0534
0535 template <TrackProxyConcept track_proxy_t>
0536 void fillGx2fSystem(
0537 const track_proxy_t track, Gx2fSystem& extendedSystem,
0538 const bool multipleScattering,
0539 const std::unordered_map<GeometryIdentifier, ScatteringProperties>&
0540 scatteringMap,
0541 std::vector<GeometryIdentifier>& geoIdVector, const Logger& logger) {
0542 std::vector<BoundMatrix> jacobianFromStart;
0543 jacobianFromStart.emplace_back(BoundMatrix::Identity());
0544
0545 for (const auto& trackState : track.trackStates()) {
0546
0547 const GeometryIdentifier geoId = trackState.referenceSurface().geometryId();
0548 ACTS_DEBUG("Start to investigate trackState on surface " << geoId);
0549 const auto typeFlags = trackState.typeFlags();
0550 const bool stateHasMeasurement =
0551 typeFlags.test(TrackStateFlag::MeasurementFlag);
0552 const bool stateHasMaterial = typeFlags.test(TrackStateFlag::MaterialFlag);
0553
0554
0555
0556
0557
0558 bool doMaterial = multipleScattering && stateHasMaterial;
0559 if (doMaterial) {
0560 const auto scatteringMapId = scatteringMap.find(geoId);
0561 assert(scatteringMapId != scatteringMap.end() &&
0562 "No scattering angles found for material surface.");
0563 doMaterial = doMaterial && scatteringMapId->second.materialIsValid();
0564 }
0565
0566
0567 if (!stateHasMeasurement && !doMaterial) {
0568 ACTS_DEBUG(" Skip state.");
0569 continue;
0570 }
0571
0572
0573 for (auto& jac : jacobianFromStart) {
0574 jac = trackState.jacobian() * jac;
0575 }
0576
0577
0578 if (stateHasMeasurement) {
0579 ACTS_DEBUG(" Handle measurement.");
0580
0581 const auto measDim = trackState.calibratedSize();
0582
0583 if (measDim < 1 || 6 < measDim) {
0584 ACTS_ERROR("Can not process state with measurement with "
0585 << measDim << " dimensions.");
0586 throw std::domain_error(
0587 "Found measurement with less than 1 or more than 6 dimension(s).");
0588 }
0589
0590 extendedSystem.ndf() += measDim;
0591
0592 visit_measurement(measDim, [&](auto N) {
0593 addMeasurementToGx2fSums<N>(extendedSystem, jacobianFromStart,
0594 trackState, logger);
0595 });
0596 }
0597
0598
0599 if (doMaterial) {
0600 ACTS_DEBUG(" Handle material");
0601
0602 jacobianFromStart.emplace_back(BoundMatrix::Identity());
0603
0604
0605 addMaterialToGx2fSums(extendedSystem, geoIdVector.size(), scatteringMap,
0606 trackState, logger);
0607
0608 geoIdVector.emplace_back(geoId);
0609 }
0610 }
0611 }
0612
0613
0614
0615
0616
0617
0618
0619
0620
0621
0622
0623
0624
0625
0626 template <TrackProxyConcept track_proxy_t>
0627 std::size_t countMaterialStates(
0628 const track_proxy_t track,
0629 const std::unordered_map<GeometryIdentifier, ScatteringProperties>&
0630 scatteringMap,
0631 const Logger& logger) {
0632 std::size_t nMaterialSurfaces = 0;
0633 ACTS_DEBUG("Count the valid material surfaces.");
0634 for (const auto& trackState : track.trackStates()) {
0635 const auto typeFlags = trackState.typeFlags();
0636 const bool stateHasMaterial = typeFlags.test(TrackStateFlag::MaterialFlag);
0637
0638 if (!stateHasMaterial) {
0639 continue;
0640 }
0641
0642
0643 const GeometryIdentifier geoId = trackState.referenceSurface().geometryId();
0644
0645 const auto scatteringMapId = scatteringMap.find(geoId);
0646 assert(scatteringMapId != scatteringMap.end() &&
0647 "No scattering angles found for material surface.");
0648 if (!scatteringMapId->second.materialIsValid()) {
0649 continue;
0650 }
0651
0652 nMaterialSurfaces++;
0653 }
0654
0655 return nMaterialSurfaces;
0656 }
0657
0658
0659
0660
0661
0662
0663
0664
0665 Eigen::VectorXd computeGx2fDeltaParams(const Gx2fSystem& extendedSystem);
0666
0667
0668
0669
0670
0671
0672
0673
0674
0675 void updateGx2fParams(
0676 BoundTrackParameters& params, const Eigen::VectorXd& deltaParamsExtended,
0677 const std::size_t nMaterialSurfaces,
0678 std::unordered_map<GeometryIdentifier, ScatteringProperties>& scatteringMap,
0679 const std::vector<GeometryIdentifier>& geoIdVector);
0680
0681
0682
0683
0684
0685
0686
0687
0688
0689
0690
0691 void updateGx2fCovarianceParams(BoundMatrix& fullCovariancePredicted,
0692 Gx2fSystem& extendedSystem);
0693
0694
0695
0696
0697
0698
0699 template <typename propagator_t, typename traj_t>
0700 class Gx2Fitter {
0701
0702 using Gx2fNavigator = typename propagator_t::Navigator;
0703
0704
0705 static constexpr bool isDirectNavigator =
0706 std::is_same_v<Gx2fNavigator, DirectNavigator>;
0707
0708 static constexpr auto kInvalid = kTrackIndexInvalid;
0709
0710 public:
0711
0712
0713
0714
0715
0716
0717
0718
0719 explicit Gx2Fitter(propagator_t pPropagator,
0720 std::unique_ptr<const Logger> _logger =
0721 getDefaultLogger("Gx2Fitter", Logging::INFO))
0722 : m_propagator(std::move(pPropagator)),
0723 m_logger{std::move(_logger)},
0724 m_actorLogger{m_logger->cloneWithSuffix("Actor")},
0725 m_addToSumLogger{m_logger->cloneWithSuffix("AddToSum")} {}
0726
0727 private:
0728
0729 propagator_t m_propagator;
0730
0731
0732 std::unique_ptr<const Logger> m_logger;
0733 std::unique_ptr<const Logger> m_actorLogger;
0734 std::unique_ptr<const Logger> m_addToSumLogger;
0735
0736 const Logger& logger() const { return *m_logger; }
0737
0738
0739
0740
0741
0742
0743
0744
0745
0746 template <typename parameters_t>
0747 class Actor {
0748 public:
0749
0750 using result_type = Gx2FitterResult<traj_t>;
0751
0752
0753 const Surface* targetSurface = nullptr;
0754
0755
0756 const std::map<GeometryIdentifier, SourceLink>* inputMeasurements = nullptr;
0757
0758
0759 bool multipleScattering = false;
0760
0761
0762 bool energyLoss = false;
0763
0764
0765
0766 FreeToBoundCorrection freeToBoundCorrection;
0767
0768
0769 std::shared_ptr<MultiTrajectory<traj_t>> outputStates;
0770
0771
0772 const Logger* actorLogger{nullptr};
0773
0774
0775 const Logger& logger() const { return *actorLogger; }
0776
0777 Gx2FitterExtensions<traj_t> extensions;
0778
0779
0780 SurfaceReached targetReached;
0781
0782
0783 const CalibrationContext* calibrationContext{nullptr};
0784
0785
0786 const parameters_t* parametersWithHypothesis = nullptr;
0787
0788
0789
0790 std::unordered_map<GeometryIdentifier, ScatteringProperties>*
0791 scatteringMap = nullptr;
0792
0793
0794
0795
0796
0797
0798
0799
0800
0801
0802
0803 template <typename propagator_state_t, typename stepper_t,
0804 typename navigator_t>
0805 Result<void> act(propagator_state_t& state, const stepper_t& stepper,
0806 const navigator_t& navigator, result_type& result,
0807 const Logger& ) const {
0808 assert(result.fittedStates && "No MultiTrajectory set");
0809
0810
0811 if (result.measurementStates == inputMeasurements->size()) {
0812 ACTS_DEBUG("Actor: finish: All measurements have been found.");
0813 result.finished = true;
0814 } else if (state.navigation.navigationBreak) {
0815 ACTS_DEBUG("Actor: finish: navigationBreak.");
0816 result.finished = true;
0817 }
0818
0819
0820 if (result.finished) {
0821
0822 if (result.measurementStates > 0) {
0823 result.missedActiveSurfaces.resize(result.measurementHoles);
0824 }
0825
0826 return Result<void>::success();
0827 }
0828
0829
0830
0831 auto surface = navigator.currentSurface(state.navigation);
0832 if (surface == nullptr) {
0833 return Result<void>::success();
0834 }
0835
0836 ++result.surfaceCount;
0837 const GeometryIdentifier geoId = surface->geometryId();
0838 ACTS_DEBUG("Surface " << geoId << " detected.");
0839
0840 const bool surfaceIsSensitive =
0841 (surface->associatedDetectorElement() != nullptr);
0842 const bool surfaceHasMaterial = (surface->surfaceMaterial() != nullptr);
0843
0844
0845
0846 bool doMaterial = multipleScattering && surfaceHasMaterial;
0847
0848
0849
0850 if (doMaterial) {
0851 ACTS_DEBUG(" The surface contains material, ...");
0852
0853 auto scatteringMapId = scatteringMap->find(geoId);
0854 if (scatteringMapId == scatteringMap->end()) {
0855 ACTS_DEBUG(" ... create entry in scattering map.");
0856
0857 Acts::detail::PointwiseMaterialInteraction interaction(surface, state,
0858 stepper);
0859
0860 const bool slabIsValid = interaction.evaluateMaterialSlab(
0861 state, navigator, MaterialUpdateStage::FullUpdate);
0862 double invSigma2 = 0.;
0863 if (slabIsValid) {
0864 const auto& particle =
0865 parametersWithHypothesis->particleHypothesis();
0866
0867 const double sigma =
0868 static_cast<double>(Acts::computeMultipleScatteringTheta0(
0869 interaction.slab, particle.absolutePdg(), particle.mass(),
0870 static_cast<float>(
0871 parametersWithHypothesis->parameters()[eBoundQOverP]),
0872 particle.absoluteCharge()));
0873 ACTS_VERBOSE(
0874 " The Highland formula gives sigma = " << sigma);
0875 invSigma2 = 1. / std::pow(sigma, 2);
0876 } else {
0877 ACTS_VERBOSE(" Material slab is not valid.");
0878 }
0879
0880 scatteringMap->emplace(
0881 geoId, ScatteringProperties{BoundVector::Zero(), invSigma2,
0882 slabIsValid});
0883 scatteringMapId = scatteringMap->find(geoId);
0884 } else {
0885 ACTS_DEBUG(" ... found entry in scattering map.");
0886 }
0887
0888 doMaterial = doMaterial && scatteringMapId->second.materialIsValid();
0889 }
0890
0891
0892 if (auto sourceLinkIt = inputMeasurements->find(geoId);
0893 sourceLinkIt != inputMeasurements->end()) {
0894 ACTS_DEBUG(" The surface contains a measurement.");
0895
0896
0897 stepper.transportCovarianceToBound(state.stepping, *surface,
0898 freeToBoundCorrection);
0899
0900
0901 auto& fittedStates = *result.fittedStates;
0902
0903
0904
0905 typename traj_t::TrackStateProxy trackStateProxy =
0906 fittedStates.makeTrackState(Gx2fConstants::trackStateMask,
0907 result.lastTrackIndex);
0908 const std::size_t currentTrackIndex = trackStateProxy.index();
0909
0910
0911
0912 {
0913 trackStateProxy.setReferenceSurface(surface->getSharedPtr());
0914
0915 auto res = stepper.boundState(state.stepping, *surface, false,
0916 freeToBoundCorrection);
0917 if (!res.ok()) {
0918 return res.error();
0919 }
0920
0921 auto& [boundParams, jacobian, pathLength] = *res;
0922
0923
0924
0925 if (doMaterial) {
0926 ACTS_DEBUG(" Update parameters with scattering angles.");
0927 const auto scatteringMapId = scatteringMap->find(geoId);
0928 ACTS_VERBOSE(
0929 " scatteringAngles: "
0930 << scatteringMapId->second.scatteringAngles().transpose());
0931 ACTS_VERBOSE(" boundParams before the update: "
0932 << boundParams.parameters().transpose());
0933 boundParams.parameters() +=
0934 scatteringMapId->second.scatteringAngles();
0935 ACTS_VERBOSE(" boundParams after the update: "
0936 << boundParams.parameters().transpose());
0937 }
0938
0939
0940 trackStateProxy.smoothed() = boundParams.parameters();
0941 trackStateProxy.smoothedCovariance() = state.stepping.cov;
0942
0943 trackStateProxy.jacobian() = jacobian;
0944 trackStateProxy.pathLength() = pathLength;
0945
0946 if (doMaterial) {
0947 stepper.update(state.stepping,
0948 transformBoundToFreeParameters(
0949 trackStateProxy.referenceSurface(),
0950 state.geoContext, trackStateProxy.smoothed()),
0951 trackStateProxy.smoothed(),
0952 trackStateProxy.smoothedCovariance(), *surface);
0953 }
0954 }
0955
0956
0957
0958 extensions.calibrator(state.geoContext, *calibrationContext,
0959 sourceLinkIt->second, trackStateProxy);
0960
0961
0962 auto typeFlags = trackStateProxy.typeFlags();
0963 typeFlags.set(TrackStateFlag::ParameterFlag);
0964 if (surfaceHasMaterial) {
0965 typeFlags.set(TrackStateFlag::MaterialFlag);
0966 }
0967
0968
0969 typeFlags.set(TrackStateFlag::MeasurementFlag);
0970
0971 ++result.processedMeasurements;
0972
0973 result.lastMeasurementIndex = currentTrackIndex;
0974 result.lastTrackIndex = currentTrackIndex;
0975
0976
0977
0978 ++result.measurementStates;
0979
0980
0981 ++result.processedStates;
0982
0983
0984
0985 result.measurementHoles = result.missedActiveSurfaces.size();
0986
0987 return Result<void>::success();
0988 }
0989
0990 if (doMaterial) {
0991
0992
0993
0994 ACTS_DEBUG(
0995 " The surface contains no measurement, but material and maybe "
0996 "a hole.");
0997
0998
0999 stepper.transportCovarianceToBound(state.stepping, *surface,
1000 freeToBoundCorrection);
1001
1002
1003 auto& fittedStates = *result.fittedStates;
1004
1005
1006
1007 typename traj_t::TrackStateProxy trackStateProxy =
1008 fittedStates.makeTrackState(Gx2fConstants::trackStateMask,
1009 result.lastTrackIndex);
1010 const std::size_t currentTrackIndex = trackStateProxy.index();
1011
1012
1013
1014 {
1015 trackStateProxy.setReferenceSurface(surface->getSharedPtr());
1016
1017 auto res = stepper.boundState(state.stepping, *surface, false,
1018 freeToBoundCorrection);
1019 if (!res.ok()) {
1020 return res.error();
1021 }
1022
1023 auto& [boundParams, jacobian, pathLength] = *res;
1024
1025
1026
1027
1028
1029 ACTS_DEBUG(" Update parameters with scattering angles.");
1030 const auto scatteringMapId = scatteringMap->find(geoId);
1031 ACTS_VERBOSE(
1032 " scatteringAngles: "
1033 << scatteringMapId->second.scatteringAngles().transpose());
1034 ACTS_VERBOSE(" boundParams before the update: "
1035 << boundParams.parameters().transpose());
1036 boundParams.parameters() +=
1037 scatteringMapId->second.scatteringAngles();
1038 ACTS_VERBOSE(" boundParams after the update: "
1039 << boundParams.parameters().transpose());
1040
1041
1042 trackStateProxy.smoothed() = boundParams.parameters();
1043 trackStateProxy.smoothedCovariance() = state.stepping.cov;
1044
1045 trackStateProxy.jacobian() = jacobian;
1046 trackStateProxy.pathLength() = pathLength;
1047
1048 stepper.update(state.stepping,
1049 transformBoundToFreeParameters(
1050 trackStateProxy.referenceSurface(),
1051 state.geoContext, trackStateProxy.smoothed()),
1052 trackStateProxy.smoothed(),
1053 trackStateProxy.smoothedCovariance(), *surface);
1054 }
1055
1056
1057 auto typeFlags = trackStateProxy.typeFlags();
1058 typeFlags.set(TrackStateFlag::ParameterFlag);
1059 typeFlags.set(TrackStateFlag::MaterialFlag);
1060
1061
1062
1063 const bool precedingMeasurementExists = (result.measurementStates > 0);
1064 if (surfaceIsSensitive && precedingMeasurementExists) {
1065 ACTS_DEBUG(" Surface is also sensitive. Marked as hole.");
1066 typeFlags.set(TrackStateFlag::HoleFlag);
1067
1068
1069 result.missedActiveSurfaces.push_back(surface);
1070 }
1071
1072 result.lastTrackIndex = currentTrackIndex;
1073
1074 ++result.processedStates;
1075
1076 return Result<void>::success();
1077 }
1078
1079 if (surfaceIsSensitive || surfaceHasMaterial) {
1080
1081
1082
1083 if (multipleScattering) {
1084 ACTS_DEBUG(
1085 " The surface contains no measurement, but maybe a hole.");
1086 } else {
1087 ACTS_DEBUG(
1088 " The surface contains no measurement, but maybe a hole "
1089 "and/or material.");
1090 }
1091
1092
1093
1094
1095 const bool precedingMeasurementExists = (result.measurementStates > 0);
1096 if (!precedingMeasurementExists && !surfaceHasMaterial) {
1097 ACTS_DEBUG(
1098 " Ignoring hole, because there are no preceding "
1099 "measurements.");
1100 return Result<void>::success();
1101 }
1102
1103 auto& fittedStates = *result.fittedStates;
1104
1105
1106
1107 typename traj_t::TrackStateProxy trackStateProxy =
1108 fittedStates.makeTrackState(Gx2fConstants::trackStateMask,
1109 result.lastTrackIndex);
1110 const std::size_t currentTrackIndex = trackStateProxy.index();
1111
1112
1113
1114 {
1115 trackStateProxy.setReferenceSurface(surface->getSharedPtr());
1116
1117 auto res = stepper.boundState(state.stepping, *surface, false,
1118 freeToBoundCorrection);
1119 if (!res.ok()) {
1120 return res.error();
1121 }
1122 const auto& [boundParams, jacobian, pathLength] = *res;
1123
1124
1125 trackStateProxy.smoothed() = boundParams.parameters();
1126 trackStateProxy.smoothedCovariance() = state.stepping.cov;
1127
1128 trackStateProxy.jacobian() = jacobian;
1129 trackStateProxy.pathLength() = pathLength;
1130 }
1131
1132
1133 auto typeFlags = trackStateProxy.typeFlags();
1134 typeFlags.set(TrackStateFlag::ParameterFlag);
1135 if (surfaceHasMaterial) {
1136 ACTS_DEBUG(" It is material.");
1137 typeFlags.set(TrackStateFlag::MaterialFlag);
1138 }
1139
1140
1141 if (surfaceIsSensitive && precedingMeasurementExists) {
1142 ACTS_DEBUG(" It is a hole.");
1143 typeFlags.set(TrackStateFlag::HoleFlag);
1144
1145 result.missedActiveSurfaces.push_back(surface);
1146 }
1147
1148 result.lastTrackIndex = currentTrackIndex;
1149
1150 ++result.processedStates;
1151
1152 return Result<void>::success();
1153 }
1154
1155 ACTS_DEBUG(" The surface contains no measurement/material/hole.");
1156 return Result<void>::success();
1157 }
1158
1159 template <typename propagator_state_t, typename stepper_t,
1160 typename navigator_t, typename result_t>
1161 bool checkAbort(propagator_state_t& , const stepper_t& ,
1162 const navigator_t& , const result_t& result,
1163 const Logger& ) const {
1164 if (result.finished) {
1165 return true;
1166 }
1167 return false;
1168 }
1169 };
1170
1171 public:
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190 template <typename source_link_iterator_t, typename start_parameters_t,
1191 typename parameters_t = BoundTrackParameters,
1192 TrackContainerFrontend track_container_t>
1193 Result<typename track_container_t::TrackProxy> fit(
1194 source_link_iterator_t it, source_link_iterator_t end,
1195 const start_parameters_t& sParameters,
1196 const Gx2FitterOptions<traj_t>& gx2fOptions,
1197 track_container_t& trackContainer) const
1198 requires(!isDirectNavigator)
1199 {
1200
1201
1202
1203 ACTS_VERBOSE("Preparing " << std::distance(it, end)
1204 << " input measurements");
1205 std::map<GeometryIdentifier, SourceLink> inputMeasurements;
1206
1207 for (; it != end; ++it) {
1208 SourceLink sl = *it;
1209 auto geoId = gx2fOptions.extensions.surfaceAccessor(sl)->geometryId();
1210 inputMeasurements.emplace(geoId, std::move(sl));
1211 }
1212
1213
1214
1215 const bool multipleScattering = gx2fOptions.multipleScattering;
1216
1217
1218 using GX2FActor = Actor<parameters_t>;
1219
1220 using GX2FResult = typename GX2FActor::result_type;
1221 using Actors = Acts::ActorList<GX2FActor>;
1222
1223 using PropagatorOptions = typename propagator_t::template Options<Actors>;
1224
1225 start_parameters_t params = sParameters;
1226 double chi2sum = 0;
1227 double oldChi2sum = std::numeric_limits<double>::max();
1228
1229
1230
1231
1232
1233 typename track_container_t::TrackContainerBackend trackContainerTempBackend;
1234 traj_t trajectoryTempBackend;
1235 TrackContainer trackContainerTemp{trackContainerTempBackend,
1236 trajectoryTempBackend};
1237
1238
1239
1240
1241 std::size_t tipIndex = kInvalid;
1242
1243
1244
1245 std::unordered_map<GeometryIdentifier, ScatteringProperties> scatteringMap;
1246
1247
1248
1249 BoundMatrix fullCovariancePredicted = BoundMatrix::Identity();
1250
1251 ACTS_VERBOSE("Initial parameters: " << params.parameters().transpose());
1252
1253
1254 ACTS_DEBUG("Start to iterate");
1255
1256
1257
1258
1259 std::size_t nUpdate = 0;
1260 for (nUpdate = 0; nUpdate < gx2fOptions.nUpdateMax; nUpdate++) {
1261 ACTS_DEBUG("nUpdate = " << nUpdate + 1 << "/" << gx2fOptions.nUpdateMax);
1262
1263
1264 PropagatorOptions propagatorOptions{gx2fOptions.propagatorPlainOptions};
1265
1266
1267
1268 for (const auto& [surfaceId, _] : inputMeasurements) {
1269 propagatorOptions.navigation.insertExternalSurface(surfaceId);
1270 }
1271
1272 auto& gx2fActor = propagatorOptions.actorList.template get<GX2FActor>();
1273 gx2fActor.inputMeasurements = &inputMeasurements;
1274 gx2fActor.multipleScattering = false;
1275 gx2fActor.extensions = gx2fOptions.extensions;
1276 gx2fActor.calibrationContext = &gx2fOptions.calibrationContext.get();
1277 gx2fActor.actorLogger = m_actorLogger.get();
1278 gx2fActor.scatteringMap = &scatteringMap;
1279 gx2fActor.parametersWithHypothesis = ¶ms;
1280
1281 auto propagatorState = m_propagator.makeState(propagatorOptions);
1282
1283 auto propagatorInitResult =
1284 m_propagator.initialize(propagatorState, params);
1285 if (!propagatorInitResult.ok()) {
1286 ACTS_ERROR("Propagation initialization failed: "
1287 << propagatorInitResult.error());
1288 return propagatorInitResult.error();
1289 }
1290
1291 auto& r = propagatorState.template get<Gx2FitterResult<traj_t>>();
1292 r.fittedStates = &trajectoryTempBackend;
1293
1294
1295
1296 trackContainerTemp.clear();
1297
1298 auto propagationResult = m_propagator.propagate(propagatorState);
1299
1300
1301 auto result =
1302 m_propagator.makeResult(std::move(propagatorState), propagationResult,
1303 propagatorOptions, false);
1304
1305 if (!result.ok()) {
1306 ACTS_ERROR("Propagation failed: " << result.error());
1307 return result.error();
1308 }
1309
1310
1311
1312 auto& propRes = *result;
1313 GX2FResult gx2fResult = std::move(propRes.template get<GX2FResult>());
1314
1315 auto track = trackContainerTemp.makeTrack();
1316 tipIndex = gx2fResult.lastMeasurementIndex;
1317
1318
1319
1320
1321 if (tipIndex == kInvalid) {
1322 ACTS_INFO("Did not find any measurements in nUpdate "
1323 << nUpdate + 1 << "/" << gx2fOptions.nUpdateMax);
1324 return Experimental::GlobalChiSquareFitterError::NotEnoughMeasurements;
1325 }
1326
1327 track.tipIndex() = tipIndex;
1328 track.linkForward();
1329
1330
1331
1332 const std::size_t nMaterialSurfaces = 0u;
1333
1334
1335
1336 const std::size_t dimsExtendedParams = eBoundSize + 2 * nMaterialSurfaces;
1337
1338
1339
1340 Gx2fSystem extendedSystem{dimsExtendedParams};
1341
1342
1343
1344
1345
1346 std::vector<GeometryIdentifier> geoIdVector;
1347
1348 fillGx2fSystem(track, extendedSystem, false, scatteringMap, geoIdVector,
1349 *m_addToSumLogger);
1350
1351 chi2sum = extendedSystem.chi2();
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365 if ((nUpdate > 0) && !extendedSystem.isWellDefined()) {
1366 ACTS_INFO("Not enough measurements. Require "
1367 << extendedSystem.findRequiredNdf() + 1 << ", but only "
1368 << extendedSystem.ndf() << " could be used.");
1369 return Experimental::GlobalChiSquareFitterError::NotEnoughMeasurements;
1370 }
1371
1372 Eigen::VectorXd deltaParamsExtended =
1373 computeGx2fDeltaParams(extendedSystem);
1374
1375 ACTS_VERBOSE("aMatrix:\n"
1376 << extendedSystem.aMatrix() << "\n"
1377 << "bVector:\n"
1378 << extendedSystem.bVector() << "\n"
1379 << "deltaParamsExtended:\n"
1380 << deltaParamsExtended << "\n"
1381 << "oldChi2sum = " << oldChi2sum << "\n"
1382 << "chi2sum = " << extendedSystem.chi2());
1383
1384 if ((gx2fOptions.relChi2changeCutOff != 0) && (nUpdate > 0) &&
1385 (std::abs(extendedSystem.chi2() / oldChi2sum - 1) <
1386 gx2fOptions.relChi2changeCutOff)) {
1387 ACTS_DEBUG("Abort with relChi2changeCutOff after "
1388 << nUpdate + 1 << "/" << gx2fOptions.nUpdateMax
1389 << " iterations.");
1390 updateGx2fCovarianceParams(fullCovariancePredicted, extendedSystem);
1391 break;
1392 }
1393
1394 if (extendedSystem.chi2() > oldChi2sum + 1e-5) {
1395 ACTS_DEBUG("chi2 not converging monotonically in update " << nUpdate);
1396 }
1397
1398
1399
1400 if (nUpdate == gx2fOptions.nUpdateMax - 1) {
1401
1402
1403
1404
1405
1406 if (gx2fOptions.nUpdateMax > 5) {
1407 ACTS_INFO("Did not converge in " << gx2fOptions.nUpdateMax
1408 << " updates.");
1409 return Experimental::GlobalChiSquareFitterError::DidNotConverge;
1410 }
1411
1412 updateGx2fCovarianceParams(fullCovariancePredicted, extendedSystem);
1413 break;
1414 }
1415
1416 updateGx2fParams(params, deltaParamsExtended, nMaterialSurfaces,
1417 scatteringMap, geoIdVector);
1418 ACTS_VERBOSE("Updated parameters: " << params.parameters().transpose());
1419
1420 oldChi2sum = extendedSystem.chi2();
1421 }
1422 ACTS_DEBUG("Finished to iterate");
1423 ACTS_VERBOSE("Final parameters: " << params.parameters().transpose());
1424
1425
1426
1427 ACTS_DEBUG("Start to evaluate material");
1428 if (multipleScattering) {
1429
1430 PropagatorOptions propagatorOptions{gx2fOptions.propagatorPlainOptions};
1431
1432
1433
1434 for (const auto& [surfaceId, _] : inputMeasurements) {
1435 propagatorOptions.navigation.insertExternalSurface(surfaceId);
1436 }
1437
1438 auto& gx2fActor = propagatorOptions.actorList.template get<GX2FActor>();
1439 gx2fActor.inputMeasurements = &inputMeasurements;
1440 gx2fActor.multipleScattering = true;
1441 gx2fActor.extensions = gx2fOptions.extensions;
1442 gx2fActor.calibrationContext = &gx2fOptions.calibrationContext.get();
1443 gx2fActor.actorLogger = m_actorLogger.get();
1444 gx2fActor.scatteringMap = &scatteringMap;
1445 gx2fActor.parametersWithHypothesis = ¶ms;
1446
1447 auto propagatorState = m_propagator.makeState(propagatorOptions);
1448
1449 auto propagatorInitResult =
1450 m_propagator.initialize(propagatorState, params);
1451 if (!propagatorInitResult.ok()) {
1452 ACTS_ERROR("Propagation initialization failed: "
1453 << propagatorInitResult.error());
1454 return propagatorInitResult.error();
1455 }
1456
1457 auto& r = propagatorState.template get<Gx2FitterResult<traj_t>>();
1458 r.fittedStates = &trajectoryTempBackend;
1459
1460
1461
1462 trackContainerTemp.clear();
1463
1464 auto propagationResult = m_propagator.propagate(propagatorState);
1465
1466
1467 auto result =
1468 m_propagator.makeResult(std::move(propagatorState), propagationResult,
1469 propagatorOptions, false);
1470
1471 if (!result.ok()) {
1472 ACTS_ERROR("Propagation failed: " << result.error());
1473 return result.error();
1474 }
1475
1476
1477
1478 auto& propRes = *result;
1479 GX2FResult gx2fResult = std::move(propRes.template get<GX2FResult>());
1480
1481 auto track = trackContainerTemp.makeTrack();
1482 tipIndex = gx2fResult.lastMeasurementIndex;
1483
1484
1485
1486
1487 if (tipIndex == kInvalid) {
1488 ACTS_INFO("Did not find any measurements in material fit.");
1489 return Experimental::GlobalChiSquareFitterError::NotEnoughMeasurements;
1490 }
1491
1492 track.tipIndex() = tipIndex;
1493 track.linkForward();
1494
1495
1496
1497 const std::size_t nMaterialSurfaces =
1498 countMaterialStates(track, scatteringMap, *m_addToSumLogger);
1499
1500
1501
1502 const std::size_t dimsExtendedParams = eBoundSize + 2 * nMaterialSurfaces;
1503
1504
1505
1506 Gx2fSystem extendedSystem{dimsExtendedParams};
1507
1508
1509
1510
1511
1512 std::vector<GeometryIdentifier> geoIdVector;
1513
1514 fillGx2fSystem(track, extendedSystem, true, scatteringMap, geoIdVector,
1515 *m_addToSumLogger);
1516
1517 chi2sum = extendedSystem.chi2();
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528 if ((nUpdate > 0) && !extendedSystem.isWellDefined()) {
1529 ACTS_INFO("Not enough measurements. Require "
1530 << extendedSystem.findRequiredNdf() + 1 << ", but only "
1531 << extendedSystem.ndf() << " could be used.");
1532 return Experimental::GlobalChiSquareFitterError::NotEnoughMeasurements;
1533 }
1534
1535 Eigen::VectorXd deltaParamsExtended =
1536 computeGx2fDeltaParams(extendedSystem);
1537
1538 ACTS_VERBOSE("aMatrix:\n"
1539 << extendedSystem.aMatrix() << "\n"
1540 << "bVector:\n"
1541 << extendedSystem.bVector() << "\n"
1542 << "deltaParamsExtended:\n"
1543 << deltaParamsExtended << "\n"
1544 << "oldChi2sum = " << oldChi2sum << "\n"
1545 << "chi2sum = " << extendedSystem.chi2());
1546
1547 chi2sum = extendedSystem.chi2();
1548
1549 updateGx2fParams(params, deltaParamsExtended, nMaterialSurfaces,
1550 scatteringMap, geoIdVector);
1551 ACTS_VERBOSE("Updated parameters: " << params.parameters().transpose());
1552
1553 updateGx2fCovarianceParams(fullCovariancePredicted, extendedSystem);
1554 }
1555 ACTS_DEBUG("Finished to evaluate material");
1556 ACTS_VERBOSE(
1557 "Final parameters after material: " << params.parameters().transpose());
1558
1559
1560 ACTS_VERBOSE("Final scattering angles:");
1561 for (const auto& [key, value] : scatteringMap) {
1562 if (!value.materialIsValid()) {
1563 continue;
1564 }
1565 const auto& angles = value.scatteringAngles();
1566 ACTS_VERBOSE(" ( " << angles[eBoundTheta] << " | " << angles[eBoundPhi]
1567 << " )");
1568 }
1569
1570 ACTS_VERBOSE("Final covariance:\n" << fullCovariancePredicted);
1571
1572
1573
1574
1575
1576
1577 if (gx2fOptions.nUpdateMax > 0) {
1578 ACTS_VERBOSE("Propagate with the final covariance.");
1579
1580 params.covariance() = fullCovariancePredicted;
1581
1582
1583 PropagatorOptions propagatorOptions{gx2fOptions.propagatorPlainOptions};
1584 auto& gx2fActor = propagatorOptions.actorList.template get<GX2FActor>();
1585 gx2fActor.inputMeasurements = &inputMeasurements;
1586 gx2fActor.multipleScattering = multipleScattering;
1587 gx2fActor.extensions = gx2fOptions.extensions;
1588 gx2fActor.calibrationContext = &gx2fOptions.calibrationContext.get();
1589 gx2fActor.actorLogger = m_actorLogger.get();
1590 gx2fActor.scatteringMap = &scatteringMap;
1591 gx2fActor.parametersWithHypothesis = ¶ms;
1592
1593 auto propagatorState = m_propagator.makeState(propagatorOptions);
1594
1595 auto propagatorInitResult =
1596 m_propagator.initialize(propagatorState, params);
1597 if (!propagatorInitResult.ok()) {
1598 ACTS_ERROR("Propagation initialization failed: "
1599 << propagatorInitResult.error());
1600 return propagatorInitResult.error();
1601 }
1602
1603 auto& r = propagatorState.template get<Gx2FitterResult<traj_t>>();
1604 r.fittedStates = &trackContainer.trackStateContainer();
1605
1606 auto propagationResult = m_propagator.propagate(propagatorState);
1607
1608
1609 auto result =
1610 m_propagator.makeResult(std::move(propagatorState), propagationResult,
1611 propagatorOptions, false);
1612
1613 if (!result.ok()) {
1614 ACTS_ERROR("Propagation failed: " << result.error());
1615 return result.error();
1616 }
1617
1618 auto& propRes = *result;
1619 GX2FResult gx2fResult = std::move(propRes.template get<GX2FResult>());
1620
1621 if (tipIndex != gx2fResult.lastMeasurementIndex) {
1622 ACTS_INFO("Final fit used unreachable measurements.");
1623 tipIndex = gx2fResult.lastMeasurementIndex;
1624
1625
1626
1627 if (tipIndex == kInvalid) {
1628 ACTS_INFO("Did not find any measurements in final propagation.");
1629 return Experimental::GlobalChiSquareFitterError::
1630 NotEnoughMeasurements;
1631 }
1632 }
1633 }
1634
1635 if (!trackContainer.hasColumn(
1636 Acts::hashString(Gx2fConstants::gx2fnUpdateColumn))) {
1637 trackContainer.template addColumn<std::uint32_t>("Gx2fnUpdateColumn");
1638 }
1639
1640
1641 auto track = trackContainer.makeTrack();
1642 track.tipIndex() = tipIndex;
1643 track.parameters() = params.parameters();
1644 track.covariance() = fullCovariancePredicted;
1645 track.setReferenceSurface(params.referenceSurface().getSharedPtr());
1646
1647 if (trackContainer.hasColumn(
1648 Acts::hashString(Gx2fConstants::gx2fnUpdateColumn))) {
1649 ACTS_DEBUG("Add nUpdate to track");
1650 track.template component<std::uint32_t>("Gx2fnUpdateColumn") =
1651 static_cast<std::uint32_t>(nUpdate);
1652 }
1653
1654
1655 calculateTrackQuantities(track);
1656
1657
1658
1659 track.chi2() = chi2sum;
1660
1661
1662 return track;
1663 }
1664 };
1665
1666 }