File indexing completed on 2025-09-16 08:27:31
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "Acts/Definitions/Direction.hpp"
0010 #include "Acts/EventData/TransformationHelpers.hpp"
0011 #include "Acts/Propagator/ConstrainedStep.hpp"
0012 #include "Acts/Propagator/EigenStepperError.hpp"
0013 #include "Acts/Propagator/detail/CovarianceEngine.hpp"
0014
0015 #include <limits>
0016
0017 template <typename E>
0018 Acts::EigenStepper<E>::EigenStepper(
0019 std::shared_ptr<const MagneticFieldProvider> bField)
0020 : m_bField(std::move(bField)) {}
0021
0022 template <typename E>
0023 auto Acts::EigenStepper<E>::makeState(const Options& options) const -> State {
0024 State state{options, m_bField->makeCache(options.magFieldContext)};
0025 return state;
0026 }
0027
0028 template <typename E>
0029 void Acts::EigenStepper<E>::initialize(State& state,
0030 const BoundTrackParameters& par) const {
0031 initialize(state, par.parameters(), par.covariance(),
0032 par.particleHypothesis(), par.referenceSurface());
0033 }
0034
0035 template <typename E>
0036 void Acts::EigenStepper<E>::initialize(State& state,
0037 const BoundVector& boundParams,
0038 const std::optional<BoundMatrix>& cov,
0039 ParticleHypothesis particleHypothesis,
0040 const Surface& surface) const {
0041 FreeVector freeParams = transformBoundToFreeParameters(
0042 surface, state.options.geoContext, boundParams);
0043
0044 state.particleHypothesis = particleHypothesis;
0045
0046 state.pathAccumulated = 0;
0047 state.nSteps = 0;
0048 state.nStepTrials = 0;
0049 state.stepSize = ConstrainedStep(state.options.maxStepSize);
0050 state.previousStepSize = 0;
0051 state.statistics = StepperStatistics();
0052
0053 state.pars = freeParams;
0054
0055
0056 state.covTransport = cov.has_value();
0057 if (state.covTransport) {
0058 state.cov = *cov;
0059 state.jacToGlobal = surface.boundToFreeJacobian(
0060 state.options.geoContext, freeParams.segment<3>(eFreePos0),
0061 freeParams.segment<3>(eFreeDir0));
0062 state.jacobian = BoundMatrix::Identity();
0063 state.jacTransport = FreeMatrix::Identity();
0064 state.derivative = FreeVector::Zero();
0065 }
0066 }
0067
0068 template <typename E>
0069 auto Acts::EigenStepper<E>::boundState(
0070 State& state, const Surface& surface, bool transportCov,
0071 const FreeToBoundCorrection& freeToBoundCorrection) const
0072 -> Result<BoundState> {
0073 return detail::boundState(
0074 state.options.geoContext, surface, state.cov, state.jacobian,
0075 state.jacTransport, state.derivative, state.jacToGlobal, state.pars,
0076 state.particleHypothesis, state.covTransport && transportCov,
0077 state.pathAccumulated, freeToBoundCorrection);
0078 }
0079
0080 template <typename E>
0081 bool Acts::EigenStepper<E>::prepareCurvilinearState(State& state) const {
0082
0083 if (state.pathAccumulated != 0) {
0084 return true;
0085 }
0086
0087
0088
0089
0090
0091 auto& sd = state.stepData;
0092 auto pos = position(state);
0093 auto fieldRes = getField(state, pos);
0094 if (!fieldRes.ok()) {
0095 return false;
0096 }
0097
0098 sd.B_first = *fieldRes;
0099 if (!state.extension.template k<0>(state, *this, nullptr, sd.k1, sd.B_first,
0100 sd.kQoP)) {
0101 return false;
0102 }
0103
0104
0105 state.derivative.template head<3>() =
0106 state.pars.template segment<3>(eFreeDir0);
0107
0108 state.derivative.template segment<3>(4) = sd.k1;
0109
0110 state.extension.finalize(state, *this, nullptr, state.pathAccumulated);
0111 return true;
0112 }
0113
0114 template <typename E>
0115 auto Acts::EigenStepper<E>::curvilinearState(
0116 State& state, bool transportCov) const -> CurvilinearState {
0117 return detail::curvilinearState(
0118 state.cov, state.jacobian, state.jacTransport, state.derivative,
0119 state.jacToGlobal, state.pars, state.particleHypothesis,
0120 state.covTransport && transportCov, state.pathAccumulated);
0121 }
0122
0123 template <typename E>
0124 void Acts::EigenStepper<E>::update(State& state, const FreeVector& freeParams,
0125 const BoundVector& ,
0126 const Covariance& covariance,
0127 const Surface& surface) const {
0128 state.pars = freeParams;
0129 state.cov = covariance;
0130 state.jacToGlobal = surface.boundToFreeJacobian(
0131 state.options.geoContext, freeParams.template segment<3>(eFreePos0),
0132 freeParams.template segment<3>(eFreeDir0));
0133 }
0134
0135 template <typename E>
0136 void Acts::EigenStepper<E>::update(State& state, const Vector3& uposition,
0137 const Vector3& udirection, double qOverP,
0138 double time) const {
0139 state.pars.template segment<3>(eFreePos0) = uposition;
0140 state.pars.template segment<3>(eFreeDir0) = udirection;
0141 state.pars[eFreeTime] = time;
0142 state.pars[eFreeQOverP] = qOverP;
0143 }
0144
0145 template <typename E>
0146 void Acts::EigenStepper<E>::transportCovarianceToCurvilinear(
0147 State& state) const {
0148 detail::transportCovarianceToCurvilinear(state.cov, state.jacobian,
0149 state.jacTransport, state.derivative,
0150 state.jacToGlobal, direction(state));
0151 }
0152
0153 template <typename E>
0154 void Acts::EigenStepper<E>::transportCovarianceToBound(
0155 State& state, const Surface& surface,
0156 const FreeToBoundCorrection& freeToBoundCorrection) const {
0157 detail::transportCovarianceToBound(
0158 state.options.geoContext, surface, state.cov, state.jacobian,
0159 state.jacTransport, state.derivative, state.jacToGlobal, state.pars,
0160 freeToBoundCorrection);
0161 }
0162
0163 template <typename E>
0164 Acts::Result<double> Acts::EigenStepper<E>::step(
0165 State& state, Direction propDir, const IVolumeMaterial* material) const {
0166
0167 auto& sd = state.stepData;
0168
0169 double errorEstimate = 0;
0170 double h2 = 0;
0171 double half_h = 0;
0172
0173 auto pos = position(state);
0174 auto dir = direction(state);
0175
0176
0177 auto fieldRes = getField(state, pos);
0178 if (!fieldRes.ok()) {
0179 return fieldRes.error();
0180 }
0181 sd.B_first = *fieldRes;
0182 if (!state.extension.template k<0>(state, *this, material, sd.k1, sd.B_first,
0183 sd.kQoP)) {
0184 return 0.;
0185 }
0186
0187 const auto calcStepSizeScaling = [&](const double errorEstimate_) -> double {
0188
0189 constexpr double lower = 0.25;
0190 constexpr double upper = 4.0;
0191
0192 constexpr double exponent = 0.25;
0193
0194 double x = state.options.stepTolerance / errorEstimate_;
0195
0196 if constexpr (exponent == 0.25) {
0197
0198 x = std::sqrt(std::sqrt(x));
0199 } else {
0200 x = std::pow(x, exponent);
0201 }
0202
0203 return std::clamp(x, lower, upper);
0204 };
0205
0206 const auto isErrorTolerable = [&](const double errorEstimate_) {
0207
0208 constexpr double marginFactor = 4.0;
0209
0210 return errorEstimate_ <= marginFactor * state.options.stepTolerance;
0211 };
0212
0213
0214
0215
0216
0217 const auto tryRungeKuttaStep = [&](const double h) -> Result<bool> {
0218
0219 constexpr auto success = &Result<bool>::success;
0220 constexpr auto failure = &Result<bool>::failure;
0221
0222
0223 h2 = h * h;
0224 half_h = h * 0.5;
0225
0226
0227 const Vector3 pos1 = pos + half_h * dir + h2 * 0.125 * sd.k1;
0228 auto field = getField(state, pos1);
0229 if (!field.ok()) {
0230 return failure(field.error());
0231 }
0232 sd.B_middle = *field;
0233
0234 if (!state.extension.template k<1>(state, *this, material, sd.k2,
0235 sd.B_middle, sd.kQoP, half_h, sd.k1)) {
0236 return success(false);
0237 }
0238
0239
0240 if (!state.extension.template k<2>(state, *this, material, sd.k3,
0241 sd.B_middle, sd.kQoP, half_h, sd.k2)) {
0242 return success(false);
0243 }
0244
0245
0246 const Vector3 pos2 = pos + h * dir + h2 * 0.5 * sd.k3;
0247 field = getField(state, pos2);
0248 if (!field.ok()) {
0249 return failure(field.error());
0250 }
0251 sd.B_last = *field;
0252 if (!state.extension.template k<3>(state, *this, material, sd.k4, sd.B_last,
0253 sd.kQoP, h, sd.k3)) {
0254 return success(false);
0255 }
0256
0257
0258 errorEstimate =
0259 h2 * ((sd.k1 - sd.k2 - sd.k3 + sd.k4).template lpNorm<1>() +
0260 std::abs(sd.kQoP[0] - sd.kQoP[1] - sd.kQoP[2] + sd.kQoP[3]));
0261
0262 errorEstimate = std::max(1e-20, errorEstimate);
0263
0264 return success(isErrorTolerable(errorEstimate));
0265 };
0266
0267 const double initialH = state.stepSize.value() * propDir;
0268 double h = initialH;
0269 std::size_t nStepTrials = 0;
0270
0271
0272 while (true) {
0273 ++nStepTrials;
0274 ++state.statistics.nAttemptedSteps;
0275
0276 auto res = tryRungeKuttaStep(h);
0277 if (!res.ok()) {
0278 return res.error();
0279 }
0280 if (!!res.value()) {
0281 break;
0282 }
0283
0284 ++state.statistics.nRejectedSteps;
0285
0286 const double stepSizeScaling = calcStepSizeScaling(errorEstimate);
0287 h *= stepSizeScaling;
0288
0289
0290
0291 if (std::abs(h) < std::abs(state.options.stepSizeCutOff)) {
0292
0293 return EigenStepperError::StepSizeStalled;
0294 }
0295
0296
0297
0298 if (nStepTrials > state.options.maxRungeKuttaStepTrials) {
0299
0300 return EigenStepperError::StepSizeAdjustmentFailed;
0301 }
0302 }
0303
0304
0305 if (state.covTransport) {
0306
0307
0308
0309 FreeMatrix D;
0310 if (!state.extension.finalize(state, *this, material, h, D)) {
0311 return EigenStepperError::StepInvalid;
0312 }
0313
0314
0315
0316
0317
0318
0319
0320
0321
0322
0323
0324
0325
0326
0327
0328
0329
0330
0331 assert((D.topLeftCorner<4, 4>().isIdentity()));
0332 assert((D.bottomLeftCorner<4, 4>().isZero()));
0333 assert((state.jacTransport.template topLeftCorner<4, 4>().isIdentity()));
0334 assert((state.jacTransport.template bottomLeftCorner<4, 4>().isZero()));
0335
0336 state.jacTransport.template topRightCorner<4, 4>() +=
0337 D.topRightCorner<4, 4>() *
0338 state.jacTransport.template bottomRightCorner<4, 4>();
0339 state.jacTransport.template bottomRightCorner<4, 4>() =
0340 (D.bottomRightCorner<4, 4>() *
0341 state.jacTransport.template bottomRightCorner<4, 4>())
0342 .eval();
0343 } else {
0344 if (!state.extension.finalize(state, *this, material, h)) {
0345 return EigenStepperError::StepInvalid;
0346 }
0347 }
0348
0349
0350 state.pars.template segment<3>(eFreePos0) +=
0351 h * dir + h2 / 6. * (sd.k1 + sd.k2 + sd.k3);
0352 state.pars.template segment<3>(eFreeDir0) +=
0353 h / 6. * (sd.k1 + 2. * (sd.k2 + sd.k3) + sd.k4);
0354 (state.pars.template segment<3>(eFreeDir0)).normalize();
0355
0356 if (state.covTransport) {
0357
0358 state.derivative.template head<3>() =
0359 state.pars.template segment<3>(eFreeDir0);
0360 state.derivative.template segment<3>(4) = sd.k4;
0361 }
0362
0363 state.pathAccumulated += h;
0364 ++state.nSteps;
0365 state.nStepTrials += nStepTrials;
0366
0367 ++state.statistics.nSuccessfulSteps;
0368 if (propDir != Direction::fromScalarZeroAsPositive(initialH)) {
0369 ++state.statistics.nReverseSteps;
0370 }
0371 state.statistics.pathLength += h;
0372 state.statistics.absolutePathLength += std::abs(h);
0373
0374 const double stepSizeScaling = calcStepSizeScaling(errorEstimate);
0375 const double nextAccuracy = std::abs(h * stepSizeScaling);
0376 const double previousAccuracy = std::abs(state.stepSize.accuracy());
0377 const double initialStepLength = std::abs(initialH);
0378 if (nextAccuracy < initialStepLength || nextAccuracy > previousAccuracy) {
0379 state.stepSize.setAccuracy(nextAccuracy);
0380 }
0381
0382 return h;
0383 }
0384
0385 template <typename E>
0386 void Acts::EigenStepper<E>::setIdentityJacobian(State& state) const {
0387 state.jacobian = BoundMatrix::Identity();
0388 }