Warning, file /include/Acts/Propagator/EigenStepper.ipp was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "Acts/EventData/TransformationHelpers.hpp"
0010 #include "Acts/Propagator/ConstrainedStep.hpp"
0011 #include "Acts/Propagator/detail/CovarianceEngine.hpp"
0012 #include "Acts/Utilities/QuickMath.hpp"
0013
0014 #include <limits>
0015
0016 template <typename E, typename A>
0017 Acts::EigenStepper<E, A>::EigenStepper(
0018 std::shared_ptr<const MagneticFieldProvider> bField,
0019 double )
0020 : m_bField(std::move(bField)) {}
0021
0022 template <typename E, typename A>
0023 auto Acts::EigenStepper<E, A>::makeState(
0024 std::reference_wrapper<const GeometryContext> gctx,
0025 std::reference_wrapper<const MagneticFieldContext> mctx,
0026 const BoundTrackParameters& par, double ssize) const -> State {
0027 return State{gctx, m_bField->makeCache(mctx), par, ssize};
0028 }
0029
0030 template <typename E, typename A>
0031 void Acts::EigenStepper<E, A>::resetState(State& state,
0032 const BoundVector& boundParams,
0033 const BoundSquareMatrix& cov,
0034 const Surface& surface,
0035 const double stepSize) const {
0036 FreeVector freeParams =
0037 transformBoundToFreeParameters(surface, state.geoContext, boundParams);
0038
0039
0040 state.pars = freeParams;
0041 state.cov = cov;
0042 state.stepSize = ConstrainedStep(stepSize);
0043 state.pathAccumulated = 0.;
0044
0045
0046 state.jacToGlobal = surface.boundToFreeJacobian(
0047 state.geoContext, freeParams.template segment<3>(eFreePos0),
0048 freeParams.template segment<3>(eFreeDir0));
0049 state.jacobian = BoundMatrix::Identity();
0050 state.jacTransport = FreeMatrix::Identity();
0051 state.derivative = FreeVector::Zero();
0052 }
0053
0054 template <typename E, typename A>
0055 auto Acts::EigenStepper<E, A>::boundState(
0056 State& state, const Surface& surface, bool transportCov,
0057 const FreeToBoundCorrection& freeToBoundCorrection) const
0058 -> Result<BoundState> {
0059 return detail::boundState(
0060 state.geoContext, surface, state.cov, state.jacobian, state.jacTransport,
0061 state.derivative, state.jacToGlobal, state.pars, state.particleHypothesis,
0062 state.covTransport && transportCov, state.pathAccumulated,
0063 freeToBoundCorrection);
0064 }
0065
0066 template <typename E, typename A>
0067 template <typename propagator_state_t, typename navigator_t>
0068 bool Acts::EigenStepper<E, A>::prepareCurvilinearState(
0069 propagator_state_t& prop_state, const navigator_t& navigator) const {
0070
0071 if (prop_state.stepping.pathAccumulated == 0.) {
0072
0073
0074
0075 if (prop_state.stepping.extension.validExtensionForStep(prop_state, *this,
0076 navigator)) {
0077
0078 auto& sd = prop_state.stepping.stepData;
0079 auto pos = position(prop_state.stepping);
0080 auto fieldRes = getField(prop_state.stepping, pos);
0081 if (fieldRes.ok()) {
0082 sd.B_first = *fieldRes;
0083 if (prop_state.stepping.extension.k1(prop_state, *this, navigator,
0084 sd.k1, sd.B_first, sd.kQoP)) {
0085
0086 prop_state.stepping.derivative.template head<3>() =
0087 prop_state.stepping.pars.template segment<3>(eFreeDir0);
0088
0089 prop_state.stepping.derivative.template segment<3>(4) = sd.k1;
0090
0091 prop_state.stepping.extension.finalize(
0092 prop_state, *this, navigator,
0093 prop_state.stepping.pathAccumulated);
0094 return true;
0095 }
0096 }
0097 }
0098 return false;
0099 }
0100 return true;
0101 }
0102
0103 template <typename E, typename A>
0104 auto Acts::EigenStepper<E, A>::curvilinearState(State& state,
0105 bool transportCov) const
0106 -> CurvilinearState {
0107 return detail::curvilinearState(
0108 state.cov, state.jacobian, state.jacTransport, state.derivative,
0109 state.jacToGlobal, state.pars, state.particleHypothesis,
0110 state.covTransport && transportCov, state.pathAccumulated);
0111 }
0112
0113 template <typename E, typename A>
0114 void Acts::EigenStepper<E, A>::update(State& state,
0115 const FreeVector& freeParams,
0116 const BoundVector& ,
0117 const Covariance& covariance,
0118 const Surface& surface) const {
0119 state.pars = freeParams;
0120 state.cov = covariance;
0121 state.jacToGlobal = surface.boundToFreeJacobian(
0122 state.geoContext, freeParams.template segment<3>(eFreePos0),
0123 freeParams.template segment<3>(eFreeDir0));
0124 }
0125
0126 template <typename E, typename A>
0127 void Acts::EigenStepper<E, A>::update(State& state, const Vector3& uposition,
0128 const Vector3& udirection, double qOverP,
0129 double time) const {
0130 state.pars.template segment<3>(eFreePos0) = uposition;
0131 state.pars.template segment<3>(eFreeDir0) = udirection;
0132 state.pars[eFreeTime] = time;
0133 state.pars[eFreeQOverP] = qOverP;
0134 }
0135
0136 template <typename E, typename A>
0137 void Acts::EigenStepper<E, A>::transportCovarianceToCurvilinear(
0138 State& state) const {
0139 detail::transportCovarianceToCurvilinear(state.cov, state.jacobian,
0140 state.jacTransport, state.derivative,
0141 state.jacToGlobal, direction(state));
0142 }
0143
0144 template <typename E, typename A>
0145 void Acts::EigenStepper<E, A>::transportCovarianceToBound(
0146 State& state, const Surface& surface,
0147 const FreeToBoundCorrection& freeToBoundCorrection) const {
0148 detail::transportCovarianceToBound(state.geoContext.get(), surface, state.cov,
0149 state.jacobian, state.jacTransport,
0150 state.derivative, state.jacToGlobal,
0151 state.pars, freeToBoundCorrection);
0152 }
0153
0154 template <typename E, typename A>
0155 template <typename propagator_state_t, typename navigator_t>
0156 Acts::Result<double> Acts::EigenStepper<E, A>::step(
0157 propagator_state_t& state, const navigator_t& navigator) const {
0158
0159 auto& sd = state.stepping.stepData;
0160
0161 double errorEstimate = 0;
0162 double h2 = 0;
0163 double half_h = 0;
0164
0165 auto pos = position(state.stepping);
0166 auto dir = direction(state.stepping);
0167
0168
0169 auto fieldRes = getField(state.stepping, pos);
0170 if (!fieldRes.ok()) {
0171 return fieldRes.error();
0172 }
0173 sd.B_first = *fieldRes;
0174 if (!state.stepping.extension.validExtensionForStep(state, *this,
0175 navigator) ||
0176 !state.stepping.extension.k1(state, *this, navigator, sd.k1, sd.B_first,
0177 sd.kQoP)) {
0178 return 0.;
0179 }
0180
0181 const auto calcStepSizeScaling = [&](const double errorEstimate_) -> double {
0182
0183 constexpr double lower = 0.25;
0184 constexpr double upper = 4.0;
0185
0186 constexpr double exponent = 0.25;
0187
0188
0189 constexpr bool tryUseFastPow{false};
0190
0191 double x = state.options.stepTolerance / errorEstimate_;
0192
0193 if constexpr (exponent == 0.25 && !tryUseFastPow) {
0194
0195 x = std::sqrt(std::sqrt(x));
0196 } else if constexpr (std::numeric_limits<double>::is_iec559 &&
0197 tryUseFastPow) {
0198 x = fastPow(x, exponent);
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.stepping, pos1);
0229 if (!field.ok()) {
0230 return failure(field.error());
0231 }
0232 sd.B_middle = *field;
0233
0234 if (!state.stepping.extension.k2(state, *this, navigator, sd.k2,
0235 sd.B_middle, sd.kQoP, half_h, sd.k1)) {
0236 return success(false);
0237 }
0238
0239
0240 if (!state.stepping.extension.k3(state, *this, navigator, 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.stepping, pos2);
0248 if (!field.ok()) {
0249 return failure(field.error());
0250 }
0251 sd.B_last = *field;
0252 if (!state.stepping.extension.k4(state, *this, navigator, 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 =
0268 state.stepping.stepSize.value() * state.options.direction;
0269 double h = initialH;
0270 std::size_t nStepTrials = 0;
0271
0272
0273 while (true) {
0274 nStepTrials++;
0275 auto res = tryRungeKuttaStep(h);
0276 if (!res.ok()) {
0277 return res.error();
0278 }
0279 if (!!res.value()) {
0280 break;
0281 }
0282
0283 const double stepSizeScaling = calcStepSizeScaling(errorEstimate);
0284 h *= stepSizeScaling;
0285
0286
0287
0288 if (std::abs(h) < std::abs(state.options.stepSizeCutOff)) {
0289
0290 return EigenStepperError::StepSizeStalled;
0291 }
0292
0293
0294
0295 if (nStepTrials > state.options.maxRungeKuttaStepTrials) {
0296
0297 return EigenStepperError::StepSizeAdjustmentFailed;
0298 }
0299 }
0300
0301
0302 if (state.stepping.covTransport) {
0303
0304
0305
0306 FreeMatrix D;
0307 if (!state.stepping.extension.finalize(state, *this, navigator, h, D)) {
0308 return EigenStepperError::StepInvalid;
0309 }
0310
0311
0312
0313
0314
0315
0316
0317
0318
0319
0320
0321
0322
0323
0324
0325
0326
0327
0328 assert((D.topLeftCorner<4, 4>().isIdentity()));
0329 assert((D.bottomLeftCorner<4, 4>().isZero()));
0330 assert((state.stepping.jacTransport.template topLeftCorner<4, 4>()
0331 .isIdentity()));
0332 assert((state.stepping.jacTransport.template bottomLeftCorner<4, 4>()
0333 .isZero()));
0334
0335 state.stepping.jacTransport.template topRightCorner<4, 4>() +=
0336 D.topRightCorner<4, 4>() *
0337 state.stepping.jacTransport.template bottomRightCorner<4, 4>();
0338 state.stepping.jacTransport.template bottomRightCorner<4, 4>() =
0339 (D.bottomRightCorner<4, 4>() *
0340 state.stepping.jacTransport.template bottomRightCorner<4, 4>())
0341 .eval();
0342 } else {
0343 if (!state.stepping.extension.finalize(state, *this, navigator, h)) {
0344 return EigenStepperError::StepInvalid;
0345 }
0346 }
0347
0348
0349 state.stepping.pars.template segment<3>(eFreePos0) +=
0350 h * dir + h2 / 6. * (sd.k1 + sd.k2 + sd.k3);
0351 state.stepping.pars.template segment<3>(eFreeDir0) +=
0352 h / 6. * (sd.k1 + 2. * (sd.k2 + sd.k3) + sd.k4);
0353 (state.stepping.pars.template segment<3>(eFreeDir0)).normalize();
0354
0355 if (state.stepping.covTransport) {
0356
0357 state.stepping.derivative.template head<3>() =
0358 state.stepping.pars.template segment<3>(eFreeDir0);
0359 state.stepping.derivative.template segment<3>(4) = sd.k4;
0360 }
0361
0362 state.stepping.pathAccumulated += h;
0363 ++state.stepping.nSteps;
0364 state.stepping.nStepTrials += nStepTrials;
0365
0366 const double stepSizeScaling = calcStepSizeScaling(errorEstimate);
0367 const double nextAccuracy = std::abs(h * stepSizeScaling);
0368 const double previousAccuracy = std::abs(state.stepping.stepSize.accuracy());
0369 const double initialStepLength = std::abs(initialH);
0370 if (nextAccuracy < initialStepLength || nextAccuracy > previousAccuracy) {
0371 state.stepping.stepSize.setAccuracy(nextAccuracy);
0372 }
0373
0374 return h;
0375 }
0376
0377 template <typename E, typename A>
0378 void Acts::EigenStepper<E, A>::setIdentityJacobian(State& state) const {
0379 state.jacobian = BoundMatrix::Identity();
0380 }