File MultiEigenStepperLoop.hpp

namespace Acts

Note

This file is foreseen for the Geometry module to replace Extent

struct MaxMomentumReducerLoop
#include <Acts/Propagator/MultiEigenStepperLoop.hpp>

Public Static Functions

template<typename stepper_state_t>
static inline ActsScalar charge(const stepper_state_t &s)
template<typename stepper_state_t>
static inline FreeVector cov(const stepper_state_t &s)
template<typename stepper_state_t>
static inline Vector3 direction(const stepper_state_t &s)
template<typename component_range_t>
static inline const auto &maxMomenutmIt(const component_range_t &cmps)
template<typename stepper_state_t>
static inline ActsScalar momentum(const stepper_state_t &s)
template<typename stepper_state_t>
static inline FreeVector pars(const stepper_state_t &s)
template<typename stepper_state_t>
static inline Vector3 position(const stepper_state_t &s)
template<typename stepper_state_t>
static inline ActsScalar time(const stepper_state_t &s)
template<typename extensionlist_t = StepperExtensionList<DefaultExtension>, typename component_reducer_t = WeightedComponentReducerLoop, typename auctioneer_t = detail::VoidAuctioneer>
class MultiEigenStepperLoop : public Acts::EigenStepper<extensionlist_t, auctioneer_t>
#include <Acts/Propagator/MultiEigenStepperLoop.hpp>

Stepper based on the EigenStepper, but handles Multi-Component Tracks (e.g., for the GSF).

Internally, this only manages a vector of EigenStepper::States. This simplifies implementation, but has several drawbacks:

  • There are certain redundancies between the global State and the component states

  • The components do not share a single magnetic-field-cache

Template Parameters
  • extensionlist_t – See EigenStepper for details

  • component_reducer_t – How to map the multi-component state to a single component

  • auctioneer_t – See EigenStepper for details

  • small_vector_size – A size-hint how much memory should be allocated by the small vector

Public Types

using BoundState = std::tuple<BoundTrackParameters, Jacobian, double>
using ConstComponentProxy = ComponentProxyBase<const typename State::Component>

A proxy struct which allows access to a single component of the multi-component state.

It has the semantics of a const reference, i.e. it requires a const reference of the single-component state it represents

using Covariance = BoundSymMatrix
using CurvilinearState = std::tuple<CurvilinearTrackParameters, Jacobian, double>
using Jacobian = BoundMatrix

Jacobian, Covariance and State defintions.

using Reducer = component_reducer_t

The reducer type.

using SingleState = typename SingleStepper::State

Typedef to the State of the single component Stepper.

using SingleStepper = EigenStepper<extensionlist_t, auctioneer_t>

Typedef to the Single-Component Eigen Stepper.

Public Functions

inline MultiEigenStepperLoop(std::shared_ptr<const MagneticFieldProvider> bField)

Constructor from a magnetic field and a optionally provided Logger.

template<typename charge_t>
inline Result<ComponentProxy> addComponent(State &state, const SingleBoundTrackParameters<charge_t> &pars, double weight) const

Add a component to the Multistepper.

Note

This function makes no garantuees about how new components are initialized, it is up to the caller to ensure that all components are valid in the end.

Note

The returned component-proxy is only garantueed to be valid until the component number is again modified

Result<BoundState> boundState(State &state, const Surface &surface, bool transportCov = true, const FreeToBoundCorrection &freeToBoundCorrection = FreeToBoundCorrection(false)) 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.

Note

This is done by combining the gaussian mixture on the specified surface. If the conversion to bound states of some components failes, these components are ignored unless all components fail. In this case an error code is returned.

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

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

  • transportCov[in] Flag steering covariance transport

  • freeToBoundCorrection[in] Flag steering non-linear correction during global to local correction

Returns

A bound state:

  • the parameters at the surface

  • the stepwise jacobian towards it (from last bound)

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

inline double charge(const State &state) const

Charge access.

Parameters

state – [in] The stepping state (thread-local cache)

inline void clearComponents(State &state) const

Reset the number of components.

inline auto componentIterable(State &state) const

Creates an iterable which can be plugged into a range-based for-loop to iterate over components.

Note

Use a for-loop with by-value semantics, since the Iterable returns a proxy internally holding a reference

inline auto constComponentIterable(const State &state) const

Creates an constant iterable which can be plugged into a range-based for-loop to iterate over components.

Note

Use a for-loop with by-value semantics, since the Iterable returns a proxy internally holding a reference

