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