Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
GenericDenseEnvironmentExtension.hpp
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file GenericDenseEnvironmentExtension.hpp
1 // This file is part of the Acts project.
2 //
3 // Copyright (C) 2018-2019 CERN for the benefit of the Acts project
4 //
5 // This Source Code Form is subject to the terms of the Mozilla Public
6 // License, v. 2.0. If a copy of the MPL was not distributed with this
7 // file, You can obtain one at http://mozilla.org/MPL/2.0/.
8 
9 #pragma once
10 
11 // Workaround for building on clang+libstdc++
13 
19 
20 #include <array>
21 #include <cmath>
22 #include <functional>
23 
24 namespace Acts {
25 namespace detail {
26 
33 template <typename scalar_t>
35  using Scalar = scalar_t;
37  using ThisVector3 = Eigen::Matrix<Scalar, 3, 1>;
38 
47  std::array<Scalar, 4> dLdl{};
49  std::array<Scalar, 4> qop{};
51  std::array<Scalar, 4> dPds{};
55  Scalar g = 0.;
57  std::array<Scalar, 4> tKi{};
59  std::array<Scalar, 4> Lambdappi{};
61  std::array<Scalar, 4> energy{};
62 
65 
77  template <typename propagator_state_t, typename stepper_t,
78  typename navigator_t>
79  int bid(const propagator_state_t& state, const stepper_t& stepper,
80  const navigator_t& navigator) const {
81  const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
82  float absQ = particleHypothesis.absoluteCharge();
83  float mass = particleHypothesis.mass();
84 
85  // Check for valid particle properties
86  if (absQ == 0. || mass == 0. ||
87  stepper.absoluteMomentum(state.stepping) <
88  state.options.momentumCutOff) {
89  return 0;
90  }
91 
92  // Check existence of a volume with material
93  if (!navigator.currentVolumeMaterial(state.navigation)) {
94  return 0;
95  }
96 
97  return 2;
98  }
99 
118  template <typename propagator_state_t, typename stepper_t,
119  typename navigator_t>
120  bool k(const propagator_state_t& state, const stepper_t& stepper,
121  const navigator_t& navigator, ThisVector3& knew, const Vector3& bField,
122  std::array<Scalar, 4>& kQoP, const int i = 0, const double h = 0.,
123  const ThisVector3& kprev = ThisVector3::Zero()) {
124  // using because of autodiff
125  using std::hypot;
126 
127  double q = stepper.charge(state.stepping);
128  const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
129  float mass = particleHypothesis.mass();
130 
131  // i = 0 is used for setup and evaluation of k
132  if (i == 0) {
133  // Set up container for energy loss
134  auto volumeMaterial = navigator.currentVolumeMaterial(state.navigation);
135  ThisVector3 position = stepper.position(state.stepping);
136  material = volumeMaterial->material(position.template cast<double>());
137  initialMomentum = stepper.absoluteMomentum(state.stepping);
139  qop[0] = stepper.qOverP(state.stepping);
140  initializeEnergyLoss(state, stepper);
141  // Evaluate k
142  knew = qop[0] * stepper.direction(state.stepping).cross(bField);
143  // Evaluate k for the time propagation
144  Lambdappi[0] = -qop[0] * qop[0] * qop[0] * g * energy[0] / (q * q);
145  //~ tKi[0] = std::hypot(1, mass / initialMomentum);
146  tKi[0] = hypot(1, mass * qop[0]);
147  kQoP[0] = Lambdappi[0];
148  } else {
149  // Update parameters and check for momentum condition
150  updateEnergyLoss(mass, h, state, stepper, i);
151  if (currentMomentum < state.options.momentumCutOff) {
152  return false;
153  }
154  // Evaluate k
155  knew = qop[i] *
156  (stepper.direction(state.stepping) + h * kprev).cross(bField);
157  // Evaluate k_i for the time propagation
158  auto qopNew = qop[0] + h * Lambdappi[i - 1];
159  Lambdappi[i] = -qopNew * qopNew * qopNew * g * energy[i] / (q * q);
160  tKi[i] = hypot(1, mass * qopNew);
161  kQoP[i] = Lambdappi[i];
162  }
163  return true;
164  }
165 
180  template <typename propagator_state_t, typename stepper_t,
181  typename navigator_t>
182  bool finalize(propagator_state_t& state, const stepper_t& stepper,
183  const navigator_t& /*navigator*/, const double h) const {
184  // using because of autodiff
185  using std::hypot;
186 
187  const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
188  float mass = particleHypothesis.mass();
189 
190  // Evaluate the new momentum
191  auto newMomentum =
192  stepper.absoluteMomentum(state.stepping) +
193  (h / 6.) * (dPds[0] + 2. * (dPds[1] + dPds[2]) + dPds[3]);
194 
195  // Break propagation if momentum becomes below cut-off
196  if (newMomentum < state.options.momentumCutOff) {
197  return false;
198  }
199 
200  // Add derivative dlambda/ds = Lambda''
201  state.stepping.derivative(7) = -hypot(mass, newMomentum) * g /
202  (newMomentum * newMomentum * newMomentum);
203 
204  // Update momentum
205  state.stepping.pars[eFreeQOverP] =
206  stepper.charge(state.stepping) / newMomentum;
207  // Add derivative dt/ds = 1/(beta * c) = sqrt(m^2 * p^{-2} + c^{-2})
208  state.stepping.derivative(3) = hypot(1, mass / newMomentum);
209  // Update time
210  state.stepping.pars[eFreeTime] +=
211  (h / 6.) * (tKi[0] + 2. * (tKi[1] + tKi[2]) + tKi[3]);
212 
213  return true;
214  }
215 
234  template <typename propagator_state_t, typename stepper_t,
235  typename navigator_t>
236  bool finalize(propagator_state_t& state, const stepper_t& stepper,
237  const navigator_t& navigator, const double h,
238  FreeMatrix& D) const {
239  return finalize(state, stepper, navigator, h) &&
240  transportMatrix(state, stepper, navigator, h, D);
241  }
242 
243  private:
255  template <typename propagator_state_t, typename stepper_t,
256  typename navigator_t>
257  bool transportMatrix(propagator_state_t& state, const stepper_t& stepper,
258  const navigator_t& /*navigator*/, const double h,
259  FreeMatrix& D) const {
279 
280  // using because of autodiff
281  using std::hypot;
282 
283  auto& sd = state.stepping.stepData;
284  auto dir = stepper.direction(state.stepping);
285  const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
286  float mass = particleHypothesis.mass();
287 
288  D = FreeMatrix::Identity();
289  const double half_h = h * 0.5;
290 
291  // This sets the reference to the sub matrices
292  // dFdx is already initialised as (3x3) zero
293  auto dFdT = D.block<3, 3>(0, 4);
294  auto dFdL = D.block<3, 1>(0, 7);
295  // dGdx is already initialised as (3x3) identity
296  auto dGdT = D.block<3, 3>(4, 4);
297  auto dGdL = D.block<3, 1>(4, 7);
298 
303 
304  Vector3 dk1dL = Vector3::Zero();
305  Vector3 dk2dL = Vector3::Zero();
306  Vector3 dk3dL = Vector3::Zero();
307  Vector3 dk4dL = Vector3::Zero();
308 
310  std::array<double, 4> jdL{};
311 
312  // Evaluation of the rightmost column without the last term.
313  jdL[0] = dLdl[0];
314  dk1dL = dir.cross(sd.B_first);
315 
316  jdL[1] = dLdl[1] * (1. + half_h * jdL[0]);
317  dk2dL = (1. + half_h * jdL[0]) * (dir + half_h * sd.k1).cross(sd.B_middle) +
318  qop[1] * half_h * dk1dL.cross(sd.B_middle);
319 
320  jdL[2] = dLdl[2] * (1. + half_h * jdL[1]);
321  dk3dL = (1. + half_h * jdL[1]) * (dir + half_h * sd.k2).cross(sd.B_middle) +
322  qop[2] * half_h * dk2dL.cross(sd.B_middle);
323 
324  jdL[3] = dLdl[3] * (1. + h * jdL[2]);
325  dk4dL = (1. + h * jdL[2]) * (dir + h * sd.k3).cross(sd.B_last) +
326  qop[3] * h * dk3dL.cross(sd.B_last);
327 
328  dk1dT(0, 1) = sd.B_first.z();
329  dk1dT(0, 2) = -sd.B_first.y();
330  dk1dT(1, 0) = -sd.B_first.z();
331  dk1dT(1, 2) = sd.B_first.x();
332  dk1dT(2, 0) = sd.B_first.y();
333  dk1dT(2, 1) = -sd.B_first.x();
334  dk1dT *= qop[0];
335 
336  dk2dT += half_h * dk1dT;
337  dk2dT = qop[1] * VectorHelpers::cross(dk2dT, sd.B_middle);
338 
339  dk3dT += half_h * dk2dT;
340  dk3dT = qop[2] * VectorHelpers::cross(dk3dT, sd.B_middle);
341 
342  dk4dT += h * dk3dT;
343  dk4dT = qop[3] * VectorHelpers::cross(dk4dT, sd.B_last);
344 
345  dFdT.setIdentity();
346  dFdT += h / 6. * (dk1dT + dk2dT + dk3dT);
347  dFdT *= h;
348 
349  dFdL = h * h / 6. * (dk1dL + dk2dL + dk3dL);
350 
351  dGdT += h / 6. * (dk1dT + 2. * (dk2dT + dk3dT) + dk4dT);
352 
353  dGdL = h / 6. * (dk1dL + 2. * (dk2dL + dk3dL) + dk4dL);
354 
355  // Evaluation of the dLambda''/dlambda term
356  D(7, 7) += (h / 6.) * (jdL[0] + 2. * (jdL[1] + jdL[2]) + jdL[3]);
357 
358  // The following comment lines refer to the application of the time being
359  // treated as a position. Since t and qop are treated independently for now,
360  // this just serves as entry point for building their relation
361  //~ double dtpp1dl = -mass * mass * qop[0] * qop[0] *
362  //~ (3. * g + qop[0] * dgdqop(energy[0], .mass,
363  //~ absPdg, meanEnergyLoss));
364 
365  double dtp1dl = qop[0] * mass * mass / hypot(1, qop[0] * mass);
366  double qopNew = qop[0] + half_h * Lambdappi[0];
367 
368  //~ double dtpp2dl = -mass * mass * qopNew *
369  //~ qopNew *
370  //~ (3. * g * (1. + half_h * jdL[0]) +
371  //~ qopNew * dgdqop(energy[1], mass, absPdgCode, meanEnergyLoss));
372 
373  double dtp2dl = qopNew * mass * mass / std::hypot(1, qopNew * mass);
374  qopNew = qop[0] + half_h * Lambdappi[1];
375 
376  //~ double dtpp3dl = -mass * mass * qopNew *
377  //~ qopNew *
378  //~ (3. * g * (1. + half_h * jdL[1]) +
379  //~ qopNew * dgdqop(energy[2], mass, absPdg, meanEnergyLoss));
380 
381  double dtp3dl = qopNew * mass * mass / hypot(1, qopNew * mass);
382  qopNew = qop[0] + half_h * Lambdappi[2];
383  double dtp4dl = qopNew * mass * mass / hypot(1, qopNew * mass);
384 
385  //~ D(3, 7) = h * mass * mass * qop[0] /
386  //~ hypot(1., mass * qop[0])
387  //~ + h * h / 6. * (dtpp1dl + dtpp2dl + dtpp3dl);
388 
389  D(3, 7) = (h / 6.) * (dtp1dl + 2. * (dtp2dl + dtp3dl) + dtp4dl);
390  return true;
391  }
392 
400  template <typename propagator_state_t, typename stepper_t>
401  void initializeEnergyLoss(const propagator_state_t& state,
402  const stepper_t& stepper) {
403  // using because of autodiff
404  using std::hypot;
405 
406  const auto& particleHypothesis = stepper.particleHypothesis(state.stepping);
407  float mass = particleHypothesis.mass();
408  PdgParticle absPdg = particleHypothesis.absolutePdg();
409  float absQ = particleHypothesis.absoluteCharge();
410 
411  energy[0] = hypot(initialMomentum, mass);
412  // use unit length as thickness to compute the energy loss per unit length
413  Acts::MaterialSlab slab(material, 1);
414  // Use the same energy loss throughout the step.
415  if (state.options.meanEnergyLoss) {
416  g = -computeEnergyLossMean(slab, absPdg, mass, static_cast<float>(qop[0]),
417  absQ);
418  } else {
419  // TODO using the unit path length is not quite right since the most
420  // probably energy loss is not independent from the path length.
421  g = -computeEnergyLossMode(slab, absPdg, mass, static_cast<float>(qop[0]),
422  absQ);
423  }
424  // Change of the momentum per path length
425  // dPds = dPdE * dEds
426  dPds[0] = g * energy[0] / initialMomentum;
427  if (state.stepping.covTransport) {
428  // Calculate the change of the energy loss per path length and
429  // inverse momentum
430  if (state.options.includeGgradient) {
431  if (state.options.meanEnergyLoss) {
433  slab, absPdg, mass, static_cast<float>(qop[0]), absQ);
434  } else {
435  // TODO path length dependence; see above
437  slab, absPdg, mass, static_cast<float>(qop[0]), absQ);
438  }
439  }
440  // Calculate term for later error propagation
441  dLdl[0] = (-qop[0] * qop[0] * g * energy[0] *
442  (3. - (initialMomentum * initialMomentum) /
443  (energy[0] * energy[0])) -
444  qop[0] * qop[0] * qop[0] * energy[0] * dgdqopValue);
445  }
446  }
447 
457  template <typename propagator_state_t, typename stepper_t>
458  void updateEnergyLoss(const double mass, const double h,
459  const propagator_state_t& state,
460  const stepper_t& stepper, const int i) {
461  // using because of autodiff
462  using std::hypot;
463 
464  // Update parameters related to a changed momentum
465  currentMomentum = initialMomentum + h * dPds[i - 1];
466  energy[i] = hypot(currentMomentum, mass);
467  dPds[i] = g * energy[i] / currentMomentum;
468  qop[i] = stepper.charge(state.stepping) / currentMomentum;
469  // Calculate term for later error propagation
470  if (state.stepping.covTransport) {
471  dLdl[i] = (-qop[i] * qop[i] * g * energy[i] *
472  (3. - (currentMomentum * currentMomentum) /
473  (energy[i] * energy[i])) -
474  qop[i] * qop[i] * qop[i] * energy[i] * dgdqopValue);
475  }
476  }
477 };
478 
479 } // namespace detail
480 } // namespace Acts