CurvilinearState curvilinearState(State &state, bool transportCov = true) 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.

Note

This is done as a simple average over the free representation and covariance of the components.

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

  • transportCov[in] Flag steering covariance transport

Returns

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)

inline Vector3 direction(const State &state) const

Momentum direction accessor.

Parameters

state – [in] The stepping state (thread-local cache)

inline Result<Vector3> getField(State &state, const Vector3 &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.

Note

This uses the cache of the first component stored in the state

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

  • pos[in] is the field position

inline double getStepSize(const State &state, ConstrainedStep::Type stype) const

Get the step size.

Note

This returns the smalles step size of all components. It uses std::abs for comparison to handle backward propagation and negative step sizes correctly.

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

  • stype – [in] The step size type to be returned

template<typename charge_t>
inline State makeState(std::reference_wrapper<const GeometryContext> gctx, std::reference_wrapper<const MagneticFieldContext> mctx, const MultiComponentBoundTrackParameters<charge_t> &par, NavigationDirection ndir = NavigationDirection::Forward, double ssize = std::numeric_limits<double>::max(), double stolerance = s_onSurfaceTolerance) const

Construct and initialize a state.

inline double momentum(const State &state) const

Absolute momentum accessor.

Parameters

state – [in] The stepping state (thread-local cache)

inline std::size_t numberComponents(const State &state) const

Get the number of components.

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

Output the Step Size of all components into one std::string.

Parameters

state – [in,out] The stepping state (thread-local cache)

inline double overstepLimit(const State &state) const

Overstep limit.

Parameters

state – [in] The stepping state (thread-local cache)

inline Vector3 position(const State &state) const

Global particle position accessor.

Parameters

state – [in] The stepping state (thread-local cache)

inline void releaseStepSize(State &state) const

Release the step-size for all components.

Parameters

state – [in,out] The stepping state (thread-local cache)

inline void removeMissedComponents(State &state) const

Remove missed components from the component state.

inline void resetState(State &state, const BoundVector &boundParams, const BoundSymMatrix &cov, const Surface &surface, const NavigationDirection navDir = NavigationDirection::Forward, const double stepSize = std::numeric_limits<double>::max()) const

Resets the state.

Parameters
  • state[inout] State of the stepper

  • boundParams[in] Parameters in bound parametrisation

  • cov[in] Covariance matrix

  • surface[in] The reference surface of the bound parameters

  • navDir[in] Navigation direction

  • stepSize[in] Step size

inline void setStepSize(State &state, double stepSize, ConstrainedStep::Type stype = ConstrainedStep::actor, bool release = true) 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

  • release – [in] Do we release the step size?

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

Perform a Runge-Kutta track parameter propagation step.

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.

Parameters

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

inline double time(const State &state) const

Time access.

Parameters

state – [in] The stepping state (thread-local cache)

inline void transportCovarianceToBound(State &state, const Surface &surface, const FreeToBoundCorrection &freeToBoundCorrection = FreeToBoundCorrection(false)) 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
  • state[inout] State of the stepper

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

  • freeToBoundCorrection[in] Flag steering non-linear correction during global to local correction to

inline void transportCovarianceToCurvilinear(State &state) const

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

Parameters

state[inout] State of the stepper

template<typename object_intersection_t>
inline 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

inline Intersection3D::Status updateSurfaceStatus(State &state, const Surface &surface, const BoundaryCheck &bcheck, LoggerWrapper logger = getDummyLogger()) 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

  • logger – [in] A LoggerWrapper instance

Public Static Attributes

static constexpr int maxComponents = std::numeric_limits<int>::max()

How many components can this stepper manage?

Private Types

using SmallVector = boost::container::small_vector<T, 16>

Small vector type for speeding up some computations where we need to accumulate stuff of components.

We think 16 is a reasonable amount here.

Private Members

std::size_t m_stepLimitAfterFirstComponentOnSurface = 50

Limits the number of steps after at least one component reached the surface.

struct ComponentProxy : public Acts::MultiEigenStepperLoop<extensionlist_t, component_reducer_t, auctioneer_t>::ComponentProxyBase<State::Component>
#include <Acts/Propagator/MultiEigenStepperLoop.hpp>

A proxy struct which allows access to a single component of the multi-component state.

It has the semantics of a mutable reference, i.e. it requires a mutable reference of the single-component state it represents

Public Types

using Base = ComponentProxyBase<typename State::Component>

Public Functions

inline ComponentProxy(typename State::Component &c, const State &s)
inline Result<BoundState> boundState(const Surface &surface, bool transportCov, const FreeToBoundCorrection &freeToBoundCorrection)
inline auto &charge()
inline auto &cov()
inline auto &derivative()
inline auto &jacobian()
inline auto &jacToGlobal()
inline auto &jacTransport()
inline auto &pars()
inline auto &pathAccumulated()
template<typename propagator_state_t>
inline auto singleState(propagator_state_t &state)
inline auto &status()
inline void update(const FreeVector &freeParams, const BoundVector &boundParams, const Covariance &covariance, const Surface &surface)
inline auto &weight()

Public Members

const State &all_state
template<typename component_t>
struct ComponentProxyBase
#include <Acts/Propagator/MultiEigenStepperLoop.hpp>

A template class which contains all const member functions, that should be available both in the mutable ComponentProxy and the ConstComponentProxy.

Template Parameters

component_t – Must be a const or mutable State::Component.

Public Functions

inline ComponentProxyBase(component_t &c)
inline auto charge() const
inline const auto &cov() const
inline const auto &derivative() const
inline const auto &jacobian() const
inline const auto &jacToGlobal() const
inline const auto &jacTransport() const
inline const auto &pars() const
inline auto pathAccumulated() const
template<typename propagator_state_t>
inline auto singleState(const propagator_state_t &state) const
inline const auto &singleStepper(const MultiEigenStepperLoop &stepper) const
inline auto status() const
inline auto weight() const

Public Members

component_t &cmp
template<typename stepping_t, typename navigation_t, typename options_t, typename geoctx_t>
struct SinglePropState
#include <Acts/Propagator/MultiEigenStepperLoop.hpp>

A helper type for providinig a propagation state which can be used with functions expecting single-component steppers and states.

Public Functions

inline SinglePropState(stepping_t &s, navigation_t &n, options_t &o, geoctx_t &g)

Public Members

geoctx_t &geoContext
navigation_t &navigation
options_t &options
stepping_t &stepping
struct State
#include <Acts/Propagator/MultiEigenStepperLoop.hpp>

Public Functions

State() = delete

No default constructor is provided.

template<typename charge_t>
inline explicit State(const GeometryContext &gctx, const MagneticFieldContext &mctx, const std::shared_ptr<const MagneticFieldProvider> &bfield, const MultiComponentBoundTrackParameters<charge_t> &multipars, NavigationDirection ndir = NavigationDirection::Forward, double ssize = std::numeric_limits<double>::max(), double stolerance = s_onSurfaceTolerance)

Constructor from the initial bound track parameters.

Note

the covariance matrix is copied when needed

Template Parameters

charge_t – Type of the bound parameter charge

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

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

  • bfield[in] the shared magnetic filed provider

  • multipars[in] The track multi-component track-parameters at start

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

  • ssize[in] is the maximum step size

  • stolerance[in] is the stepping tolerance

Public Members

SmallVector<Component> components

The components of which the state consists.

bool covTransport = false
std::reference_wrapper<const GeometryContext> geoContext

geoContext

std::reference_wrapper<const MagneticFieldContext> magContext

MagneticFieldContext.

NavigationDirection navDir
double pathAccumulated = 0.
std::optional<std::size_t> stepCounterAfterFirstComponentOnSurface

Step-limit counter which limits the number of steps when one component reached a surface.

int steps = 0
struct Component
#include <Acts/Propagator/MultiEigenStepperLoop.hpp>

The struct that stores the individual components.

Public Members

SingleState state
Intersection3D::Status status
ActsScalar weight
struct WeightedComponentReducerLoop
#include <Acts/Propagator/MultiEigenStepperLoop.hpp>

Reducer struct for the Loop MultiEigenStepper which reduces the multicomponent state to simply by summing the weighted values.

Public Static Functions

template<typename stepper_state_t>
static inline ActsScalar charge(const stepper_state_t &s)
template<typename stepper_state_t>
static inline FreeVector cov(const stepper_state_t &s)
template<typename stepper_state_t>
static inline Vector3 direction(const stepper_state_t &s)
template<typename stepper_state_t>
static inline ActsScalar momentum(const stepper_state_t &s)
template<typename stepper_state_t>
static inline FreeVector pars(const stepper_state_t &s)
template<typename stepper_state_t>
static inline Vector3 position(const stepper_state_t &s)
template<typename stepper_state_t>
static inline ActsScalar time(const stepper_state_t &s)
template<typename component_range_t>
static inline Vector3 toVector3(const component_range_t &comps, const FreeIndices i)