Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-07-01 07:52:29

0001 // This file is part of the ACTS project.
0002 //
0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/.
0008 
0009 #pragma once
0010 
0011 #include "Acts/Propagator/EigenStepper.hpp"
0012 
0013 #include "Acts/Definitions/Direction.hpp"
0014 #include "Acts/EventData/TransformationHelpers.hpp"
0015 #include "Acts/Propagator/ConstrainedStep.hpp"
0016 #include "Acts/Propagator/EigenStepperError.hpp"
0017 #include "Acts/Propagator/detail/CovarianceEngine.hpp"
0018 
0019 template <typename E>
0020 Acts::EigenStepper<E>::EigenStepper(
0021     std::shared_ptr<const MagneticFieldProvider> bField)
0022     : m_bField(std::move(bField)) {}
0023 
0024 template <typename E>
0025 auto Acts::EigenStepper<E>::makeState(const Options& options) const -> State {
0026   State state{options, m_bField->makeCache(options.magFieldContext)};
0027   return state;
0028 }
0029 
0030 template <typename E>
0031 void Acts::EigenStepper<E>::initialize(State& state,
0032                                        const BoundTrackParameters& par) const {
0033   initialize(state, par.parameters(), par.covariance(),
0034              par.particleHypothesis(), par.referenceSurface());
0035 }
0036 
0037 template <typename E>
0038 void Acts::EigenStepper<E>::initialize(State& state,
0039                                        const BoundVector& boundParams,
0040                                        const std::optional<BoundMatrix>& cov,
0041                                        ParticleHypothesis particleHypothesis,
0042                                        const Surface& surface) const {
0043   FreeVector freeParams = transformBoundToFreeParameters(
0044       surface, state.options.geoContext, boundParams);
0045 
0046   state.particleHypothesis = particleHypothesis;
0047 
0048   state.pathAccumulated = 0;
0049   state.nSteps = 0;
0050   state.nStepTrials = 0;
0051   state.stepSize = ConstrainedStep();
0052   state.stepSize.setAccuracy(state.options.initialStepSize);
0053   state.stepSize.setUser(state.options.maxStepSize);
0054   state.previousStepSize = 0;
0055   state.statistics = StepperStatistics();
0056 
0057   state.pars = freeParams;
0058 
0059   // Init the jacobian matrix if needed
0060   state.covTransport = cov.has_value();
0061   if (state.covTransport) {
0062     state.cov = *cov;
0063     state.jacToGlobal = surface.boundToFreeJacobian(
0064         state.options.geoContext, freeParams.segment<3>(eFreePos0),
0065         freeParams.segment<3>(eFreeDir0));
0066     state.jacobian = BoundMatrix::Identity();
0067     state.jacTransport = FreeMatrix::Identity();
0068     state.derivative = FreeVector::Zero();
0069   }
0070 }
0071 
0072 template <typename E>
0073 auto Acts::EigenStepper<E>::boundState(
0074     State& state, const Surface& surface, bool transportCov,
0075     const FreeToBoundCorrection& freeToBoundCorrection) const
0076     -> Result<BoundState> {
0077   return detail::boundState(
0078       state.options.geoContext, surface, state.cov, state.jacobian,
0079       state.jacTransport, state.derivative, state.jacToGlobal, state.pars,
0080       state.particleHypothesis, state.covTransport && transportCov,
0081       state.pathAccumulated, freeToBoundCorrection);
0082 }
0083 
0084 template <typename E>
0085 bool Acts::EigenStepper<E>::prepareCurvilinearState(State& state) const {
0086   // test whether the accumulated path has still its initial value.
0087   if (state.pathAccumulated != 0) {
0088     return true;
0089   }
0090 
0091   // if no step was executed the path length derivates have not been
0092   // computed but are needed to compute the curvilinear covariance. The
0093   // derivates are given by k1 for a zero step width.
0094   // First Runge-Kutta point (at current position)
0095   auto& sd = state.stepData;
0096   auto pos = position(state);
0097   auto fieldRes = getField(state, pos);
0098   if (!fieldRes.ok()) {
0099     return false;
0100   }
0101 
0102   sd.B_first = *fieldRes;
0103   if (!state.extension.template k<0>(state, *this, nullptr, sd.k1, sd.B_first,
0104                                      sd.kQoP)) {
0105     return false;
0106   }
0107 
0108   // dr/ds :
0109   state.derivative.template head<3>() =
0110       state.pars.template segment<3>(eFreeDir0);
0111   // d (dr/ds) / ds :
0112   state.derivative.template segment<3>(4) = sd.k1;
0113   // to set dt/ds :
0114   state.extension.finalize(state, *this, nullptr, state.pathAccumulated);
0115   return true;
0116 }
0117 
0118 template <typename E>
0119 auto Acts::EigenStepper<E>::curvilinearState(
0120     State& state, bool transportCov) const -> BoundState {
0121   return detail::curvilinearState(
0122       state.cov, state.jacobian, state.jacTransport, state.derivative,
0123       state.jacToGlobal, state.pars, state.particleHypothesis,
0124       state.covTransport && transportCov, state.pathAccumulated);
0125 }
0126 
0127 template <typename E>
0128 void Acts::EigenStepper<E>::update(State& state, const FreeVector& freeParams,
0129                                    const BoundVector& /*boundParams*/,
0130                                    const Covariance& covariance,
0131                                    const Surface& surface) const {
0132   state.pars = freeParams;
0133   state.cov = covariance;
0134   state.jacToGlobal = surface.boundToFreeJacobian(
0135       state.options.geoContext, freeParams.template segment<3>(eFreePos0),
0136       freeParams.template segment<3>(eFreeDir0));
0137 }
0138 
0139 template <typename E>
0140 void Acts::EigenStepper<E>::update(State& state, const Vector3& uposition,
0141                                    const Vector3& udirection, double qOverP,
0142                                    double time) const {
0143   state.pars.template segment<3>(eFreePos0) = uposition;
0144   state.pars.template segment<3>(eFreeDir0) = udirection;
0145   state.pars[eFreeTime] = time;
0146   state.pars[eFreeQOverP] = qOverP;
0147 }
0148 
0149 template <typename E>
0150 void Acts::EigenStepper<E>::transportCovarianceToCurvilinear(
0151     State& state) const {
0152   detail::transportCovarianceToCurvilinear(state.cov, state.jacobian,
0153                                            state.jacTransport, state.derivative,
0154                                            state.jacToGlobal, direction(state));
0155 }
0156 
0157 template <typename E>
0158 void Acts::EigenStepper<E>::transportCovarianceToBound(
0159     State& state, const Surface& surface,
0160     const FreeToBoundCorrection& freeToBoundCorrection) const {
0161   detail::transportCovarianceToBound(
0162       state.options.geoContext, surface, state.cov, state.jacobian,
0163       state.jacTransport, state.derivative, state.jacToGlobal, state.pars,
0164       freeToBoundCorrection);
0165 }
0166 
0167 template <typename E>
0168 Acts::Result<double> Acts::EigenStepper<E>::step(
0169     State& state, Direction propDir, const IVolumeMaterial* material) const {
0170   // Runge-Kutta integrator state
0171   auto& sd = state.stepData;
0172 
0173   double errorEstimate = 0;
0174   double h2 = 0;
0175   double half_h = 0;
0176 
0177   auto pos = position(state);
0178   auto dir = direction(state);
0179 
0180   // First Runge-Kutta point (at current position)
0181   auto fieldRes = getField(state, pos);
0182   if (!fieldRes.ok()) {
0183     return fieldRes.error();
0184   }
0185   sd.B_first = *fieldRes;
0186   if (!state.extension.template k<0>(state, *this, material, sd.k1, sd.B_first,
0187                                      sd.kQoP)) {
0188     return 0.;
0189   }
0190 
0191   const auto calcStepSizeScaling = [&](const double errorEstimate_) -> double {
0192     // For details about these values see ATL-SOFT-PUB-2009-001
0193     constexpr double lower = 0.25;
0194     constexpr double upper = 4.0;
0195     // This is given by the order of the Runge-Kutta method
0196     constexpr double exponent = 0.25;
0197 
0198     double x = state.options.stepTolerance / errorEstimate_;
0199 
0200     if constexpr (exponent == 0.25) {
0201       // This is 3x faster than std::pow
0202       x = std::sqrt(std::sqrt(x));
0203     } else {
0204       x = std::pow(x, exponent);
0205     }
0206 
0207     return std::clamp(x, lower, upper);
0208   };
0209 
0210   const auto isErrorTolerable = [&](const double errorEstimate_) {
0211     // For details about these values see ATL-SOFT-PUB-2009-001
0212     constexpr double marginFactor = 4.0;
0213 
0214     return errorEstimate_ <= marginFactor * state.options.stepTolerance;
0215   };
0216 
0217   // The following functor starts to perform a Runge-Kutta step of a certain
0218   // size, going up to the point where it can return an estimate of the local
0219   // integration error. The results are stated in the local variables above,
0220   // allowing integration to continue once the error is deemed satisfactory
0221   const auto tryRungeKuttaStep = [&](const double h) -> Result<bool> {
0222     // helpers because bool and std::error_code are ambiguous
0223     constexpr auto success = &Result<bool>::success;
0224     constexpr auto failure = &Result<bool>::failure;
0225 
0226     // State the square and half of the step size
0227     h2 = h * h;
0228     half_h = h * 0.5;
0229 
0230     // Second Runge-Kutta point
0231     const Vector3 pos1 = pos + half_h * dir + h2 * 0.125 * sd.k1;
0232     auto field = getField(state, pos1);
0233     if (!field.ok()) {
0234       return failure(field.error());
0235     }
0236     sd.B_middle = *field;
0237 
0238     if (!state.extension.template k<1>(state, *this, material, sd.k2,
0239                                        sd.B_middle, sd.kQoP, half_h, sd.k1)) {
0240       return success(false);
0241     }
0242 
0243     // Third Runge-Kutta point
0244     if (!state.extension.template k<2>(state, *this, material, sd.k3,
0245                                        sd.B_middle, sd.kQoP, half_h, sd.k2)) {
0246       return success(false);
0247     }
0248 
0249     // Last Runge-Kutta point
0250     const Vector3 pos2 = pos + h * dir + h2 * 0.5 * sd.k3;
0251     field = getField(state, pos2);
0252     if (!field.ok()) {
0253       return failure(field.error());
0254     }
0255     sd.B_last = *field;
0256     if (!state.extension.template k<3>(state, *this, material, sd.k4, sd.B_last,
0257                                        sd.kQoP, h, sd.k3)) {
0258       return success(false);
0259     }
0260 
0261     // Compute and check the local integration error estimate
0262     errorEstimate =
0263         h2 * ((sd.k1 - sd.k2 - sd.k3 + sd.k4).template lpNorm<1>() +
0264               std::abs(sd.kQoP[0] - sd.kQoP[1] - sd.kQoP[2] + sd.kQoP[3]));
0265     // Protect against division by zero
0266     errorEstimate = std::max(1e-20, errorEstimate);
0267 
0268     return success(isErrorTolerable(errorEstimate));
0269   };
0270 
0271   const double initialH = state.stepSize.value() * propDir;
0272   double h = initialH;
0273   std::size_t nStepTrials = 0;
0274   // Select and adjust the appropriate Runge-Kutta step size as given
0275   // ATL-SOFT-PUB-2009-001
0276   while (true) {
0277     ++nStepTrials;
0278     ++state.statistics.nAttemptedSteps;
0279 
0280     auto res = tryRungeKuttaStep(h);
0281     if (!res.ok()) {
0282       return res.error();
0283     }
0284     if (!!res.value()) {
0285       break;
0286     }
0287 
0288     ++state.statistics.nRejectedSteps;
0289 
0290     const double stepSizeScaling = calcStepSizeScaling(errorEstimate);
0291     h *= stepSizeScaling;
0292 
0293     // If step size becomes too small the particle remains at the initial
0294     // place
0295     if (std::abs(h) < std::abs(state.options.stepSizeCutOff)) {
0296       // Not moving due to too low momentum needs an aborter
0297       return EigenStepperError::StepSizeStalled;
0298     }
0299 
0300     // If the parameter is off track too much or given stepSize is not
0301     // appropriate
0302     if (nStepTrials > state.options.maxRungeKuttaStepTrials) {
0303       // Too many trials, have to abort
0304       return EigenStepperError::StepSizeAdjustmentFailed;
0305     }
0306   }
0307 
0308   // When doing error propagation, update the associated Jacobian matrix
0309   if (state.covTransport) {
0310     // using the direction before updated below
0311 
0312     // The step transport matrix in global coordinates
0313     FreeMatrix D;
0314     if (!state.extension.finalize(state, *this, material, h, D)) {
0315       return EigenStepperError::StepInvalid;
0316     }
0317 
0318     // See the documentation of Acts::blockedMult for a description of blocked
0319     // matrix multiplication. However, we can go one further. Let's assume that
0320     // some of these sub-matrices are zero matrices 0₈ and identity matrices
0321     // I₈, namely:
0322     //
0323     // D₁₁ = I₈, J₁₁ = I₈, D₂₁ = 0₈, J₂₁ = 0₈
0324     //
0325     // Which gives:
0326     //
0327     // K₁₁ = I₈  * I₈  + D₁₂ * 0₈  = I₈
0328     // K₁₂ = I₈  * J₁₂ + D₁₂ * J₂₂ = J₁₂ + D₁₂ * J₂₂
0329     // K₂₁ = 0₈  * I₈  + D₂₂ * 0₈  = 0₈
0330     // K₂₂ = 0₈  * J₁₂ + D₂₂ * J₂₂ = D₂₂ * J₂₂
0331     //
0332     // Furthermore, we're constructing K in place of J, and since
0333     // K₁₁ = I₈ = J₁₁ and K₂₁ = 0₈ = D₂₁, we don't actually need to touch those
0334     // sub-matrices at all!
0335     assert((D.topLeftCorner<4, 4>().isIdentity()));
0336     assert((D.bottomLeftCorner<4, 4>().isZero()));
0337     assert((state.jacTransport.template topLeftCorner<4, 4>().isIdentity()));
0338     assert((state.jacTransport.template bottomLeftCorner<4, 4>().isZero()));
0339 
0340     state.jacTransport.template topRightCorner<4, 4>() +=
0341         D.topRightCorner<4, 4>() *
0342         state.jacTransport.template bottomRightCorner<4, 4>();
0343     state.jacTransport.template bottomRightCorner<4, 4>() =
0344         (D.bottomRightCorner<4, 4>() *
0345          state.jacTransport.template bottomRightCorner<4, 4>())
0346             .eval();
0347   } else {
0348     if (!state.extension.finalize(state, *this, material, h)) {
0349       return EigenStepperError::StepInvalid;
0350     }
0351   }
0352 
0353   // Update the track parameters according to the equations of motion
0354   state.pars.template segment<3>(eFreePos0) +=
0355       h * dir + h2 / 6. * (sd.k1 + sd.k2 + sd.k3);
0356   state.pars.template segment<3>(eFreeDir0) +=
0357       h / 6. * (sd.k1 + 2. * (sd.k2 + sd.k3) + sd.k4);
0358   (state.pars.template segment<3>(eFreeDir0)).normalize();
0359 
0360   if (state.covTransport) {
0361     // using the updated direction
0362     state.derivative.template head<3>() =
0363         state.pars.template segment<3>(eFreeDir0);
0364     state.derivative.template segment<3>(4) = sd.k4;
0365   }
0366 
0367   state.pathAccumulated += h;
0368   ++state.nSteps;
0369   state.nStepTrials += nStepTrials;
0370 
0371   ++state.statistics.nSuccessfulSteps;
0372   if (propDir != Direction::fromScalarZeroAsPositive(initialH)) {
0373     ++state.statistics.nReverseSteps;
0374   }
0375   state.statistics.pathLength += h;
0376   state.statistics.absolutePathLength += std::abs(h);
0377 
0378   const double stepSizeScaling = calcStepSizeScaling(errorEstimate);
0379   const double nextAccuracy = std::abs(h * stepSizeScaling);
0380   const double previousAccuracy = std::abs(state.stepSize.accuracy());
0381   const double initialStepLength = std::abs(initialH);
0382   if (nextAccuracy < initialStepLength || nextAccuracy > previousAccuracy) {
0383     state.stepSize.setAccuracy(nextAccuracy);
0384   }
0385 
0386   return h;
0387 }
0388 
0389 template <typename E>
0390 void Acts::EigenStepper<E>::setIdentityJacobian(State& state) const {
0391   state.jacobian = BoundMatrix::Identity();
0392 }