Template Class EigenStepper

Nested Relationships

Class Documentation

template<typename bfield_t, typename extensionlist_t = StepperExtensionList<DefaultExtension>, typename auctioneer_t = detail::VoidAuctioneer>
class Acts::EigenStepper

Runge-Kutta-Nystroem stepper based on Eigen implementation for the following ODE:

r = (x,y,z) … global position T = (Ax,Ay,Az) … momentum direction (normalized)

dr/ds = T dT/ds = q/p * (T x B)

with s being the arc length of the track, q the charge of the particle, p its momentum and B the magnetic field

Public Types

using BField = bfield_t
using BoundState = std::tuple<BoundParameters, Jacobian, double>
using Covariance = BoundSymMatrix
using CurvilinearState = std::tuple<CurvilinearParameters, Jacobian, double>
using Jacobian = BoundMatrix

Jacobian, Covariance and State defintions.

Public Functions

EigenStepper(BField bField)

Constructor requires knowledge of the detector’s magnetic field.

BoundState boundState(State &state, const Surface &surface) const

Create and return the bound state at the current position.

This transports (if necessary) the covariance to the surface and creates a bound state. It does not check if the transported state is at the surface, this needs to be guaranteed by the propagator

Return

A bound state:

  • the parameters at the surface

  • the stepwise jacobian towards it (from last bound)

  • and the path length (from start - for ordering)

Parameters
  • [in] state: State that will be presented as BoundState

  • [in] surface: The surface to which we bind the state

double charge(const State &state) const

Charge access.

Parameters
  • state: [in] The stepping state (thread-local cache)

void covarianceTransport(State &state) const

Method for on-demand transport of the covariance to a new curvilinear frame at current position, or direction of the state.

Return

the full transport jacobian

Parameters
  • [inout] state: State of the stepper

void covarianceTransport(State &state, const Surface &surface) const

Method for on-demand transport of the covariance to a new curvilinear frame at current position, or direction of the state.

Note

no check is done if the position is actually on the surface

Template Parameters
  • surface_t: the Surface type

Parameters
  • [inout] state: State of the stepper

  • [in] surface: is the surface to which the covariance is forwarded to

CurvilinearState curvilinearState(State &state) const

Create and return a curvilinear state at the current position.

This transports (if necessary) the covariance to the current position and creates a curvilinear state.

Return

A curvilinear state:

  • the curvilinear parameters at given position

  • the stepweise jacobian towards it (from last bound)

  • and the path length (from start - for ordering)

Parameters
  • [in] state: State that will be presented as CurvilinearState

Vector3D direction(const State &state) const

Momentum direction accessor.

Parameters
  • state: [in] The stepping state (thread-local cache)

Vector3D getField(State &state, const Vector3D &pos) const

Get the field for the stepping, it checks first if the access is still within the Cell, and updates the cell if necessary.

Parameters
  • [inout] state: is the propagation state associated with the track the magnetic field cell is used (and potentially updated)

  • [in] pos: is the field position

double momentum(const State &state) const

Actual momentum accessor.

Parameters
  • state: [in] The stepping state (thread-local cache)

std::string outputStepSize(const State &state) const

Output the Step Size - single component.

Parameters
  • state: [in,out] The stepping state (thread-local cache)

double overstepLimit(const State&) const

Overstep limit.

Parameters
  • state: [in] The stepping state (thread-local cache)

Vector3D position(const State &state) const

Global particle position accessor.

Parameters
  • state: [in] The stepping state (thread-local cache)

void releaseStepSize(State &state) const

Release the Step size.

Parameters
  • state: [in,out] The stepping state (thread-local cache)

void setStepSize(State &state, double stepSize, ConstrainedStep::Type stype = ConstrainedStep::actor) const

Set Step size - explicitely with a double.

Parameters
  • state: [in,out] The stepping state (thread-local cache)

  • stepSize: [in] The step size value

  • stype: [in] The step size type to be set

template<typename propagator_state_t>
Acts::Result<double> step(propagator_state_t &state) const
template<typename propagator_state_t>
Result<double> step(propagator_state_t &state) const

Perform a Runge-Kutta track parameter propagation step.

Parameters
  • [inout] state: is the propagation state associated with the track parameters that are being propagated.

                     the state contains the desired step size.
                     It can be negative during backwards track
                     propagation,
                     and since we're using an adaptive algorithm, it can
                     be modified by the stepper class during propagation.
    

