Warning, file /include/Acts/Propagator/DenseEnvironmentExtension.hpp 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 #pragma once
0010
0011
0012 #include "Acts/Utilities/detail/ReferenceWrapperAnyCompat.hpp"
0013
0014 #include "Acts/Definitions/Algebra.hpp"
0015 #include "Acts/Definitions/PdgParticle.hpp"
0016 #include "Acts/Definitions/TrackParametrization.hpp"
0017 #include "Acts/Material/Interactions.hpp"
0018 #include "Acts/Propagator/AbortList.hpp"
0019 #include "Acts/Propagator/ActionList.hpp"
0020 #include "Acts/Propagator/Propagator.hpp"
0021 #include "Acts/Utilities/VectorHelpers.hpp"
0022
0023 #include <array>
0024 #include <cmath>
0025
0026 namespace Acts {
0027
0028
0029
0030
0031
0032
0033 struct DenseEnvironmentExtension {
0034 using Scalar = ActsScalar;
0035
0036 using ThisVector3 = Eigen::Matrix<Scalar, 3, 1>;
0037
0038
0039 Scalar currentMomentum = 0.;
0040
0041 Scalar initialMomentum = 0.;
0042
0043
0044 Material material;
0045
0046 std::array<Scalar, 4> dLdl{};
0047
0048 std::array<Scalar, 4> qop{};
0049
0050 std::array<Scalar, 4> dPds{};
0051
0052 Scalar dgdqopValue = 0.;
0053
0054 Scalar g = 0.;
0055
0056 std::array<Scalar, 4> tKi{};
0057
0058 std::array<Scalar, 4> Lambdappi{};
0059
0060 std::array<Scalar, 4> energy{};
0061
0062
0063
0064
0065
0066
0067
0068
0069
0070
0071
0072
0073 template <typename propagator_state_t, typename stepper_t,
0074 typename navigator_t>
0075 int bid(const propagator_state_t& state, const stepper_t& stepper,
0076 const navigator_t& navigator) const {
0077 const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
0078 float absQ = particleHypothesis.absoluteCharge();
0079 float mass = particleHypothesis.mass();
0080
0081
0082 if (absQ == 0. || mass == 0. ||
0083 stepper.absoluteMomentum(state.stepping) <
0084 state.options.momentumCutOff) {
0085 return 0;
0086 }
0087
0088
0089 if (!navigator.currentVolumeMaterial(state.navigation)) {
0090 return 0;
0091 }
0092
0093 return 2;
0094 }
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104
0105
0106
0107
0108
0109
0110
0111
0112
0113
0114 template <typename propagator_state_t, typename stepper_t,
0115 typename navigator_t>
0116 bool k(const propagator_state_t& state, const stepper_t& stepper,
0117 const navigator_t& navigator, ThisVector3& knew, const Vector3& bField,
0118 std::array<Scalar, 4>& kQoP, const int i = 0, const double h = 0.,
0119 const ThisVector3& kprev = ThisVector3::Zero()) {
0120 double q = stepper.charge(state.stepping);
0121 const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
0122 float mass = particleHypothesis.mass();
0123
0124
0125 if (i == 0) {
0126
0127 auto volumeMaterial = navigator.currentVolumeMaterial(state.navigation);
0128 ThisVector3 position = stepper.position(state.stepping);
0129 material = volumeMaterial->material(position.template cast<double>());
0130 initialMomentum = stepper.absoluteMomentum(state.stepping);
0131 currentMomentum = initialMomentum;
0132 qop[0] = stepper.qOverP(state.stepping);
0133 initializeEnergyLoss(state, stepper);
0134
0135 knew = qop[0] * stepper.direction(state.stepping).cross(bField);
0136
0137 Lambdappi[0] = -qop[0] * qop[0] * qop[0] * g * energy[0] / (q * q);
0138
0139 tKi[0] = std::hypot(1, mass * qop[0]);
0140 kQoP[0] = Lambdappi[0];
0141 } else {
0142
0143 updateEnergyLoss(mass, h, state, stepper, i);
0144 if (currentMomentum < state.options.momentumCutOff) {
0145 return false;
0146 }
0147
0148 knew = qop[i] *
0149 (stepper.direction(state.stepping) + h * kprev).cross(bField);
0150
0151 auto qopNew = qop[0] + h * Lambdappi[i - 1];
0152 Lambdappi[i] = -qopNew * qopNew * qopNew * g * energy[i] / (q * q);
0153 tKi[i] = std::hypot(1, mass * qopNew);
0154 kQoP[i] = Lambdappi[i];
0155 }
0156 return true;
0157 }
0158
0159
0160
0161
0162
0163
0164
0165
0166
0167
0168
0169
0170
0171
0172
0173 template <typename propagator_state_t, typename stepper_t,
0174 typename navigator_t>
0175 bool finalize(propagator_state_t& state, const stepper_t& stepper,
0176 const navigator_t& , const double h) const {
0177 const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
0178 float mass = particleHypothesis.mass();
0179
0180
0181 auto newMomentum =
0182 stepper.absoluteMomentum(state.stepping) +
0183 (h / 6.) * (dPds[0] + 2. * (dPds[1] + dPds[2]) + dPds[3]);
0184
0185
0186 if (newMomentum < state.options.momentumCutOff) {
0187 return false;
0188 }
0189
0190
0191 state.stepping.derivative(7) = -std::hypot(mass, newMomentum) * g /
0192 (newMomentum * newMomentum * newMomentum);
0193
0194
0195 state.stepping.pars[eFreeQOverP] =
0196 stepper.charge(state.stepping) / newMomentum;
0197
0198 state.stepping.derivative(3) = std::hypot(1, mass / newMomentum);
0199
0200 state.stepping.pars[eFreeTime] +=
0201 (h / 6.) * (tKi[0] + 2. * (tKi[1] + tKi[2]) + tKi[3]);
0202
0203 return true;
0204 }
0205
0206
0207
0208
0209
0210
0211
0212
0213
0214
0215
0216
0217
0218
0219
0220
0221
0222
0223
0224 template <typename propagator_state_t, typename stepper_t,
0225 typename navigator_t>
0226 bool finalize(propagator_state_t& state, const stepper_t& stepper,
0227 const navigator_t& navigator, const double h,
0228 FreeMatrix& D) const {
0229 return finalize(state, stepper, navigator, h) &&
0230 transportMatrix(state, stepper, h, D);
0231 }
0232
0233 private:
0234
0235
0236
0237
0238
0239
0240
0241
0242
0243
0244
0245 template <typename propagator_state_t, typename stepper_t>
0246 bool transportMatrix(propagator_state_t& state, const stepper_t& stepper,
0247 const double h, FreeMatrix& D) const {
0248
0249
0250
0251
0252
0253
0254
0255
0256
0257
0258
0259
0260
0261
0262
0263
0264
0265
0266
0267
0268 auto& sd = state.stepping.stepData;
0269 auto dir = stepper.direction(state.stepping);
0270 const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
0271 float mass = particleHypothesis.mass();
0272
0273 D = FreeMatrix::Identity();
0274 const double half_h = h * 0.5;
0275
0276
0277
0278 auto dFdT = D.block<3, 3>(0, 4);
0279 auto dFdL = D.block<3, 1>(0, 7);
0280
0281 auto dGdT = D.block<3, 3>(4, 4);
0282 auto dGdL = D.block<3, 1>(4, 7);
0283
0284 ActsMatrix<3, 3> dk1dT = ActsMatrix<3, 3>::Zero();
0285 ActsMatrix<3, 3> dk2dT = ActsMatrix<3, 3>::Identity();
0286 ActsMatrix<3, 3> dk3dT = ActsMatrix<3, 3>::Identity();
0287 ActsMatrix<3, 3> dk4dT = ActsMatrix<3, 3>::Identity();
0288
0289 Vector3 dk1dL = Vector3::Zero();
0290 Vector3 dk2dL = Vector3::Zero();
0291 Vector3 dk3dL = Vector3::Zero();
0292 Vector3 dk4dL = Vector3::Zero();
0293
0294
0295 std::array<double, 4> jdL{};
0296
0297
0298 jdL[0] = dLdl[0];
0299 dk1dL = dir.cross(sd.B_first);
0300
0301 jdL[1] = dLdl[1] * (1. + half_h * jdL[0]);
0302 dk2dL = (1. + half_h * jdL[0]) * (dir + half_h * sd.k1).cross(sd.B_middle) +
0303 qop[1] * half_h * dk1dL.cross(sd.B_middle);
0304
0305 jdL[2] = dLdl[2] * (1. + half_h * jdL[1]);
0306 dk3dL = (1. + half_h * jdL[1]) * (dir + half_h * sd.k2).cross(sd.B_middle) +
0307 qop[2] * half_h * dk2dL.cross(sd.B_middle);
0308
0309 jdL[3] = dLdl[3] * (1. + h * jdL[2]);
0310 dk4dL = (1. + h * jdL[2]) * (dir + h * sd.k3).cross(sd.B_last) +
0311 qop[3] * h * dk3dL.cross(sd.B_last);
0312
0313 dk1dT(0, 1) = sd.B_first.z();
0314 dk1dT(0, 2) = -sd.B_first.y();
0315 dk1dT(1, 0) = -sd.B_first.z();
0316 dk1dT(1, 2) = sd.B_first.x();
0317 dk1dT(2, 0) = sd.B_first.y();
0318 dk1dT(2, 1) = -sd.B_first.x();
0319 dk1dT *= qop[0];
0320
0321 dk2dT += half_h * dk1dT;
0322 dk2dT = qop[1] * VectorHelpers::cross(dk2dT, sd.B_middle);
0323
0324 dk3dT += half_h * dk2dT;
0325 dk3dT = qop[2] * VectorHelpers::cross(dk3dT, sd.B_middle);
0326
0327 dk4dT += h * dk3dT;
0328 dk4dT = qop[3] * VectorHelpers::cross(dk4dT, sd.B_last);
0329
0330 dFdT.setIdentity();
0331 dFdT += h / 6. * (dk1dT + dk2dT + dk3dT);
0332 dFdT *= h;
0333
0334 dFdL = h * h / 6. * (dk1dL + dk2dL + dk3dL);
0335
0336 dGdT += h / 6. * (dk1dT + 2. * (dk2dT + dk3dT) + dk4dT);
0337
0338 dGdL = h / 6. * (dk1dL + 2. * (dk2dL + dk3dL) + dk4dL);
0339
0340
0341 D(7, 7) += (h / 6.) * (jdL[0] + 2. * (jdL[1] + jdL[2]) + jdL[3]);
0342
0343
0344
0345
0346
0347
0348
0349
0350 double dtp1dl = qop[0] * mass * mass / std::hypot(1, qop[0] * mass);
0351 double qopNew = qop[0] + half_h * Lambdappi[0];
0352
0353
0354
0355
0356
0357
0358 double dtp2dl = qopNew * mass * mass / std::hypot(1, qopNew * mass);
0359 qopNew = qop[0] + half_h * Lambdappi[1];
0360
0361
0362
0363
0364
0365
0366 double dtp3dl = qopNew * mass * mass / std::hypot(1, qopNew * mass);
0367 qopNew = qop[0] + half_h * Lambdappi[2];
0368 double dtp4dl = qopNew * mass * mass / std::hypot(1, qopNew * mass);
0369
0370
0371
0372
0373
0374 D(3, 7) = (h / 6.) * (dtp1dl + 2. * (dtp2dl + dtp3dl) + dtp4dl);
0375 return true;
0376 }
0377
0378
0379
0380
0381
0382
0383
0384
0385
0386 template <typename propagator_state_t, typename stepper_t>
0387 void initializeEnergyLoss(const propagator_state_t& state,
0388 const stepper_t& stepper) {
0389 const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
0390 float mass = particleHypothesis.mass();
0391 PdgParticle absPdg = particleHypothesis.absolutePdg();
0392 float absQ = particleHypothesis.absoluteCharge();
0393
0394 energy[0] = std::hypot(initialMomentum, mass);
0395
0396 Acts::MaterialSlab slab(material, 1);
0397
0398 if (state.options.meanEnergyLoss) {
0399 g = -computeEnergyLossMean(slab, absPdg, mass, static_cast<float>(qop[0]),
0400 absQ);
0401 } else {
0402
0403
0404 g = -computeEnergyLossMode(slab, absPdg, mass, static_cast<float>(qop[0]),
0405 absQ);
0406 }
0407
0408
0409 dPds[0] = g * energy[0] / initialMomentum;
0410 if (state.stepping.covTransport) {
0411
0412
0413 if (state.options.includeGradient) {
0414 if (state.options.meanEnergyLoss) {
0415 dgdqopValue = deriveEnergyLossMeanQOverP(
0416 slab, absPdg, mass, static_cast<float>(qop[0]), absQ);
0417 } else {
0418
0419 dgdqopValue = deriveEnergyLossModeQOverP(
0420 slab, absPdg, mass, static_cast<float>(qop[0]), absQ);
0421 }
0422 }
0423
0424 dLdl[0] = (-qop[0] * qop[0] * g * energy[0] *
0425 (3. - (initialMomentum * initialMomentum) /
0426 (energy[0] * energy[0])) -
0427 qop[0] * qop[0] * qop[0] * energy[0] * dgdqopValue);
0428 }
0429 }
0430
0431
0432
0433
0434
0435
0436
0437
0438
0439
0440
0441
0442 template <typename propagator_state_t, typename stepper_t>
0443 void updateEnergyLoss(const double mass, const double h,
0444 const propagator_state_t& state,
0445 const stepper_t& stepper, const int i) {
0446
0447 currentMomentum = initialMomentum + h * dPds[i - 1];
0448 energy[i] = std::hypot(currentMomentum, mass);
0449 dPds[i] = g * energy[i] / currentMomentum;
0450 qop[i] = stepper.charge(state.stepping) / currentMomentum;
0451
0452 if (state.stepping.covTransport) {
0453 dLdl[i] = (-qop[i] * qop[i] * g * energy[i] *
0454 (3. - (currentMomentum * currentMomentum) /
0455 (energy[i] * energy[i])) -
0456 qop[i] * qop[i] * qop[i] * energy[i] * dgdqopValue);
0457 }
0458 }
0459 };
0460
0461 template <typename action_list_t = ActionList<>,
0462 typename aborter_list_t = AbortList<>>
0463 struct DenseStepperPropagatorOptions
0464 : public PropagatorOptions<action_list_t, aborter_list_t> {
0465
0466 DenseStepperPropagatorOptions(
0467 const DenseStepperPropagatorOptions<action_list_t, aborter_list_t>&
0468 dspo) = default;
0469
0470
0471
0472
0473
0474 DenseStepperPropagatorOptions(const GeometryContext& gctx,
0475 const MagneticFieldContext& mctx)
0476 : PropagatorOptions<action_list_t, aborter_list_t>(gctx, mctx) {}
0477
0478
0479 bool meanEnergyLoss = true;
0480
0481
0482 bool includeGradient = true;
0483
0484
0485 double momentumCutOff = 0.;
0486
0487
0488
0489
0490
0491
0492 template <typename extended_aborter_list_t>
0493 DenseStepperPropagatorOptions<action_list_t, extended_aborter_list_t> extend(
0494 extended_aborter_list_t aborters) const {
0495 DenseStepperPropagatorOptions<action_list_t, extended_aborter_list_t>
0496 eoptions(this->geoContext, this->magFieldContext);
0497
0498
0499 eoptions.setPlainOptions(*this);
0500
0501
0502 eoptions.actionList = std::move(this->actionList);
0503 eoptions.abortList = std::move(aborters);
0504
0505
0506 eoptions.meanEnergyLoss = meanEnergyLoss;
0507 eoptions.includeGradient = includeGradient;
0508 eoptions.momentumCutOff = momentumCutOff;
0509
0510
0511 return eoptions;
0512 }
0513 };
0514
0515 }