Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-06-08 07:47:30

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/Definitions/Algebra.hpp"
0012 #include "Acts/Definitions/Common.hpp"
0013 #include "Acts/Definitions/TrackParametrization.hpp"
0014 #include "Acts/Definitions/Units.hpp"
0015 #include "Acts/EventData/BoundTrackParameters.hpp"
0016 #include "Acts/EventData/TransformationHelpers.hpp"
0017 #include "Acts/EventData/detail/CorrectedTransformationFreeToBound.hpp"
0018 #include "Acts/Geometry/GeometryContext.hpp"
0019 #include "Acts/MagneticField/MagneticFieldProvider.hpp"
0020 #include "Acts/Propagator/ConstrainedStep.hpp"
0021 #include "Acts/Propagator/NavigationTarget.hpp"
0022 #include "Acts/Propagator/StepperOptions.hpp"
0023 #include "Acts/Propagator/StepperStatistics.hpp"
0024 #include "Acts/Propagator/detail/SteppingHelper.hpp"
0025 #include "Acts/Surfaces/Surface.hpp"
0026 #include "Acts/Utilities/Intersection.hpp"
0027 #include "Acts/Utilities/Result.hpp"
0028 
0029 #include <cmath>
0030 
0031 namespace Acts {
0032 
0033 class IVolumeMaterial;
0034 
0035 /// @brief the AtlasStepper implementation for the
0036 ///
0037 /// This is based original stepper code from the ATLAS RungeKuttaPropagator
0038 class AtlasStepper {
0039  public:
0040   /// Type alias for bound track parameters
0041   using BoundParameters = BoundTrackParameters;
0042   /// Type alias for Jacobian matrix
0043   using Jacobian = BoundMatrix;
0044   /// Type alias for covariance matrix
0045   using Covariance = BoundMatrix;
0046   /// Type alias for bound state (parameters, jacobian, path length)
0047   using BoundState = std::tuple<BoundParameters, Jacobian, double>;
0048 
0049   /// Configuration for constructing an AtlasStepper.
0050   struct Config {
0051     /// Magnetic field provider
0052     std::shared_ptr<const MagneticFieldProvider> bField;
0053   };
0054 
0055   /// Stepper options extending plain stepper settings.
0056   struct Options : public StepperPlainOptions {
0057     /// Constructor from context objects
0058     /// @param gctx Geometry context
0059     /// @param mctx Magnetic field context
0060     Options(const GeometryContext& gctx, const MagneticFieldContext& mctx)
0061         : StepperPlainOptions(gctx, mctx) {}
0062 
0063     /// Set plain options
0064     /// @param options Plain stepper options to set
0065     void setPlainOptions(const StepperPlainOptions& options) {
0066       static_cast<StepperPlainOptions&>(*this) = options;
0067     }
0068   };
0069 
0070   /// @brief Nested State struct for the local caching
0071   struct State {
0072     /// Constructor
0073     ///
0074     /// @tparam Type of TrackParameters
0075     ///
0076     /// @param [in] optionsIn The options for the stepper
0077     /// @param [in] fieldCacheIn The magnetic field cache
0078     State(const Options& optionsIn, MagneticFieldProvider::Cache fieldCacheIn)
0079         : options(optionsIn), fieldCache(std::move(fieldCacheIn)) {}
0080 
0081     /// Stepper options for this propagation
0082     Options options;
0083 
0084     /// Particle hypothesis used for energy loss calculations
0085     ParticleHypothesis particleHypothesis = ParticleHypothesis::pion();
0086 
0087     // optimisation that init is not called twice
0088     /// Flag indicating whether the state is ready for stepping
0089     bool state_ready = false;
0090     // configuration
0091     /// Flag indicating whether Jacobian calculation is enabled
0092     bool useJacobian = false;
0093     /// Current step size
0094     double step = 0;
0095     /// Maximum path length for this propagation
0096     double maxPathLength = 0;
0097     /// Flag indicating magnetic field boundary conditions
0098     bool mcondition = false;
0099     /// Flag indicating whether magnetic field gradients are needed
0100     bool needgradient = false;
0101     /// Flag indicating new magnetic field evaluation needed
0102     bool newfield = true;
0103     // internal parameters to be used
0104     /// Current magnetic field vector
0105     Vector3 field = Vector3::Zero();
0106     /// Parameter vector for track state and derivatives
0107     std::array<double, 60> pVector{};
0108 
0109     /// Storage pattern of pVector
0110     ///                   /dL0    /dL1    /dPhi   /the   /dCM   /dT
0111     /// X  ->P[0]  dX /   P[ 8]   P[16]   P[24]   P[32]   P[40]  P[48]
0112     /// Y  ->P[1]  dY /   P[ 9]   P[17]   P[25]   P[33]   P[41]  P[49]
0113     /// Z  ->P[2]  dZ /   P[10]   P[18]   P[26]   P[34]   P[42]  P[50]
0114     /// T  ->P[3]  dT/    P[11]   P[19]   P[27]   P[35]   P[43]  P[51]
0115     /// Ax ->P[4]  dAx/   P[12]   P[20]   P[28]   P[36]   P[44]  P[52]
0116     /// Ay ->P[5]  dAy/   P[13]   P[21]   P[29]   P[37]   P[45]  P[53]
0117     /// Az ->P[6]  dAz/   P[14]   P[22]   P[30]   P[38]   P[46]  P[54]
0118     /// CM ->P[7]  dCM/   P[15]   P[23]   P[31]   P[39]   P[47]  P[55]
0119     /// Cache: P[56] - P[59]
0120 
0121     // result
0122     double parameters[eBoundSize] = {0., 0., 0., 0., 0., 0.};
0123     /// Pointer to external covariance matrix
0124     const Covariance* covariance = nullptr;
0125     /// Local covariance matrix storage
0126     Covariance cov = Covariance::Zero();
0127     /// Flag indicating whether covariance transport is enabled
0128     bool covTransport = false;
0129     /// Jacobian matrix storage for parameter derivatives
0130     double jacobian[eBoundSize * eBoundSize] = {};
0131 
0132     // accumulated path length cache
0133     /// Accumulated path length during propagation
0134     double pathAccumulated = 0.;
0135 
0136     /// Total number of performed steps
0137     std::size_t nSteps = 0;
0138 
0139     /// Totoal number of attempted steps
0140     std::size_t nStepTrials = 0;
0141 
0142     // Adaptive step size of the runge-kutta integration
0143     /// Constrained step size for adaptive integration
0144     ConstrainedStep stepSize;
0145 
0146     // Previous step size for overstep estimation
0147     /// Previous step size for overstep estimation
0148     double previousStepSize = 0.;
0149 
0150     /// It caches the current magnetic field cell and stays (and interpolates)
0151     ///  within as long as this is valid. See step() code for details.
0152     MagneticFieldProvider::Cache fieldCache;
0153 
0154     /// Debug output
0155     /// the string where debug messages are stored (optionally)
0156     bool debug = false;
0157     /// String buffer for debug output messages
0158     std::string debugString = "";
0159     /// buffer & formatting for consistent output
0160     std::size_t debugPfxWidth = 30;
0161     /// Width for debug message formatting
0162     std::size_t debugMsgWidth = 50;
0163 
0164     /// The statistics of the stepper
0165     StepperStatistics statistics;
0166   };
0167 
0168   /// Construct AtlasStepper with magnetic field provider
0169   /// @param bField Shared pointer to magnetic field provider
0170   explicit AtlasStepper(std::shared_ptr<const MagneticFieldProvider> bField)
0171       : m_bField(std::move(bField)) {}
0172 
0173   /// Construct AtlasStepper with configuration
0174   /// @param config Configuration containing magnetic field provider and other settings
0175   explicit AtlasStepper(const Config& config) : m_bField(config.bField) {}
0176 
0177   /// Create stepper state with given options
0178   /// @param options Stepper options including magnetic field context
0179   /// @return Stepper state initialized with options and magnetic field cache
0180   State makeState(const Options& options) const {
0181     State state{options, m_bField->makeCache(options.magFieldContext)};
0182     return state;
0183   }
0184 
0185   /// Initialize stepper state from bound track parameters
0186   /// @param state Stepper state to initialize
0187   /// @param par Bound track parameters containing initial conditions
0188   void initialize(State& state, const BoundParameters& par) const {
0189     initialize(state, par.parameters(), par.covariance(),
0190                par.particleHypothesis(), par.referenceSurface());
0191   }
0192 
0193   /// Initialize stepper state from individual parameters
0194   /// @param state Stepper state to initialize
0195   /// @param boundParams Bound track parameter vector
0196   /// @param cov Optional covariance matrix
0197   /// @param particleHypothesis Particle hypothesis (mass, charge)
0198   /// @param surface Reference surface for bound parameters
0199   void initialize(State& state, const BoundVector& boundParams,
0200                   const std::optional<BoundMatrix>& cov,
0201                   ParticleHypothesis particleHypothesis,
0202                   const Surface& surface) const {
0203     state.particleHypothesis = particleHypothesis;
0204 
0205     state.pathAccumulated = 0;
0206     state.nSteps = 0;
0207     state.nStepTrials = 0;
0208     state.stepSize = ConstrainedStep();
0209     state.stepSize.setAccuracy(state.options.initialStepSize);
0210     state.stepSize.setUser(state.options.maxStepSize);
0211     state.previousStepSize = 0;
0212     state.statistics = StepperStatistics();
0213 
0214     // The rest of this constructor is copy&paste of AtlasStepper::update() -
0215     // this is a nasty but working solution for the stepper state without
0216     // functions
0217 
0218     const auto& Vp = boundParams;
0219 
0220     double Sf = std::sin(Vp[eBoundPhi]);
0221     double Cf = std::cos(Vp[eBoundPhi]);
0222     double Se = std::sin(Vp[eBoundTheta]);
0223     double Ce = std::cos(Vp[eBoundTheta]);
0224 
0225     const Vector3 dir = {Cf * Se, Sf * Se, Ce};
0226     const auto pos = surface.localToGlobal(
0227         state.options.geoContext, boundParams.segment<2>(eBoundLoc0), dir);
0228 
0229     double* pVector = state.pVector.data();
0230 
0231     pVector[0] = pos[ePos0];
0232     pVector[1] = pos[ePos1];
0233     pVector[2] = pos[ePos2];
0234     pVector[3] = boundParams[eBoundTime];
0235     pVector[4] = dir[ePos0];
0236     pVector[5] = dir[ePos1];
0237     pVector[6] = dir[ePos2];
0238     pVector[7] = Vp[eBoundQOverP];
0239 
0240     // @todo: remove magic numbers - is that the charge ?
0241     if (std::abs(pVector[7]) < .000000000000001) {
0242       pVector[7] < 0. ? pVector[7] = -.000000000000001
0243                       : pVector[7] = .000000000000001;
0244     }
0245 
0246     // prepare the jacobian if we have a covariance
0247     state.covTransport = cov.has_value();
0248     if (state.covTransport) {
0249       // copy the covariance matrix
0250       state.covariance = new BoundMatrix(*cov);
0251       state.useJacobian = true;
0252       const auto transform =
0253           surface.referenceFrame(state.options.geoContext, pos, dir);
0254 
0255       pVector[8] = transform(0, eBoundLoc0);
0256       pVector[16] = transform(0, eBoundLoc1);
0257       pVector[24] = 0.;
0258       pVector[32] = 0.;
0259       pVector[40] = 0.;
0260       pVector[48] = 0.;  // dX /
0261 
0262       pVector[9] = transform(1, eBoundLoc0);
0263       pVector[17] = transform(1, eBoundLoc1);
0264       pVector[25] = 0.;
0265       pVector[33] = 0.;
0266       pVector[41] = 0.;
0267       pVector[49] = 0.;  // dY /
0268 
0269       pVector[10] = transform(2, eBoundLoc0);
0270       pVector[18] = transform(2, eBoundLoc1);
0271       pVector[26] = 0.;
0272       pVector[34] = 0.;
0273       pVector[42] = 0.;
0274       pVector[50] = 0.;  // dZ /
0275 
0276       pVector[11] = 0.;
0277       pVector[19] = 0.;
0278       pVector[27] = 0.;
0279       pVector[35] = 0.;
0280       pVector[43] = 0.;
0281       pVector[51] = 1.;  // dT/
0282 
0283       pVector[12] = 0.;
0284       pVector[20] = 0.;
0285       pVector[28] = -Sf * Se;  // - sin(phi) * cos(theta)
0286       pVector[36] = Cf * Ce;   // cos(phi) * cos(theta)
0287       pVector[44] = 0.;
0288       pVector[52] = 0.;  // dAx/
0289 
0290       pVector[13] = 0.;
0291       pVector[21] = 0.;
0292       pVector[29] = Cf * Se;  // cos(phi) * sin(theta)
0293       pVector[37] = Sf * Ce;  // sin(phi) * cos(theta)
0294       pVector[45] = 0.;
0295       pVector[53] = 0.;  // dAy/
0296 
0297       pVector[14] = 0.;
0298       pVector[22] = 0.;
0299       pVector[30] = 0.;
0300       pVector[38] = -Se;  // - sin(theta)
0301       pVector[46] = 0.;
0302       pVector[54] = 0.;  // dAz/
0303 
0304       pVector[15] = 0.;
0305       pVector[23] = 0.;
0306       pVector[31] = 0.;
0307       pVector[39] = 0.;
0308       pVector[47] = 1.;
0309       pVector[55] = 0.;  // dCM/
0310 
0311       pVector[56] = 0.;
0312       pVector[57] = 0.;
0313       pVector[58] = 0.;
0314       pVector[59] = 0.;
0315 
0316       // special treatment for surface types
0317       // the disc needs polar coordinate adaptations
0318       if (surface.type() == Surface::Disc) {
0319         double lCf = std::cos(Vp[1]);
0320         double lSf = std::sin(Vp[1]);
0321         double Ax[3] = {transform(0, 0), transform(1, 0), transform(2, 0)};
0322         double Ay[3] = {transform(0, 1), transform(1, 1), transform(2, 1)};
0323         double d0 = lCf * Ax[0] + lSf * Ay[0];
0324         double d1 = lCf * Ax[1] + lSf * Ay[1];
0325         double d2 = lCf * Ax[2] + lSf * Ay[2];
0326         pVector[8] = d0;
0327         pVector[9] = d1;
0328         pVector[10] = d2;
0329         pVector[16] = Vp[0] * (lCf * Ay[0] - lSf * Ax[0]);
0330         pVector[17] = Vp[0] * (lCf * Ay[1] - lSf * Ax[1]);
0331         pVector[18] = Vp[0] * (lCf * Ay[2] - lSf * Ax[2]);
0332       }
0333       // the line needs components that relate direction change
0334       // with global frame change
0335       if (surface.type() == Surface::Perigee ||
0336           surface.type() == Surface::Straw) {
0337         // sticking to the nomenclature of the original RkPropagator
0338         // - axis pointing along the drift/transverse direction
0339         double B[3] = {transform(0, 0), transform(1, 0), transform(2, 0)};
0340         // - axis along the straw
0341         double A[3] = {transform(0, 1), transform(1, 1), transform(2, 1)};
0342         // - normal vector of the reference frame
0343         double C[3] = {transform(0, 2), transform(1, 2), transform(2, 2)};
0344 
0345         // projection of direction onto normal vector of reference frame
0346         double PC = pVector[4] * C[0] + pVector[5] * C[1] + pVector[6] * C[2];
0347         double Bn = 1. / PC;
0348 
0349         double Bx2 = -A[2] * pVector[29];
0350         double Bx3 = A[1] * pVector[38] - A[2] * pVector[37];
0351 
0352         double By2 = A[2] * pVector[28];
0353         double By3 = A[2] * pVector[36] - A[0] * pVector[38];
0354 
0355         double Bz2 = A[0] * pVector[29] - A[1] * pVector[28];
0356         double Bz3 = A[0] * pVector[37] - A[1] * pVector[36];
0357 
0358         double B2 = B[0] * Bx2 + B[1] * By2 + B[2] * Bz2;
0359         double B3 = B[0] * Bx3 + B[1] * By3 + B[2] * Bz3;
0360 
0361         Bx2 = (Bx2 - B[0] * B2) * Bn;
0362         Bx3 = (Bx3 - B[0] * B3) * Bn;
0363         By2 = (By2 - B[1] * B2) * Bn;
0364         By3 = (By3 - B[1] * B3) * Bn;
0365         Bz2 = (Bz2 - B[2] * B2) * Bn;
0366         Bz3 = (Bz3 - B[2] * B3) * Bn;
0367 
0368         //  /dPhi      |     /dThe       |
0369         pVector[24] = Bx2 * Vp[0];
0370         pVector[32] = Bx3 * Vp[0];  // dX/
0371         pVector[25] = By2 * Vp[0];
0372         pVector[33] = By3 * Vp[0];  // dY/
0373         pVector[26] = Bz2 * Vp[0];
0374         pVector[34] = Bz3 * Vp[0];  // dZ/
0375       }
0376     }
0377 
0378     // now declare the state as ready
0379     state.state_ready = true;
0380   }
0381 
0382   /// Get the field for the stepping
0383   /// It checks first if the access is still within the Cell,
0384   /// and updates the cell if necessary, then it takes the field
0385   /// from the cell
0386   /// @param [in,out] state is the stepper state associated with the track
0387   ///                 the magnetic field cell is used (and potentially updated)
0388   /// @param [in] pos is the field position
0389   /// @return Magnetic field vector at the given position or error
0390   Result<Vector3> getField(State& state, const Vector3& pos) const {
0391     // get the field from the cell
0392     auto res = m_bField->getField(pos, state.fieldCache);
0393     if (res.ok()) {
0394       state.field = *res;
0395     }
0396     return res;
0397   }
0398 
0399   /// Get particle position from stepper state
0400   /// @param state Stepper state containing track information
0401   /// @return 3D position vector
0402   Vector3 position(const State& state) const {
0403     return Vector3(state.pVector[0], state.pVector[1], state.pVector[2]);
0404   }
0405 
0406   /// Get particle direction from stepper state
0407   /// @param state Stepper state containing track information
0408   /// @return 3D unit direction vector
0409   Vector3 direction(const State& state) const {
0410     return Vector3(state.pVector[4], state.pVector[5], state.pVector[6]);
0411   }
0412 
0413   /// Get charge over momentum from stepper state
0414   /// @param state Stepper state containing track information
0415   /// @return Charge over momentum (q/p)
0416   double qOverP(const State& state) const { return state.pVector[7]; }
0417 
0418   /// Absolute momentum accessor
0419   ///
0420   /// @param state [in] The stepping state (thread-local cache)
0421   /// @return Absolute momentum magnitude
0422   double absoluteMomentum(const State& state) const {
0423     return particleHypothesis(state).extractMomentum(qOverP(state));
0424   }
0425 
0426   /// Get 3D momentum vector from stepper state
0427   /// @param state Stepper state containing track information
0428   /// @return 3D momentum vector
0429   Vector3 momentum(const State& state) const {
0430     return absoluteMomentum(state) * direction(state);
0431   }
0432 
0433   /// Charge access
0434   ///
0435   /// @param state [in] The stepping state (thread-local cache)
0436   /// @return Electric charge of the particle
0437   double charge(const State& state) const {
0438     return particleHypothesis(state).extractCharge(qOverP(state));
0439   }
0440 
0441   /// Particle hypothesis
0442   ///
0443   /// @param state [in] The stepping state (thread-local cache)
0444   /// @return Reference to the particle hypothesis used
0445   const ParticleHypothesis& particleHypothesis(const State& state) const {
0446     return state.particleHypothesis;
0447   }
0448 
0449   /// Overstep limit
0450   /// @return Overstep limit for this stepper (negative value)
0451   double overstepLimit(const State& /*state*/) const {
0452     return -m_overstepLimit;
0453   }
0454 
0455   /// Time access
0456   /// @param state The stepping state to get time from
0457   /// @return Current time coordinate
0458   double time(const State& state) const { return state.pVector[3]; }
0459 
0460   /// Update surface status
0461   ///
0462   /// This method intersect the provided surface and update the navigation
0463   /// step estimation accordingly (hence it changes the state). It also
0464   /// returns the status of the intersection to trigger onSurface in case
0465   /// the surface is reached.
0466   ///
0467   /// @param [in,out] state The stepping state (thread-local cache)
0468   /// @param [in] surface The surface provided
0469   /// @param [in] index The surface intersection index
0470   /// @param [in] navDir The navigation direction
0471   /// @param [in] boundaryTolerance The boundary check for this status update
0472   /// @param [in] surfaceTolerance Surface tolerance used for intersection
0473   /// @param [in] stype The step size type to be set
0474   /// @param [in] logger Logger instance to use
0475   /// @return IntersectionStatus indicating how the surface was encountered
0476   IntersectionStatus updateSurfaceStatus(
0477       State& state, const Surface& surface, std::uint8_t index,
0478       Direction navDir, const BoundaryTolerance& boundaryTolerance,
0479       double surfaceTolerance, ConstrainedStep::Type stype,
0480       const Logger& logger = getDummyLogger()) const {
0481     return detail::updateSingleSurfaceStatus<AtlasStepper>(
0482         *this, state, surface, index, navDir, boundaryTolerance,
0483         surfaceTolerance, stype, logger);
0484   }
0485 
0486   /// Update step size
0487   ///
0488   /// It checks the status to the reference surface & updates
0489   /// the step size accordingly
0490   ///
0491   /// @param state [in,out] The stepping state (thread-local cache)
0492   /// @param target [in] The NavigationTarget
0493   /// @param direction [in] The propagation direction
0494   /// @param stype [in] The step size type to be set
0495   void updateStepSize(State& state, const NavigationTarget& target,
0496                       Direction direction, ConstrainedStep::Type stype) const {
0497     static_cast<void>(direction);
0498     double stepSize = target.pathLength();
0499     updateStepSize(state, stepSize, stype);
0500   }
0501 
0502   /// Update step size - explicitly with a double
0503   ///
0504   /// @param [in,out] state The stepping state (thread-local cache)
0505   /// @param [in] stepSize The step size value
0506   /// @param [in] stype The step size type to be set
0507   void updateStepSize(State& state, double stepSize,
0508                       ConstrainedStep::Type stype) const {
0509     state.previousStepSize = state.stepSize.value();
0510     state.stepSize.update(stepSize, stype);
0511   }
0512 
0513   /// Release the Step size
0514   ///
0515   /// @param [in,out] state The stepping state (thread-local cache)
0516   /// @param [in] stype The step size type to be released
0517   void releaseStepSize(State& state, ConstrainedStep::Type stype) const {
0518     state.stepSize.release(stype);
0519   }
0520 
0521   /// Get the step size
0522   ///
0523   /// @param state [in] The stepping state (thread-local cache)
0524   /// @param stype [in] The step size type to be returned
0525   /// @return Current step size for the specified constraint type
0526   double getStepSize(const State& state, ConstrainedStep::Type stype) const {
0527     return state.stepSize.value(stype);
0528   }
0529 
0530   /// Output the Step Size - single component
0531   ///
0532   /// @param [in,out] state The stepping state (thread-local cache)
0533   /// @return String representation of the current step size
0534   std::string outputStepSize(const State& state) const {
0535     return state.stepSize.toString();
0536   }
0537 
0538   /// Create and return the bound state at the current position
0539   ///
0540   ///
0541   /// @param [in] state State that will be presented as @c BoundState
0542   /// @param [in] surface The surface to which we bind the state
0543   /// @param [in] transportCov Flag steering covariance transport
0544   /// @param [in] freeToBoundCorrection Correction for non-linearity effect during transform from free to bound
0545   ///
0546   /// @return A bound state:
0547   ///   - the parameters at the surface
0548   ///   - the stepwise jacobian towards it
0549   ///   - and the path length (from start - for ordering)
0550   Result<BoundState> boundState(
0551       State& state, const Surface& surface, bool transportCov = true,
0552       const FreeToBoundCorrection& freeToBoundCorrection =
0553           FreeToBoundCorrection(false)) const {
0554     // the convert method invalidates the state (in case it's reused)
0555     state.state_ready = false;
0556     // extract state information
0557     Acts::Vector4 pos4;
0558     pos4[ePos0] = state.pVector[0];
0559     pos4[ePos1] = state.pVector[1];
0560     pos4[ePos2] = state.pVector[2];
0561     pos4[eTime] = state.pVector[3];
0562     Acts::Vector3 dir;
0563     dir[eMom0] = state.pVector[4];
0564     dir[eMom1] = state.pVector[5];
0565     dir[eMom2] = state.pVector[6];
0566     const auto qOverP = state.pVector[7];
0567 
0568     // The transport of the covariance
0569     std::optional<Covariance> covOpt = std::nullopt;
0570     if (state.covTransport && transportCov) {
0571       transportCovarianceToBound(state, surface, freeToBoundCorrection);
0572     }
0573     if (state.cov != Covariance::Zero()) {
0574       covOpt = state.cov;
0575     }
0576 
0577     // Fill the end parameters
0578     auto parameters = BoundTrackParameters::create(
0579         state.options.geoContext, surface.getSharedPtr(), pos4, dir, qOverP,
0580         std::move(covOpt), state.particleHypothesis);
0581     if (!parameters.ok()) {
0582       return parameters.error();
0583     }
0584 
0585     Jacobian jacobian(state.jacobian);
0586 
0587     return BoundState(std::move(*parameters), jacobian.transpose(),
0588                       state.pathAccumulated);
0589   }
0590 
0591   /// @brief If necessary fill additional members needed for curvilinearState
0592   ///
0593   /// Compute path length derivatives in case they have not been computed
0594   /// yet, which is the case if no step has been executed yet.
0595   ///
0596   /// @param [in, out] state The stepping state (thread-local cache)
0597   /// @return true if nothing is missing after this call, false otherwise.
0598   bool prepareCurvilinearState(State& state) const {
0599     static_cast<void>(state);
0600     return true;
0601   }
0602 
0603   /// Create and return a curvilinear state at the current position
0604   ///
0605   ///
0606   /// @param [in] state State that will be presented as @c CurvilinearState
0607   /// @param [in] transportCov Flag steering covariance transport
0608   ///
0609   /// @return A curvilinear state:
0610   ///   - the curvilinear parameters at given position
0611   ///   - the stepweise jacobian towards it
0612   ///   - and the path length (from start - for ordering)
0613   BoundState curvilinearState(State& state, bool transportCov = true) const {
0614     // the convert method invalidates the state (in case it's reused)
0615     state.state_ready = false;
0616     // extract state information
0617     Acts::Vector4 pos4;
0618     pos4[ePos0] = state.pVector[0];
0619     pos4[ePos1] = state.pVector[1];
0620     pos4[ePos2] = state.pVector[2];
0621     pos4[eTime] = state.pVector[3];
0622     Acts::Vector3 dir;
0623     dir[eMom0] = state.pVector[4];
0624     dir[eMom1] = state.pVector[5];
0625     dir[eMom2] = state.pVector[6];
0626     const auto qOverP = state.pVector[7];
0627 
0628     std::optional<Covariance> covOpt = std::nullopt;
0629     if (state.covTransport && transportCov) {
0630       transportCovarianceToCurvilinear(state);
0631     }
0632     if (state.cov != Covariance::Zero()) {
0633       covOpt = state.cov;
0634     }
0635 
0636     BoundTrackParameters parameters = BoundTrackParameters::createCurvilinear(
0637         pos4, dir, qOverP, std::move(covOpt), state.particleHypothesis);
0638 
0639     Jacobian jacobian(state.jacobian);
0640 
0641     return BoundState(std::move(parameters), jacobian.transpose(),
0642                       state.pathAccumulated);
0643   }
0644 
0645   /// The state update method
0646   ///
0647   /// @param [in,out] state The stepper state for
0648   /// @param [in] parameters The new free track parameters at start
0649   /// @param [in] boundParams Corresponding bound parameters
0650   /// @param [in] covariance The updated covariance matrix
0651   /// @param [in] surface The surface used to update the pVector
0652   void update(State& state, const FreeVector& parameters,
0653               const BoundVector& boundParams, const Covariance& covariance,
0654               const Surface& surface) const {
0655     Vector3 direction = parameters.template segment<3>(eFreeDir0).normalized();
0656     state.pVector[0] = parameters[eFreePos0];
0657     state.pVector[1] = parameters[eFreePos1];
0658     state.pVector[2] = parameters[eFreePos2];
0659     state.pVector[3] = parameters[eFreeTime];
0660     state.pVector[4] = direction.x();
0661     state.pVector[5] = direction.y();
0662     state.pVector[6] = direction.z();
0663     state.pVector[7] = std::copysign(parameters[eFreeQOverP], state.pVector[7]);
0664 
0665     // @todo: remove magic numbers - is that the charge ?
0666     if (std::abs(state.pVector[7]) < .000000000000001) {
0667       state.pVector[7] < 0. ? state.pVector[7] = -.000000000000001
0668                             : state.pVector[7] = .000000000000001;
0669     }
0670 
0671     // prepare the jacobian if we have a covariance
0672     // copy the covariance matrix
0673 
0674     Vector3 pos(state.pVector[0], state.pVector[1], state.pVector[2]);
0675     Vector3 mom(state.pVector[4], state.pVector[5], state.pVector[6]);
0676 
0677     double Sf = std::sin(boundParams[eBoundPhi]);
0678     double Cf = std::cos(boundParams[eBoundPhi]);
0679     double Se = std::sin(boundParams[eBoundTheta]);
0680     double Ce = std::cos(boundParams[eBoundTheta]);
0681 
0682     const auto transform =
0683         surface.referenceFrame(state.options.geoContext, pos, mom);
0684 
0685     state.pVector[8] = transform(0, eBoundLoc0);
0686     state.pVector[16] = transform(0, eBoundLoc1);
0687     state.pVector[24] = 0.;
0688     state.pVector[32] = 0.;
0689     state.pVector[40] = 0.;
0690     state.pVector[48] = 0.;  // dX /
0691 
0692     state.pVector[9] = transform(1, eBoundLoc0);
0693     state.pVector[17] = transform(1, eBoundLoc1);
0694     state.pVector[25] = 0.;
0695     state.pVector[33] = 0.;
0696     state.pVector[41] = 0.;
0697     state.pVector[49] = 0.;  // dY /
0698 
0699     state.pVector[10] = transform(2, eBoundLoc0);
0700     state.pVector[18] = transform(2, eBoundLoc1);
0701     state.pVector[26] = 0.;
0702     state.pVector[34] = 0.;
0703     state.pVector[42] = 0.;
0704     state.pVector[50] = 0.;  // dZ /
0705 
0706     state.pVector[11] = 0.;
0707     state.pVector[19] = 0.;
0708     state.pVector[27] = 0.;
0709     state.pVector[35] = 0.;
0710     state.pVector[43] = 0.;
0711     state.pVector[51] = 1.;  // dT/
0712 
0713     state.pVector[12] = 0.;
0714     state.pVector[20] = 0.;
0715     state.pVector[28] = -Sf * Se;  // - sin(phi) * cos(theta)
0716     state.pVector[36] = Cf * Ce;   // cos(phi) * cos(theta)
0717     state.pVector[44] = 0.;
0718     state.pVector[52] = 0.;  // dAx/
0719 
0720     state.pVector[13] = 0.;
0721     state.pVector[21] = 0.;
0722     state.pVector[29] = Cf * Se;  // cos(phi) * sin(theta)
0723     state.pVector[37] = Sf * Ce;  // sin(phi) * cos(theta)
0724     state.pVector[45] = 0.;
0725     state.pVector[53] = 0.;  // dAy/
0726 
0727     state.pVector[14] = 0.;
0728     state.pVector[22] = 0.;
0729     state.pVector[30] = 0.;
0730     state.pVector[38] = -Se;  // - sin(theta)
0731     state.pVector[46] = 0.;
0732     state.pVector[54] = 0.;  // dAz/
0733 
0734     state.pVector[15] = 0.;
0735     state.pVector[23] = 0.;
0736     state.pVector[31] = 0.;
0737     state.pVector[39] = 0.;
0738     state.pVector[47] = 1.;
0739     state.pVector[55] = 0.;  // dCM/
0740 
0741     state.pVector[56] = 0.;
0742     state.pVector[57] = 0.;
0743     state.pVector[58] = 0.;
0744     state.pVector[59] = 0.;
0745 
0746     // special treatment for surface types
0747     // the disc needs polar coordinate adaptations
0748     if (surface.type() == Surface::Disc) {
0749       double lCf = std::cos(boundParams[eBoundLoc1]);
0750       double lSf = std::sin(boundParams[eBoundLoc1]);
0751       double Ax[3] = {transform(0, 0), transform(1, 0), transform(2, 0)};
0752       double Ay[3] = {transform(0, 1), transform(1, 1), transform(2, 1)};
0753       double d0 = lCf * Ax[0] + lSf * Ay[0];
0754       double d1 = lCf * Ax[1] + lSf * Ay[1];
0755       double d2 = lCf * Ax[2] + lSf * Ay[2];
0756       state.pVector[8] = d0;
0757       state.pVector[9] = d1;
0758       state.pVector[10] = d2;
0759       state.pVector[16] = boundParams[eBoundLoc0] * (lCf * Ay[0] - lSf * Ax[0]);
0760       state.pVector[17] = boundParams[eBoundLoc0] * (lCf * Ay[1] - lSf * Ax[1]);
0761       state.pVector[18] = boundParams[eBoundLoc0] * (lCf * Ay[2] - lSf * Ax[2]);
0762     }
0763     // the line needs components that relate direction change
0764     // with global frame change
0765     if (surface.type() == Surface::Perigee ||
0766         surface.type() == Surface::Straw) {
0767       // sticking to the nomenclature of the original RkPropagator
0768       // - axis pointing along the drift/transverse direction
0769       double B[3] = {transform(0, 0), transform(1, 0), transform(2, 0)};
0770       // - axis along the straw
0771       double A[3] = {transform(0, 1), transform(1, 1), transform(2, 1)};
0772       // - normal vector of the reference frame
0773       double C[3] = {transform(0, 2), transform(1, 2), transform(2, 2)};
0774 
0775       // projection of direction onto normal vector of reference frame
0776       double PC = state.pVector[4] * C[0] + state.pVector[5] * C[1] +
0777                   state.pVector[6] * C[2];
0778       double Bn = 1. / PC;
0779 
0780       double Bx2 = -A[2] * state.pVector[29];
0781       double Bx3 = A[1] * state.pVector[38] - A[2] * state.pVector[37];
0782 
0783       double By2 = A[2] * state.pVector[28];
0784       double By3 = A[2] * state.pVector[36] - A[0] * state.pVector[38];
0785 
0786       double Bz2 = A[0] * state.pVector[29] - A[1] * state.pVector[28];
0787       double Bz3 = A[0] * state.pVector[37] - A[1] * state.pVector[36];
0788 
0789       double B2 = B[0] * Bx2 + B[1] * By2 + B[2] * Bz2;
0790       double B3 = B[0] * Bx3 + B[1] * By3 + B[2] * Bz3;
0791 
0792       Bx2 = (Bx2 - B[0] * B2) * Bn;
0793       Bx3 = (Bx3 - B[0] * B3) * Bn;
0794       By2 = (By2 - B[1] * B2) * Bn;
0795       By3 = (By3 - B[1] * B3) * Bn;
0796       Bz2 = (Bz2 - B[2] * B2) * Bn;
0797       Bz3 = (Bz3 - B[2] * B3) * Bn;
0798 
0799       //  /dPhi      |     /dThe       |
0800       state.pVector[24] = Bx2 * boundParams[eBoundLoc0];
0801       state.pVector[32] = Bx3 * boundParams[eBoundLoc0];  // dX/
0802       state.pVector[25] = By2 * boundParams[eBoundLoc0];
0803       state.pVector[33] = By3 * boundParams[eBoundLoc0];  // dY/
0804       state.pVector[26] = Bz2 * boundParams[eBoundLoc0];
0805       state.pVector[34] = Bz3 * boundParams[eBoundLoc0];  // dZ/
0806     }
0807 
0808     state.covariance = new BoundMatrix(covariance);
0809     state.covTransport = true;
0810     state.useJacobian = true;
0811 
0812     // declare the state as ready
0813     state.state_ready = true;
0814   }
0815 
0816   /// Method to update momentum, direction and p
0817   ///
0818   /// @param state The state object
0819   /// @param uposition the updated position
0820   /// @param udirection the updated direction
0821   /// @param qop the updated momentum value
0822   /// @param time the update time
0823   void update(State& state, const Vector3& uposition, const Vector3& udirection,
0824               double qop, double time) const {
0825     // update the vector
0826     state.pVector[0] = uposition[0];
0827     state.pVector[1] = uposition[1];
0828     state.pVector[2] = uposition[2];
0829     state.pVector[3] = time;
0830     state.pVector[4] = udirection[0];
0831     state.pVector[5] = udirection[1];
0832     state.pVector[6] = udirection[2];
0833     state.pVector[7] = qop;
0834   }
0835 
0836   /// Method for on-demand transport of the covariance
0837   /// to a new curvilinear frame at current  position,
0838   /// or direction of the state
0839   ///
0840   /// @param [in,out] state State of the stepper
0841   void transportCovarianceToCurvilinear(State& state) const {
0842     double P[60];
0843     for (unsigned int i = 0; i < 60; ++i) {
0844       P[i] = state.pVector[i];
0845     }
0846 
0847     double p = 1. / P[7];
0848     P[40] *= p;
0849     P[41] *= p;
0850     P[42] *= p;
0851     P[44] *= p;
0852     P[45] *= p;
0853     P[46] *= p;
0854 
0855     double An = std::sqrt(P[4] * P[4] + P[5] * P[5]);
0856     double Ax[3];
0857     if (An != 0.) {
0858       Ax[0] = -P[5] / An;
0859       Ax[1] = P[4] / An;
0860       Ax[2] = 0.;
0861     } else {
0862       Ax[0] = 1.;
0863       Ax[1] = 0.;
0864       Ax[2] = 0.;
0865     }
0866 
0867     double Ay[3] = {-Ax[1] * P[6], Ax[0] * P[6], An};
0868     double S[3] = {P[4], P[5], P[6]};
0869 
0870     double A = P[4] * S[0] + P[5] * S[1] + P[6] * S[2];
0871     if (A != 0.) {
0872       A = 1. / A;
0873     }
0874     S[0] *= A;
0875     S[1] *= A;
0876     S[2] *= A;
0877 
0878     double s0 = P[8] * S[0] + P[9] * S[1] + P[10] * S[2];
0879     double s1 = P[16] * S[0] + P[17] * S[1] + P[18] * S[2];
0880     double s2 = P[24] * S[0] + P[25] * S[1] + P[26] * S[2];
0881     double s3 = P[32] * S[0] + P[33] * S[1] + P[34] * S[2];
0882     double s4 = P[40] * S[0] + P[41] * S[1] + P[42] * S[2];
0883 
0884     P[8] -= (s0 * P[4]);
0885     P[9] -= (s0 * P[5]);
0886     P[10] -= (s0 * P[6]);
0887     P[11] -= (s0 * P[59]);
0888     P[12] -= (s0 * P[56]);
0889     P[13] -= (s0 * P[57]);
0890     P[14] -= (s0 * P[58]);
0891     P[16] -= (s1 * P[4]);
0892     P[17] -= (s1 * P[5]);
0893     P[18] -= (s1 * P[6]);
0894     P[19] -= (s1 * P[59]);
0895     P[20] -= (s1 * P[56]);
0896     P[21] -= (s1 * P[57]);
0897     P[22] -= (s1 * P[58]);
0898     P[24] -= (s2 * P[4]);
0899     P[25] -= (s2 * P[5]);
0900     P[26] -= (s2 * P[6]);
0901     P[27] -= (s2 * P[59]);
0902     P[28] -= (s2 * P[56]);
0903     P[29] -= (s2 * P[57]);
0904     P[30] -= (s2 * P[58]);
0905     P[32] -= (s3 * P[4]);
0906     P[33] -= (s3 * P[5]);
0907     P[34] -= (s3 * P[6]);
0908     P[35] -= (s3 * P[59]);
0909     P[36] -= (s3 * P[56]);
0910     P[37] -= (s3 * P[57]);
0911     P[38] -= (s3 * P[58]);
0912     P[40] -= (s4 * P[4]);
0913     P[41] -= (s4 * P[5]);
0914     P[42] -= (s4 * P[6]);
0915     P[43] -= (s4 * P[59]);
0916     P[44] -= (s4 * P[56]);
0917     P[45] -= (s4 * P[57]);
0918     P[46] -= (s4 * P[58]);
0919 
0920     double P3 = 0, P4 = 0, C = P[4] * P[4] + P[5] * P[5];
0921     if (C > 1.e-20) {
0922       C = 1. / C;
0923       P3 = P[4] * C;
0924       P4 = P[5] * C;
0925       C = -sqrt(C);
0926     } else {
0927       C = -1.e10;
0928       P3 = 1.;
0929       P4 = 0.;
0930     }
0931 
0932     // Jacobian production
0933     //
0934     state.jacobian[0] = Ax[0] * P[8] + Ax[1] * P[9];    // dL0/dL0
0935     state.jacobian[1] = Ax[0] * P[16] + Ax[1] * P[17];  // dL0/dL1
0936     state.jacobian[2] = Ax[0] * P[24] + Ax[1] * P[25];  // dL0/dPhi
0937     state.jacobian[3] = Ax[0] * P[32] + Ax[1] * P[33];  // dL0/dThe
0938     state.jacobian[4] = Ax[0] * P[40] + Ax[1] * P[41];  // dL0/dCM
0939     state.jacobian[5] = 0.;                             // dL0/dT
0940 
0941     state.jacobian[6] = Ay[0] * P[8] + Ay[1] * P[9] + Ay[2] * P[10];  // dL1/dL0
0942     state.jacobian[7] =
0943         Ay[0] * P[16] + Ay[1] * P[17] + Ay[2] * P[18];  // dL1/dL1
0944     state.jacobian[8] =
0945         Ay[0] * P[24] + Ay[1] * P[25] + Ay[2] * P[26];  // dL1/dPhi
0946     state.jacobian[9] =
0947         Ay[0] * P[32] + Ay[1] * P[33] + Ay[2] * P[34];  // dL1/dThe
0948     state.jacobian[10] =
0949         Ay[0] * P[40] + Ay[1] * P[41] + Ay[2] * P[42];  // dL1/dCM
0950     state.jacobian[11] = 0.;                            // dL1/dT
0951 
0952     state.jacobian[12] = P3 * P[13] - P4 * P[12];  // dPhi/dL0
0953     state.jacobian[13] = P3 * P[21] - P4 * P[20];  // dPhi/dL1
0954     state.jacobian[14] = P3 * P[29] - P4 * P[28];  // dPhi/dPhi
0955     state.jacobian[15] = P3 * P[37] - P4 * P[36];  // dPhi/dThe
0956     state.jacobian[16] = P3 * P[45] - P4 * P[44];  // dPhi/dCM
0957     state.jacobian[17] = 0.;                       // dPhi/dT
0958 
0959     state.jacobian[18] = C * P[14];  // dThe/dL0
0960     state.jacobian[19] = C * P[22];  // dThe/dL1
0961     state.jacobian[20] = C * P[30];  // dThe/dPhi
0962     state.jacobian[21] = C * P[38];  // dThe/dThe
0963     state.jacobian[22] = C * P[46];  // dThe/dCM
0964     state.jacobian[23] = 0.;         // dThe/dT
0965 
0966     state.jacobian[24] = 0.;     // dCM /dL0
0967     state.jacobian[25] = 0.;     // dCM /dL1
0968     state.jacobian[26] = 0.;     // dCM /dPhi
0969     state.jacobian[27] = 0.;     // dCM /dTheta
0970     state.jacobian[28] = P[47];  // dCM /dCM
0971     state.jacobian[29] = 0.;     // dCM/dT
0972 
0973     state.jacobian[30] = P[11];  // dT/dL0
0974     state.jacobian[31] = P[19];  // dT/dL1
0975     state.jacobian[32] = P[27];  // dT/dPhi
0976     state.jacobian[33] = P[35];  // dT/dThe
0977     state.jacobian[34] = P[43];  // dT/dCM
0978     state.jacobian[35] = P[51];  // dT/dT
0979 
0980     Eigen::Map<Eigen::Matrix<double, eBoundSize, eBoundSize, Eigen::RowMajor>>
0981         J(state.jacobian);
0982     state.cov = J * (*state.covariance) * J.transpose();
0983   }
0984 
0985   /// Method for on-demand transport of the covariance
0986   /// to a new curvilinear frame at current position,
0987   /// or direction of the state
0988   ///
0989   /// @param [in,out] state State of the stepper
0990   /// @param [in] surface is the surface to which the covariance is forwarded to
0991   void transportCovarianceToBound(
0992       State& state, const Surface& surface,
0993       const FreeToBoundCorrection& /*freeToBoundCorrection*/ =
0994           FreeToBoundCorrection(false)) const {
0995     Acts::Vector3 gp(state.pVector[0], state.pVector[1], state.pVector[2]);
0996     Acts::Vector3 mom(state.pVector[4], state.pVector[5], state.pVector[6]);
0997 
0998     double P[60];
0999     for (unsigned int i = 0; i < 60; ++i) {
1000       P[i] = state.pVector[i];
1001     }
1002 
1003     mom /= std::abs(state.pVector[7]);
1004 
1005     double p = 1. / state.pVector[7];
1006     P[40] *= p;
1007     P[41] *= p;
1008     P[42] *= p;
1009     P[44] *= p;
1010     P[45] *= p;
1011     P[46] *= p;
1012 
1013     const auto fFrame =
1014         surface.referenceFrame(state.options.geoContext, gp, mom);
1015 
1016     double Ax[3] = {fFrame(0, 0), fFrame(1, 0), fFrame(2, 0)};
1017     double Ay[3] = {fFrame(0, 1), fFrame(1, 1), fFrame(2, 1)};
1018     double S[3] = {fFrame(0, 2), fFrame(1, 2), fFrame(2, 2)};
1019 
1020     // this is the projection of direction onto the local normal vector
1021     double A = P[4] * S[0] + P[5] * S[1] + P[6] * S[2];
1022 
1023     if (A != 0.) {
1024       A = 1. / A;
1025     }
1026 
1027     S[0] *= A;
1028     S[1] *= A;
1029     S[2] *= A;
1030 
1031     double s0 = P[8] * S[0] + P[9] * S[1] + P[10] * S[2];
1032     double s1 = P[16] * S[0] + P[17] * S[1] + P[18] * S[2];
1033     double s2 = P[24] * S[0] + P[25] * S[1] + P[26] * S[2];
1034     double s3 = P[32] * S[0] + P[33] * S[1] + P[34] * S[2];
1035     double s4 = P[40] * S[0] + P[41] * S[1] + P[42] * S[2];
1036 
1037     // in case of line-type surfaces - we need to take into account that
1038     // the reference frame changes with variations of all local
1039     // parameters
1040     if (surface.type() == Surface::Straw ||
1041         surface.type() == Surface::Perigee) {
1042       // vector from position to center
1043       double x = P[0] - surface.center(state.options.geoContext).x();
1044       double y = P[1] - surface.center(state.options.geoContext).y();
1045       double z = P[2] - surface.center(state.options.geoContext).z();
1046 
1047       // this is the projection of the direction onto the local y axis
1048       double d = P[4] * Ay[0] + P[5] * Ay[1] + P[6] * Ay[2];
1049 
1050       // this is cos(beta)
1051       double a = (1. - d) * (1. + d);
1052       if (a != 0.) {
1053         a = 1. / a;  // i.e. 1./(1-d^2)
1054       }
1055 
1056       // that's the modified norm vector
1057       double X = d * Ay[0] - P[4];  //
1058       double Y = d * Ay[1] - P[5];  //
1059       double Z = d * Ay[2] - P[6];  //
1060 
1061       // d0 to d1
1062       double d0 = P[12] * Ay[0] + P[13] * Ay[1] + P[14] * Ay[2];
1063       double d1 = P[20] * Ay[0] + P[21] * Ay[1] + P[22] * Ay[2];
1064       double d2 = P[28] * Ay[0] + P[29] * Ay[1] + P[30] * Ay[2];
1065       double d3 = P[36] * Ay[0] + P[37] * Ay[1] + P[38] * Ay[2];
1066       double d4 = P[44] * Ay[0] + P[45] * Ay[1] + P[46] * Ay[2];
1067 
1068       s0 = (((P[8] * X + P[9] * Y + P[10] * Z) + x * (d0 * Ay[0] - P[12])) +
1069             (y * (d0 * Ay[1] - P[13]) + z * (d0 * Ay[2] - P[14]))) *
1070            (-a);
1071 
1072       s1 = (((P[16] * X + P[17] * Y + P[18] * Z) + x * (d1 * Ay[0] - P[20])) +
1073             (y * (d1 * Ay[1] - P[21]) + z * (d1 * Ay[2] - P[22]))) *
1074            (-a);
1075       s2 = (((P[24] * X + P[25] * Y + P[26] * Z) + x * (d2 * Ay[0] - P[28])) +
1076             (y * (d2 * Ay[1] - P[29]) + z * (d2 * Ay[2] - P[30]))) *
1077            (-a);
1078       s3 = (((P[32] * X + P[33] * Y + P[34] * Z) + x * (d3 * Ay[0] - P[36])) +
1079             (y * (d3 * Ay[1] - P[37]) + z * (d3 * Ay[2] - P[38]))) *
1080            (-a);
1081       s4 = (((P[40] * X + P[41] * Y + P[42] * Z) + x * (d4 * Ay[0] - P[44])) +
1082             (y * (d4 * Ay[1] - P[45]) + z * (d4 * Ay[2] - P[46]))) *
1083            (-a);
1084     }
1085 
1086     P[8] -= (s0 * P[4]);
1087     P[9] -= (s0 * P[5]);
1088     P[10] -= (s0 * P[6]);
1089     P[11] -= (s0 * P[59]);
1090     P[12] -= (s0 * P[56]);
1091     P[13] -= (s0 * P[57]);
1092     P[14] -= (s0 * P[58]);
1093 
1094     P[16] -= (s1 * P[4]);
1095     P[17] -= (s1 * P[5]);
1096     P[18] -= (s1 * P[6]);
1097     P[19] -= (s1 * P[59]);
1098     P[20] -= (s1 * P[56]);
1099     P[21] -= (s1 * P[57]);
1100     P[22] -= (s1 * P[58]);
1101 
1102     P[24] -= (s2 * P[4]);
1103     P[25] -= (s2 * P[5]);
1104     P[26] -= (s2 * P[6]);
1105     P[27] -= (s2 * P[59]);
1106     P[28] -= (s2 * P[56]);
1107     P[29] -= (s2 * P[57]);
1108     P[30] -= (s2 * P[58]);
1109 
1110     P[32] -= (s3 * P[4]);
1111     P[33] -= (s3 * P[5]);
1112     P[34] -= (s3 * P[6]);
1113     P[35] -= (s3 * P[59]);
1114     P[36] -= (s3 * P[56]);
1115     P[37] -= (s3 * P[57]);
1116     P[38] -= (s3 * P[58]);
1117 
1118     P[40] -= (s4 * P[4]);
1119     P[41] -= (s4 * P[5]);
1120     P[42] -= (s4 * P[6]);
1121     P[43] -= (s4 * P[59]);
1122     P[44] -= (s4 * P[56]);
1123     P[45] -= (s4 * P[57]);
1124     P[46] -= (s4 * P[58]);
1125 
1126     double P3 = 0, P4 = 0, C = P[4] * P[4] + P[5] * P[5];
1127     if (C > 1.e-20) {
1128       C = 1. / C;
1129       P3 = P[4] * C;
1130       P4 = P[5] * C;
1131       C = -sqrt(C);
1132     } else {
1133       C = -1.e10;
1134       P3 = 1.;
1135       P4 = 0.;
1136     }
1137 
1138     double MA[3] = {Ax[0], Ax[1], Ax[2]};
1139     double MB[3] = {Ay[0], Ay[1], Ay[2]};
1140     // Jacobian production of transport and to_local
1141     if (surface.type() == Surface::Disc) {
1142       // the vector from the disc surface to the p
1143       const auto& sfc = surface.center(state.options.geoContext);
1144       double d[3] = {P[0] - sfc(0), P[1] - sfc(1), P[2] - sfc(2)};
1145       // this needs the transformation to polar coordinates
1146       double RC = d[0] * Ax[0] + d[1] * Ax[1] + d[2] * Ax[2];
1147       double RS = d[0] * Ay[0] + d[1] * Ay[1] + d[2] * Ay[2];
1148       double R2 = RC * RC + RS * RS;
1149 
1150       // inverse radius
1151       double Ri = 1. / std::sqrt(R2);
1152       MA[0] = (RC * Ax[0] + RS * Ay[0]) * Ri;
1153       MA[1] = (RC * Ax[1] + RS * Ay[1]) * Ri;
1154       MA[2] = (RC * Ax[2] + RS * Ay[2]) * Ri;
1155       MB[0] = (RC * Ay[0] - RS * Ax[0]) * (Ri = 1. / R2);
1156       MB[1] = (RC * Ay[1] - RS * Ax[1]) * Ri;
1157       MB[2] = (RC * Ay[2] - RS * Ax[2]) * Ri;
1158     }
1159 
1160     state.jacobian[0] = MA[0] * P[8] + MA[1] * P[9] + MA[2] * P[10];  // dL0/dL0
1161     state.jacobian[1] =
1162         MA[0] * P[16] + MA[1] * P[17] + MA[2] * P[18];  // dL0/dL1
1163     state.jacobian[2] =
1164         MA[0] * P[24] + MA[1] * P[25] + MA[2] * P[26];  // dL0/dPhi
1165     state.jacobian[3] =
1166         MA[0] * P[32] + MA[1] * P[33] + MA[2] * P[34];  // dL0/dThe
1167     state.jacobian[4] =
1168         MA[0] * P[40] + MA[1] * P[41] + MA[2] * P[42];  // dL0/dCM
1169     state.jacobian[5] = 0.;                             // dL0/dT
1170 
1171     state.jacobian[6] = MB[0] * P[8] + MB[1] * P[9] + MB[2] * P[10];  // dL1/dL0
1172     state.jacobian[7] =
1173         MB[0] * P[16] + MB[1] * P[17] + MB[2] * P[18];  // dL1/dL1
1174     state.jacobian[8] =
1175         MB[0] * P[24] + MB[1] * P[25] + MB[2] * P[26];  // dL1/dPhi
1176     state.jacobian[9] =
1177         MB[0] * P[32] + MB[1] * P[33] + MB[2] * P[34];  // dL1/dThe
1178     state.jacobian[10] =
1179         MB[0] * P[40] + MB[1] * P[41] + MB[2] * P[42];  // dL1/dCM
1180     state.jacobian[11] = 0.;                            // dL1/dT
1181 
1182     state.jacobian[12] = P3 * P[13] - P4 * P[12];  // dPhi/dL0
1183     state.jacobian[13] = P3 * P[21] - P4 * P[20];  // dPhi/dL1
1184     state.jacobian[14] = P3 * P[29] - P4 * P[28];  // dPhi/dPhi
1185     state.jacobian[15] = P3 * P[37] - P4 * P[36];  // dPhi/dThe
1186     state.jacobian[16] = P3 * P[45] - P4 * P[44];  // dPhi/dCM
1187     state.jacobian[17] = 0.;                       // dPhi/dT
1188 
1189     state.jacobian[18] = C * P[14];  // dThe/dL0
1190     state.jacobian[19] = C * P[22];  // dThe/dL1
1191     state.jacobian[20] = C * P[30];  // dThe/dPhi
1192     state.jacobian[21] = C * P[38];  // dThe/dThe
1193     state.jacobian[22] = C * P[46];  // dThe/dCM
1194     state.jacobian[23] = 0.;         // dThe/dT
1195 
1196     state.jacobian[24] = 0.;     // dCM /dL0
1197     state.jacobian[25] = 0.;     // dCM /dL1
1198     state.jacobian[26] = 0.;     // dCM /dPhi
1199     state.jacobian[27] = 0.;     // dCM /dTheta
1200     state.jacobian[28] = P[47];  // dCM /dCM
1201     state.jacobian[29] = 0.;     // dCM/dT
1202 
1203     state.jacobian[30] = P[11];  // dT/dL0
1204     state.jacobian[31] = P[19];  // dT/dL1
1205     state.jacobian[32] = P[27];  // dT/dPhi
1206     state.jacobian[33] = P[35];  // dT/dThe
1207     state.jacobian[34] = P[43];  // dT/dCM
1208     state.jacobian[35] = P[51];  // dT/dT
1209 
1210     Eigen::Map<Eigen::Matrix<double, eBoundSize, eBoundSize, Eigen::RowMajor>>
1211         J(state.jacobian);
1212     state.cov = J * (*state.covariance) * J.transpose();
1213   }
1214 
1215   /// Perform the actual step on the state
1216   ///
1217   /// @param [in,out] state State of the stepper
1218   /// @param propDir is the direction of propagation
1219   /// @param material is the optional volume material we are stepping through.
1220   //         This is simply ignored if `nullptr`.
1221   /// @return the result of the step
1222   ///
1223   /// @note The state contains the desired step size. It can be negative during
1224   ///       backwards track propagation, and since we're using an adaptive
1225   ///       algorithm, it can be modified by the stepper class during
1226   ///       propagation.
1227   Result<double> step(State& state, Direction propDir,
1228                       const IVolumeMaterial* material) const {
1229     static_cast<void>(material);
1230 
1231     // we use h for keeping the nominclature with the original atlas code
1232     auto h = state.stepSize.value() * propDir;
1233     bool Jac = state.useJacobian;
1234 
1235     double* R = &(state.pVector[0]);  // Coordinates
1236     double* A = &(state.pVector[4]);  // Directions
1237     double* sA = &(state.pVector[56]);
1238     // Invert mometum/2.
1239     double Pi = 0.5 * state.pVector[7];
1240     //    double dltm = 0.0002 * .03;
1241     Vector3 f0, f;
1242 
1243     // if new field is required get it
1244     if (state.newfield) {
1245       const Vector3 pos(R[0], R[1], R[2]);
1246       // This is sd.B_first in EigenStepper
1247       auto fRes = getField(state, pos);
1248       if (!fRes.ok()) {
1249         return fRes.error();
1250       }
1251       f0 = *fRes;
1252     } else {
1253       f0 = state.field;
1254     }
1255 
1256     bool Helix = false;
1257     // if (std::abs(S) < m_cfg.helixStep) Helix = true;
1258 
1259     const auto calcStepSizeScaling =
1260         [&](const double errorEstimate_) -> double {
1261       // For details about these values see ATL-SOFT-PUB-2009-001
1262       constexpr double lower = 0.25;
1263       constexpr double upper = 4.0;
1264       // This is given by the order of the Runge-Kutta method
1265       constexpr double exponent = 0.25;
1266 
1267       double x = state.options.stepTolerance / errorEstimate_;
1268 
1269       if constexpr (exponent == 0.25) {
1270         // This is 3x faster than std::pow
1271         x = std::sqrt(std::sqrt(x));
1272       } else {
1273         x = std::pow(x, exponent);
1274       }
1275 
1276       return std::clamp(x, lower, upper);
1277     };
1278 
1279     const auto isErrorTolerable = [&](const double errorEstimate_) {
1280       // For details about these values see ATL-SOFT-PUB-2009-001
1281       constexpr double marginFactor = 4.0;
1282 
1283       return errorEstimate_ <= marginFactor * state.options.stepTolerance;
1284     };
1285 
1286     double EST = 0;
1287     double initialH = h;
1288 
1289     std::size_t nStepTrials = 0;
1290     while (h != 0.) {
1291       nStepTrials++;
1292 
1293       // PS2 is h/(2*momentum) in EigenStepper
1294       double S3 = (1. / 3.) * h, S4 = .25 * h, PS2 = Pi * h;
1295 
1296       // First point
1297       //
1298       // H0 is (h/(2*momentum) * sd.B_first) in EigenStepper
1299       double H0[3] = {f0[0] * PS2, f0[1] * PS2, f0[2] * PS2};
1300       // { A0, B0, C0 } is (h/2 * sd.k1) in EigenStepper
1301       double A0 = A[1] * H0[2] - A[2] * H0[1];
1302       double B0 = A[2] * H0[0] - A[0] * H0[2];
1303       double C0 = A[0] * H0[1] - A[1] * H0[0];
1304       // { A2, B2, C2 } is (h/2 * sd.k1 + direction) in EigenStepper
1305       double A2 = A0 + A[0];
1306       double B2 = B0 + A[1];
1307       double C2 = C0 + A[2];
1308       // { A1, B1, C1 } is (h/2 * sd.k1 + 2*direction) in EigenStepper
1309       double A1 = A2 + A[0];
1310       double B1 = B2 + A[1];
1311       double C1 = C2 + A[2];
1312 
1313       // Second point
1314       //
1315       if (!Helix) {
1316         // This is pos1 in EigenStepper
1317         const Vector3 pos(R[0] + A1 * S4, R[1] + B1 * S4, R[2] + C1 * S4);
1318         // This is sd.B_middle in EigenStepper
1319         auto fRes = getField(state, pos);
1320         if (!fRes.ok()) {
1321           return fRes.error();
1322         }
1323         f = *fRes;
1324       } else {
1325         f = f0;
1326       }
1327 
1328       // H1 is (h/(2*momentum) * sd.B_middle) in EigenStepper
1329       double H1[3] = {f[0] * PS2, f[1] * PS2, f[2] * PS2};
1330       // { A3, B3, C3 } is (direction + h/2 * sd.k2) in EigenStepper
1331       double A3 = (A[0] + B2 * H1[2]) - C2 * H1[1];
1332       double B3 = (A[1] + C2 * H1[0]) - A2 * H1[2];
1333       double C3 = (A[2] + A2 * H1[1]) - B2 * H1[0];
1334       // { A4, B4, C4 } is (direction + h/2 * sd.k3) in EigenStepper
1335       double A4 = (A[0] + B3 * H1[2]) - C3 * H1[1];
1336       double B4 = (A[1] + C3 * H1[0]) - A3 * H1[2];
1337       double C4 = (A[2] + A3 * H1[1]) - B3 * H1[0];
1338       // { A5, B5, C5 } is (direction + h * sd.k3) in EigenStepper
1339       double A5 = 2. * A4 - A[0];
1340       double B5 = 2. * B4 - A[1];
1341       double C5 = 2. * C4 - A[2];
1342 
1343       // Last point
1344       //
1345       if (!Helix) {
1346         // This is pos2 in EigenStepper
1347         const Vector3 pos(R[0] + h * A4, R[1] + h * B4, R[2] + h * C4);
1348         // This is sd.B_last in Eigen stepper
1349         auto fRes = getField(state, pos);
1350         if (!fRes.ok()) {
1351           return fRes.error();
1352         }
1353         f = *fRes;
1354       } else {
1355         f = f0;
1356       }
1357 
1358       // H2 is (h/(2*momentum) * sd.B_last) in EigenStepper
1359       double H2[3] = {f[0] * PS2, f[1] * PS2, f[2] * PS2};
1360       // { A6, B6, C6 } is (h/2 * sd.k4) in EigenStepper
1361       double A6 = B5 * H2[2] - C5 * H2[1];
1362       double B6 = C5 * H2[0] - A5 * H2[2];
1363       double C6 = A5 * H2[1] - B5 * H2[0];
1364 
1365       // Test approximation quality on give step and possible step reduction
1366       //
1367       // This is (h2 * (sd.k1 - sd.k2 - sd.k3 + sd.k4).template lpNorm<1>())
1368       // in EigenStepper
1369       EST = 2. * std::abs(h) *
1370             (std::abs((A1 + A6) - (A3 + A4)) + std::abs((B1 + B6) - (B3 + B4)) +
1371              std::abs((C1 + C6) - (C3 + C4)));
1372       EST = std::max(1e-20, EST);
1373       if (!isErrorTolerable(EST)) {
1374         const double stepSizeScaling = calcStepSizeScaling(EST);
1375         h *= stepSizeScaling;
1376         // neutralize the sign of h again
1377         state.stepSize.setAccuracy(h * propDir);
1378         //        dltm = 0.;
1379         continue;
1380       }
1381 
1382       //      if (EST < dltm) h *= 2.;
1383 
1384       // Parameters calculation
1385       //
1386       double A00 = A[0], A11 = A[1], A22 = A[2];
1387 
1388       A[0] = 2. * A3 + (A0 + A5 + A6);
1389       A[1] = 2. * B3 + (B0 + B5 + B6);
1390       A[2] = 2. * C3 + (C0 + C5 + C6);
1391 
1392       double D = (A[0] * A[0] + A[1] * A[1]) + (A[2] * A[2] - 9.);
1393       double Sl = 2. / h;
1394       D = (1. / 3.) - ((1. / 648.) * D) * (12. - D);
1395 
1396       R[0] += (A2 + A3 + A4) * S3;
1397       R[1] += (B2 + B3 + B4) * S3;
1398       R[2] += (C2 + C3 + C4) * S3;
1399       A[0] *= D;
1400       A[1] *= D;
1401       A[2] *= D;
1402       sA[0] = A6 * Sl;
1403       sA[1] = B6 * Sl;
1404       sA[2] = C6 * Sl;
1405 
1406       double mass = particleHypothesis(state).mass();
1407       double momentum = absoluteMomentum(state);
1408 
1409       // Evaluate the time propagation
1410       double dtds = std::sqrt(1 + mass * mass / (momentum * momentum));
1411       state.pVector[3] += h * dtds;
1412       state.pVector[59] = dtds;
1413       state.field = f;
1414       state.newfield = false;
1415 
1416       if (Jac) {
1417         double dtdl = h * mass * mass * qOverP(state) / dtds;
1418         state.pVector[43] += dtdl;
1419 
1420         // Jacobian calculation
1421         //
1422         double* d2A = &state.pVector[28];
1423         double* d3A = &state.pVector[36];
1424         double* d4A = &state.pVector[44];
1425         double d2A0 = H0[2] * d2A[1] - H0[1] * d2A[2];
1426         double d2B0 = H0[0] * d2A[2] - H0[2] * d2A[0];
1427         double d2C0 = H0[1] * d2A[0] - H0[0] * d2A[1];
1428         double d3A0 = H0[2] * d3A[1] - H0[1] * d3A[2];
1429         double d3B0 = H0[0] * d3A[2] - H0[2] * d3A[0];
1430         double d3C0 = H0[1] * d3A[0] - H0[0] * d3A[1];
1431         double d4A0 = (A0 + H0[2] * d4A[1]) - H0[1] * d4A[2];
1432         double d4B0 = (B0 + H0[0] * d4A[2]) - H0[2] * d4A[0];
1433         double d4C0 = (C0 + H0[1] * d4A[0]) - H0[0] * d4A[1];
1434         double d2A2 = d2A0 + d2A[0];
1435         double d2B2 = d2B0 + d2A[1];
1436         double d2C2 = d2C0 + d2A[2];
1437         double d3A2 = d3A0 + d3A[0];
1438         double d3B2 = d3B0 + d3A[1];
1439         double d3C2 = d3C0 + d3A[2];
1440         double d4A2 = d4A0 + d4A[0];
1441         double d4B2 = d4B0 + d4A[1];
1442         double d4C2 = d4C0 + d4A[2];
1443         double d0 = d4A[0] - A00;
1444         double d1 = d4A[1] - A11;
1445         double d2 = d4A[2] - A22;
1446         double d2A3 = (d2A[0] + d2B2 * H1[2]) - d2C2 * H1[1];
1447         double d2B3 = (d2A[1] + d2C2 * H1[0]) - d2A2 * H1[2];
1448         double d2C3 = (d2A[2] + d2A2 * H1[1]) - d2B2 * H1[0];
1449         double d3A3 = (d3A[0] + d3B2 * H1[2]) - d3C2 * H1[1];
1450         double d3B3 = (d3A[1] + d3C2 * H1[0]) - d3A2 * H1[2];
1451         double d3C3 = (d3A[2] + d3A2 * H1[1]) - d3B2 * H1[0];
1452         double d4A3 = ((A3 + d0) + d4B2 * H1[2]) - d4C2 * H1[1];
1453         double d4B3 = ((B3 + d1) + d4C2 * H1[0]) - d4A2 * H1[2];
1454         double d4C3 = ((C3 + d2) + d4A2 * H1[1]) - d4B2 * H1[0];
1455         double d2A4 = (d2A[0] + d2B3 * H1[2]) - d2C3 * H1[1];
1456         double d2B4 = (d2A[1] + d2C3 * H1[0]) - d2A3 * H1[2];
1457         double d2C4 = (d2A[2] + d2A3 * H1[1]) - d2B3 * H1[0];
1458         double d3A4 = (d3A[0] + d3B3 * H1[2]) - d3C3 * H1[1];
1459         double d3B4 = (d3A[1] + d3C3 * H1[0]) - d3A3 * H1[2];
1460         double d3C4 = (d3A[2] + d3A3 * H1[1]) - d3B3 * H1[0];
1461         double d4A4 = ((A4 + d0) + d4B3 * H1[2]) - d4C3 * H1[1];
1462         double d4B4 = ((B4 + d1) + d4C3 * H1[0]) - d4A3 * H1[2];
1463         double d4C4 = ((C4 + d2) + d4A3 * H1[1]) - d4B3 * H1[0];
1464         double d2A5 = 2. * d2A4 - d2A[0];
1465         double d2B5 = 2. * d2B4 - d2A[1];
1466         double d2C5 = 2. * d2C4 - d2A[2];
1467         double d3A5 = 2. * d3A4 - d3A[0];
1468         double d3B5 = 2. * d3B4 - d3A[1];
1469         double d3C5 = 2. * d3C4 - d3A[2];
1470         double d4A5 = 2. * d4A4 - d4A[0];
1471         double d4B5 = 2. * d4B4 - d4A[1];
1472         double d4C5 = 2. * d4C4 - d4A[2];
1473         double d2A6 = d2B5 * H2[2] - d2C5 * H2[1];
1474         double d2B6 = d2C5 * H2[0] - d2A5 * H2[2];
1475         double d2C6 = d2A5 * H2[1] - d2B5 * H2[0];
1476         double d3A6 = d3B5 * H2[2] - d3C5 * H2[1];
1477         double d3B6 = d3C5 * H2[0] - d3A5 * H2[2];
1478         double d3C6 = d3A5 * H2[1] - d3B5 * H2[0];
1479         double d4A6 = d4B5 * H2[2] - d4C5 * H2[1];
1480         double d4B6 = d4C5 * H2[0] - d4A5 * H2[2];
1481         double d4C6 = d4A5 * H2[1] - d4B5 * H2[0];
1482 
1483         double* dR = &state.pVector[24];
1484         dR[0] += (d2A2 + d2A3 + d2A4) * S3;
1485         dR[1] += (d2B2 + d2B3 + d2B4) * S3;
1486         dR[2] += (d2C2 + d2C3 + d2C4) * S3;
1487         d2A[0] = ((d2A0 + 2. * d2A3) + (d2A5 + d2A6)) * (1. / 3.);
1488         d2A[1] = ((d2B0 + 2. * d2B3) + (d2B5 + d2B6)) * (1. / 3.);
1489         d2A[2] = ((d2C0 + 2. * d2C3) + (d2C5 + d2C6)) * (1. / 3.);
1490 
1491         dR = &state.pVector[32];
1492         dR[0] += (d3A2 + d3A3 + d3A4) * S3;
1493         dR[1] += (d3B2 + d3B3 + d3B4) * S3;
1494         dR[2] += (d3C2 + d3C3 + d3C4) * S3;
1495         d3A[0] = ((d3A0 + 2. * d3A3) + (d3A5 + d3A6)) * (1. / 3.);
1496         d3A[1] = ((d3B0 + 2. * d3B3) + (d3B5 + d3B6)) * (1. / 3.);
1497         d3A[2] = ((d3C0 + 2. * d3C3) + (d3C5 + d3C6)) * (1. / 3.);
1498 
1499         dR = &state.pVector[40];
1500         dR[0] += (d4A2 + d4A3 + d4A4) * S3;
1501         dR[1] += (d4B2 + d4B3 + d4B4) * S3;
1502         dR[2] += (d4C2 + d4C3 + d4C4) * S3;
1503         d4A[0] = ((d4A0 + 2. * d4A3) + (d4A5 + d4A6 + A6)) * (1. / 3.);
1504         d4A[1] = ((d4B0 + 2. * d4B3) + (d4B5 + d4B6 + B6)) * (1. / 3.);
1505         d4A[2] = ((d4C0 + 2. * d4C3) + (d4C5 + d4C6 + C6)) * (1. / 3.);
1506       }
1507 
1508       break;
1509     }
1510 
1511     state.pathAccumulated += h;
1512     ++state.nSteps;
1513     state.nStepTrials += nStepTrials;
1514 
1515     const double stepSizeScaling = calcStepSizeScaling(EST);
1516     const double nextAccuracy = std::abs(h * stepSizeScaling);
1517     const double previousAccuracy = std::abs(state.stepSize.accuracy());
1518     const double initialStepLength = std::abs(initialH);
1519     if (nextAccuracy < initialStepLength || nextAccuracy > previousAccuracy) {
1520       state.stepSize.setAccuracy(nextAccuracy);
1521     }
1522 
1523     return h;
1524   }
1525 
1526  private:
1527   std::shared_ptr<const MagneticFieldProvider> m_bField;
1528 
1529   /// Overstep limit: could/should be dynamic
1530   double m_overstepLimit = 100 * UnitConstants::um;
1531 };
1532 
1533 }  // namespace Acts