double time(const State &state) const

Time access.

Parameters
  • state: [in] The stepping state (thread-local cache)

void update(State &state, const BoundParameters &pars) const

Method to update a stepper state to the some parameters.

Parameters
  • [inout] state: State object that will be updated

  • [in] pars: Parameters that will be written into state

void update(State &state, const Vector3D &uposition, const Vector3D &udirection, double up, double time) const

Method to update momentum, direction and p.

Parameters
  • [inout] state: State object that will be updated

  • [in] uposition: the updated position

  • [in] udirection: the updated direction

  • [in] up: the updated momentum value

template<typename object_intersection_t>
void updateStepSize(State &state, const object_intersection_t &oIntersection, bool release = true) const

Update step size.

This method intersects the provided surface and update the navigation step estimation accordingly (hence it changes the state). It also returns the status of the intersection to trigger onSurface in case the surface is reached.

Parameters
  • state: [in,out] The stepping state (thread-local cache)

  • oIntersection: [in] The ObjectIntersection to layer, boundary, etc

  • release: [in] boolean to trigger step size release

Intersection::Status updateSurfaceStatus(State &state, const Surface &surface, const BoundaryCheck &bcheck) const

Update surface status.

It checks the status to the reference surface & updates the step size accordingly

Parameters
  • state: [in,out] The stepping state (thread-local cache)

  • surface: [in] The surface provided

  • bcheck: [in] The boundary check for this status update

struct State

State for track parameter propagation.

It contains the stepping information and is provided thread local by the propagator

Public Functions

State() = delete

Default constructor - deleted.

template<typename parameters_t>
State(std::reference_wrapper<const GeometryContext> gctx, std::reference_wrapper<const MagneticFieldContext> mctx, const parameters_t &par, NavigationDirection ndir = forward, double ssize = std::numeric_limits<double>::max(), double stolerance = s_onSurfaceTolerance)

Constructor from the initial track parameters.

Note

the covariance matrix is copied when needed

Parameters
  • [in] gctx: is the context object for the geometry

  • [in] mctx: is the context object for the magnetic field

  • [in] par: The track parameters at start

  • [in] ndir: The navigation direciton w.r.t momentum

  • [in] ssize: is the maximum step size

  • [in] stolerance: is the stepping tolerance

Public Members

auctioneer_t auctioneer

Auctioneer for choosing the extension.

Vector3D B_first

Magnetic field evaulations.

Vector3D B_last
Vector3D B_middle
Covariance cov = Covariance::Zero()
bool covTransport = false

Covariance matrix (and indicator) / associated with the initial error on track parameters.

FreeVector derivative = FreeVector::Zero()

The propagation derivative.

Vector3D dir = Vector3D(1., 0., 0.)

Momentum direction (normalized)

extensionlist_t extension

List of algorithmic extensions.

BField::Cache fieldCache

This caches the current magnetic field cell and stays (and interpolates) within it as long as this is valid.

See step() code for details.

std::reference_wrapper<const GeometryContext> geoContext

The geometry context.

Jacobian jacobian = Jacobian::Identity()

The full jacobian of the transport entire transport.

BoundToFreeMatrix jacToGlobal = BoundToFreeMatrix::Zero()

Jacobian from local to the global frame.

FreeMatrix jacTransport = FreeMatrix::Identity()

Pure transport jacobian part from runge kutta integration.

Vector3D k1

k_i of the RKN4 algorithm

Vector3D k2
Vector3D k3
Vector3D k4
std::array<double, 4> kQoP

k_i elements of the momenta

NavigationDirection navDir

Navigation direction, this is needed for searching.

double p = 0.

Momentum.

double pathAccumulated = 0.

Accummulated path length state.

Vector3D pos = Vector3D(0., 0., 0.)

Global particle position.

double previousStepSize = 0.

Last performed step (for overstep limit calculation)

double q = 1.

The charge.

struct Acts::EigenStepper::State::[anonymous] stepData

Storage of magnetic field and the sub steps during a RKN4 step.

ConstrainedStep stepSize = {std::numeric_limits<double>::max()}

Adaptive step size of the runge-kutta integration.

double t = 0.

Propagated time.

double tolerance = s_onSurfaceTolerance

The tolerance for the stepping.